Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2012-2020 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 "PathMSDBase.h" 23 : #include "core/PlumedMain.h" 24 : 25 : using namespace std; 26 : 27 : namespace PLMD { 28 : namespace colvar { 29 : 30 : //+PLUMEDOC COLVAR PROPERTYMAP 31 : /* 32 : Calculate generic property maps. 33 : 34 : This Colvar calculates the property maps according to the work of Spiwok \cite Spiwok:2011ce. 35 : 36 : 37 : Basically it calculates 38 : \f{eqnarray*}{ 39 : X=\frac{\sum_i X_i*\exp(-\lambda D_i(x))}{\sum_i \exp(-\lambda D_i(x))} \\ 40 : Y=\frac{\sum_i Y_i*\exp(-\lambda D_i(x))}{\sum_i \exp(-\lambda D_i(x))} \\ 41 : \cdots\\ 42 : zzz=-\frac{1}{\lambda}\log(\sum_i \exp(-\lambda D_i(x))) 43 : \f} 44 : 45 : where the parameters \f$X_i\f$ and \f$Y_i\f$ are provided in the input pdb (allv.pdb in this case) and 46 : \f$D_i(x)\f$ is the mean squared displacement after optimal alignment calculated on the pdb frames you input (see Kearsley). 47 : 48 : 49 : When running with periodic boundary conditions, the atoms should be 50 : in the proper periodic image. This is done automatically since PLUMED 2.5, 51 : by considering the ordered list of atoms and rebuilding molecules using a procedure 52 : that is equivalent to that done in \ref WHOLEMOLECULES . Notice that 53 : rebuilding is local to this action. This is different from \ref WHOLEMOLECULES 54 : which actually modifies the coordinates stored in PLUMED. 55 : 56 : In case you want to recover the old behavior you should use the NOPBC flag. 57 : In that case you need to take care that atoms are in the correct 58 : periodic image. 59 : 60 : \par Examples 61 : 62 : \plumedfile 63 : p3: PROPERTYMAP REFERENCE=allv.pdb PROPERTY=X,Y LAMBDA=69087 NEIGH_SIZE=8 NEIGH_STRIDE=4 64 : PRINT ARG=p3.X,p3.Y,p3.zzz STRIDE=1 FILE=colvar FMT=%8.4f 65 : \endplumedfile 66 : 67 : note that NEIGH_STRIDE=4 NEIGH_SIZE=8 control the neighbor list parameter (optional but 68 : recommended for performance) and states that the neighbor list will be calculated every 4 69 : steps and consider only the closest 8 member to the actual md snapshots. 70 : 71 : In this case the input line instructs plumed to look for two properties X and Y with attached values in the REMARK 72 : line of the reference pdb (Note: No spaces from X and = and 1 !!!!). 73 : e.g. 74 : 75 : \auxfile{allv.pdb} 76 : REMARK X=1 Y=2 77 : ATOM 1 CL ALA 1 -3.171 0.295 2.045 1.00 1.00 78 : ATOM 5 CLP ALA 1 -1.819 -0.143 1.679 1.00 1.00 79 : END 80 : REMARK X=2 Y=3 81 : ATOM 1 CL ALA 1 -3.175 0.365 2.024 1.00 1.00 82 : ATOM 5 CLP ALA 1 -1.814 -0.106 1.685 1.00 1.00 83 : END 84 : \endauxfile 85 : 86 : \note 87 : The implementation of this collective variable and of \ref PATHMSD 88 : is shared, as well as most input options. 89 : 90 : */ 91 : //+ENDPLUMEDOC 92 : 93 22 : class PropertyMap : public PathMSDBase { 94 : public: 95 : explicit PropertyMap(const ActionOptions&); 96 : static void registerKeywords(Keywords& keys); 97 : }; 98 : 99 7378 : PLUMED_REGISTER_ACTION(PropertyMap,"PROPERTYMAP") 100 : 101 12 : void PropertyMap::registerKeywords(Keywords& keys) { 102 12 : PathMSDBase::registerKeywords(keys); 103 48 : keys.add("compulsory","PROPERTY","the property to be used in the indexing: this goes in the REMARK field of the reference"); 104 12 : ActionWithValue::useCustomisableComponents(keys); 105 48 : keys.addOutputComponent("zzz","default","the minimum distance from the reference points"); 106 12 : } 107 : 108 11 : PropertyMap::PropertyMap(const ActionOptions&ao): 109 : Action(ao), 110 11 : PathMSDBase(ao) 111 : { 112 : // this is the only additional keyword needed 113 22 : parseVector("PROPERTY",labels); 114 11 : checkRead(); 115 11 : log<<" Bibliography " 116 33 : <<plumed.cite("Spiwok V, Kralova B J. Chem. Phys. 135, 224504 (2011)") 117 22 : <<"\n"; 118 11 : if(labels.size()==0) { 119 : char buf[500]; 120 : sprintf(buf,"Need to specify PROPERTY with this action\n"); 121 0 : plumed_merror(buf); 122 : } else { 123 88 : for(unsigned i=0; i<labels.size(); i++) { 124 44 : log<<" found custom propety to be found in the REMARK line: "<<labels[i].c_str()<<"\n"; 125 66 : addComponentWithDerivatives(labels[i]); componentIsNotPeriodic(labels[i]); 126 : } 127 : // add distance anyhow 128 33 : addComponentWithDerivatives("zzz"); componentIsNotPeriodic("zzz"); 129 : //reparse the REMARK field and pick the index 130 1408 : for(unsigned i=0; i<pdbv.size(); i++) { 131 : // now look for X=1.34555 Y=5.6677 132 : vector<double> labelvals; 133 3696 : for(unsigned j=0; j<labels.size(); j++) { 134 : double val; 135 924 : if( pdbv[i].getArgumentValue(labels[j],val) ) {labelvals.push_back(val);} 136 : else { 137 : char buf[500]; 138 : sprintf(buf,"PROPERTY LABEL \" %s \" NOT FOUND IN REMARK FOR FRAME %u \n",labels[j].c_str(),i); 139 0 : plumed_merror(buf); 140 : }; 141 : } 142 462 : indexvec.push_back(labelvals); 143 : } 144 : } 145 11 : requestAtoms(pdbv[0].getAtomNumbers()); 146 : 147 11 : } 148 : 149 : } 150 5517 : } 151 : 152 :