Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2012-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 "PathMSDBase.h" 23 : #include "core/PlumedMain.h" 24 : 25 : namespace PLMD { 26 : namespace colvar { 27 : 28 : //+PLUMEDOC COLVAR PATHMSD 29 : /* 30 : This Colvar calculates path collective variables. 31 : 32 : This is the Path Collective Variables implementation 33 : ( see \cite brand07 ). 34 : This variable computes the progress along a given set of frames that is provided 35 : in input ("sss" component) and the distance from them ("zzz" component). 36 : (see below). 37 : 38 : When running with periodic boundary conditions, the atoms should be 39 : in the proper periodic image. This is done automatically since PLUMED 2.5, 40 : by considering the ordered list of atoms and rebuilding molecules with a procedure 41 : that is equivalent to that done in \ref WHOLEMOLECULES . Notice that 42 : rebuilding is local to this action. This is different from \ref WHOLEMOLECULES 43 : which actually modifies the coordinates stored in PLUMED. 44 : 45 : In case you want to recover the old behavior you should use the NOPBC flag. 46 : In that case you need to take care that atoms are in the correct 47 : periodic image. 48 : 49 : \par Examples 50 : 51 : Here below is a case where you have defined three frames and you want to 52 : calculate the progress along the path and the distance from it in p1 53 : 54 : \plumedfile 55 : p1: PATHMSD REFERENCE=file.pdb LAMBDA=500.0 NEIGH_STRIDE=4 NEIGH_SIZE=8 56 : PRINT ARG=p1.sss,p1.zzz STRIDE=1 FILE=colvar FMT=%8.4f 57 : \endplumedfile 58 : 59 : note that NEIGH_STRIDE=4 NEIGH_SIZE=8 control the neighbor list parameter (optional but 60 : recommended for performance) and states that the neighbor list will be calculated every 4 61 : steps and consider only the closest 8 member to the actual md snapshots. 62 : 63 : This input must be accompanied by a REFERENCE PDB file in which the positions of each of the frames are specified 64 : separated using either END or ENDMDL as shown below: 65 : 66 : \auxfile{file.pdb} 67 : ATOM 1 CL ALA 1 -3.171 0.295 2.045 1.00 1.00 68 : ATOM 5 CLP ALA 1 -1.819 -0.143 1.679 1.00 1.00 69 : ATOM 6 OL ALA 1 -1.177 -0.889 2.401 1.00 1.00 70 : ATOM 7 NL ALA 1 -1.313 0.341 0.529 1.00 1.00 71 : END 72 : ATOM 1 CL ALA 1 -3.175 0.365 2.024 1.00 1.00 73 : ATOM 5 CLP ALA 1 -1.814 -0.106 1.685 1.00 1.00 74 : ATOM 6 OL ALA 1 -1.201 -0.849 2.425 1.00 1.00 75 : ATOM 7 NL ALA 1 -1.296 0.337 0.534 1.00 1.00 76 : END 77 : ATOM 1 CL ALA 1 -2.990 0.383 2.277 1.00 1.00 78 : ATOM 5 CLP ALA 1 -1.664 -0.085 1.831 1.00 1.00 79 : ATOM 6 OL ALA 1 -0.987 -0.835 2.533 1.00 1.00 80 : ATOM 7 NL ALA 1 -1.227 0.364 0.646 1.00 1.00 81 : END 82 : \endauxfile 83 : 84 : \note 85 : The implementation of this collective variable and of \ref PROPERTYMAP 86 : is shared, as well as most input options. 87 : 88 : 89 : */ 90 : //+ENDPLUMEDOC 91 : 92 : class PathMSD : public PathMSDBase { 93 : public: 94 : explicit PathMSD(const ActionOptions&); 95 : static void registerKeywords(Keywords& keys); 96 : }; 97 : 98 : PLUMED_REGISTER_ACTION(PathMSD,"PATHMSD") 99 : 100 16 : void PathMSD::registerKeywords(Keywords& keys) { 101 16 : PathMSDBase::registerKeywords(keys); 102 32 : keys.addOutputComponent("sss","default","scalar","the position on the path"); 103 32 : keys.addOutputComponent("zzz","default","scalar","the distance from the path"); 104 16 : } 105 : 106 14 : PathMSD::PathMSD(const ActionOptions&ao): 107 14 : Action(ao),PathMSDBase(ao) 108 : { 109 14 : checkRead(); 110 : 111 14 : log<<" Bibliography " 112 28 : <<plumed.cite("Branduardi, Gervasio, Parrinello J. Chem. Phys. 126, 054103 (2007)") 113 28 : <<"\n"; 114 : // no need to read anything 115 42 : addComponentWithDerivatives("sss"); componentIsNotPeriodic("sss"); 116 28 : addComponentWithDerivatives("zzz"); componentIsNotPeriodic("zzz"); 117 14 : requestAtoms(pdbv[0].getAtomNumbers()); 118 : 119 14 : double i=1.; 120 634 : for(unsigned it=0 ; it<nframes ; ++it) { 121 620 : std::vector<double> v; v.push_back(i); 122 620 : indexvec.push_back(v); i+=1.; 123 : } 124 14 : } 125 : 126 : } 127 : 128 : }