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 : Path collective variables were introduced in the paper cited in the bibliography. 33 : This CV defines a a curvilinear path using a set of reference configurations. Variables 34 : that measure the instantaneous position, sss, of the system on this path and the distance from 35 : this path, zzz, are then computed. The following input illustrates the syntax that is used 36 : with this variable 37 : 38 : ```plumed 39 : #SETTINGS INPUTFILES=regtest/basic/rt39/all.pdb 40 : p1: PATHMSD REFERENCE=regtest/basic/rt39/all.pdb LAMBDA=500.0 NEIGH_STRIDE=4 NEIGH_SIZE=8 41 : PRINT ARG=p1.sss,p1.zzz STRIDE=1 FILE=colvar 42 : ``` 43 : 44 : The `NEIGH_STRIDE=4` and `NEIGH_SIZE=8` keywords here control the neighbor list parameter (optional but 45 : recommended for performance) and states that the neighbor list will be calculated every 4 46 : steps and consider only the closest 8 member to the actual md snapshots. 47 : 48 : When running with periodic boundary conditions, the atoms should be 49 : in the proper periodic image. This is done automatically since PLUMED 2.5, 50 : by considering the ordered list of atoms and rebuilding molecules with a procedure 51 : that is equivalent to that done in [WHOLEMOLECULES](WHOLEMOLECULES.md). Notice that 52 : rebuilding is local to this action. This is different from [WHOLEMOLECULES](WHOLEMOLECULES.md) 53 : which actually modifies the coordinates stored in PLUMED. 54 : 55 : In case you want to recover the old behavior you should use the NOPBC flag. 56 : In that case you need to take care that atoms are in the correct 57 : periodic image. 58 : 59 : The implementation of this collective variable and of [PROPERTYMAP](PROPERTYMAP.md) 60 : is shared, as well as most input options. 61 : 62 : 63 : */ 64 : //+ENDPLUMEDOC 65 : 66 : class PathMSD : public PathMSDBase { 67 : public: 68 : explicit PathMSD(const ActionOptions&); 69 : static void registerKeywords(Keywords& keys); 70 : }; 71 : 72 : PLUMED_REGISTER_ACTION(PathMSD,"PATHMSD") 73 : 74 16 : void PathMSD::registerKeywords(Keywords& keys) { 75 16 : PathMSDBase::registerKeywords(keys); 76 32 : keys.addOutputComponent("sss","default","scalar","the position on the path"); 77 32 : keys.addOutputComponent("zzz","default","scalar","the distance from the path"); 78 16 : keys.addDOI("10.1063/1.2432340"); 79 16 : } 80 : 81 14 : PathMSD::PathMSD(const ActionOptions&ao): 82 14 : Action(ao),PathMSDBase(ao) { 83 14 : checkRead(); 84 : 85 14 : log<<" Bibliography " 86 28 : <<plumed.cite("Branduardi, Gervasio, Parrinello J. Chem. Phys. 126, 054103 (2007)") 87 28 : <<"\n"; 88 : // no need to read anything 89 28 : addComponentWithDerivatives("sss"); 90 28 : componentIsNotPeriodic("sss"); 91 28 : addComponentWithDerivatives("zzz"); 92 14 : componentIsNotPeriodic("zzz"); 93 14 : requestAtoms(pdbv[0].getAtomNumbers()); 94 : 95 14 : double i=1.; 96 634 : for(unsigned it=0 ; it<nframes ; ++it) { 97 : std::vector<double> v; 98 : v.push_back(i); 99 620 : indexvec.push_back(v); 100 620 : i+=1.; 101 : } 102 14 : } 103 : 104 : } 105 : 106 : }