Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2015-2023 The plumed team 3 : (see the PEOPLE file at the root of the distribution for a list of names) 4 : 5 : See http://www.plumed.org for more information. 6 : 7 : This file is part of plumed, version 2. 8 : 9 : plumed is free software: you can redistribute it and/or modify 10 : it under the terms of the GNU Lesser General Public License as published by 11 : the Free Software Foundation, either version 3 of the License, or 12 : (at your option) any later version. 13 : 14 : plumed is distributed in the hope that it will be useful, 15 : but WITHOUT ANY WARRANTY; without even the implied warranty of 16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the 17 : GNU Lesser General Public License for more details. 18 : 19 : You should have received a copy of the GNU Lesser General Public License 20 : along with plumed. If not, see <http://www.gnu.org/licenses/>. 21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */ 22 : #include "LandmarkSelectionBase.h" 23 : #include "core/ActionRegister.h" 24 : #include "tools/Random.h" 25 : #include <iostream> 26 : 27 : //+PLUMEDOC LANDMARKS LANDMARK_SELECT_STAGED 28 : /* 29 : Select a set of landmarks using the staged algorithm. 30 : 31 : \par Examples 32 : 33 : */ 34 : //+ENDPLUMEDOC 35 : 36 : 37 : namespace PLMD { 38 : namespace analysis { 39 : 40 : class LandmarkStaged : public LandmarkSelectionBase { 41 : private: 42 : unsigned seed; 43 : double gamma; 44 : public: 45 : static void registerKeywords( Keywords& keys ); 46 : explicit LandmarkStaged( const ActionOptions& ao ); 47 : void selectLandmarks() override; 48 : }; 49 : 50 10419 : PLUMED_REGISTER_ACTION(LandmarkStaged,"LANDMARK_SELECT_STAGED") 51 : 52 1 : void LandmarkStaged::registerKeywords( Keywords& keys ) { 53 1 : LandmarkSelectionBase::registerKeywords(keys); 54 2 : keys.add("compulsory","GAMMA","the gamma parameter to be used in weights"); 55 2 : keys.add("compulsory","SEED","1234","a random number seed"); 56 1 : } 57 : 58 0 : LandmarkStaged::LandmarkStaged( const ActionOptions& ao ): 59 : Action(ao), 60 0 : LandmarkSelectionBase(ao) 61 : { 62 0 : parse("SEED",seed); parse("GAMMA",gamma); 63 0 : log.printf(" probability of selecting voronoi polyhedra equal to exp(-weight/%f) \n", gamma ); 64 0 : } 65 : 66 0 : void LandmarkStaged::selectLandmarks() { 67 0 : unsigned int n = getNumberOfDataPoints(); // The number of landmarks to pick 68 0 : unsigned int N = my_input_data->getNumberOfDataPoints(); // The total number of frames we can choose from 69 0 : unsigned int m = static_cast<int>( std::sqrt(n*N) ); 70 0 : std::vector<unsigned> fpslandmarks(m); 71 : // Select first point at random 72 0 : Random random; random.setSeed(-seed); double rand=random.RandU01(); 73 0 : fpslandmarks[0] = std::floor( N*rand ); 74 : 75 : // using FPS we want to find m landmarks where m = sqrt(nN) 76 : // Now find distance to all other points 77 : Matrix<double> distances( m, N ); 78 0 : for(unsigned int i=0; i<N; ++i) { 79 0 : distances(0,i) = my_input_data->getDissimilarity( fpslandmarks[0], i ); 80 : } 81 : 82 : // Now find all other landmarks 83 0 : for(unsigned i=1; i<m; ++i) { 84 : // Find point that has the largest minimum distance from the landmarks selected thus far 85 : double maxd=0; 86 0 : for(unsigned j=0; j<N; ++j) { 87 0 : double mind=distances(0,j); 88 0 : for(unsigned k=1; k<i; ++k) { 89 0 : if( distances(k,j)<mind ) { mind=distances(k,j); } 90 : } 91 0 : if( mind>maxd ) { maxd=mind; fpslandmarks[i]=j; } 92 : } 93 0 : for(unsigned k=0; k<N; ++k) distances(i,k) = my_input_data->getDissimilarity( fpslandmarks[i], k ); 94 : } 95 : 96 : // Initial FPS selection of m landmarks completed 97 : // Now find voronoi weights of these m points 98 0 : std::vector<unsigned> poly_assign( N ); 99 0 : std::vector<double> weights( m, 0 ); 100 0 : voronoiAnalysis( fpslandmarks, weights, poly_assign ); 101 : 102 : //Calculate total weight of voronoi polyhedras 103 0 : double vweight=0; for(unsigned i=0; i<m; i++) vweight += std::exp( -weights[i] / gamma ); 104 : 105 0 : std::vector<bool> selected(N, false); unsigned ncount=0; 106 0 : while ( ncount<n) { 107 : // generate random number and check which point it belongs to. select only it was not selected before 108 0 : double rand = vweight*random.RandU01(); 109 : double running_vweight=0; 110 0 : for(unsigned jpoly=0; jpoly<m; ++jpoly) { 111 0 : running_vweight+=std::exp( -weights[jpoly] / gamma ); 112 0 : if( running_vweight>=rand ) { 113 : double tweight=0; 114 0 : for(unsigned i=0; i<poly_assign.size(); ++i) { 115 0 : if( poly_assign[i]==jpoly ) tweight += getWeight( i ); 116 : } 117 0 : double rand_poly = tweight*random.RandU01(); 118 : double running_tweight=0; 119 0 : for(unsigned i=0; i<N; ++i) { 120 0 : if( poly_assign[i]==jpoly ) { 121 0 : running_tweight += getWeight( i ); 122 0 : if( running_tweight>=rand_poly && !selected[i] ) { 123 0 : selectFrame(i); selected[i]=true; ncount++; break; 124 0 : } else if( running_tweight>=rand_poly ) { 125 : break; 126 : } 127 : } 128 : } 129 : } 130 : } 131 : } 132 0 : } 133 : 134 : } 135 : }