Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2012-2019 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 PATHMSD
31 : /*
32 : This Colvar calculates path collective variables.
33 :
34 : This is the Path Collective Variables implementation
35 : ( see \cite brand07 ).
36 : This variable computes the progress along a given set of frames that is provided
37 : in input ("sss" component) and the distance from them ("zzz" component).
38 : (see below).
39 :
40 : \warning
41 : The molecule used for \ref PATHMSD calculation should be whole (both atoms used in alignment and in displacement calculation).
42 : In case it is broken by the host MD code, please use \ref WHOLEMOLECULES to reconstruct it before \ref PATHMSD calculation.
43 :
44 : \par Examples
45 :
46 : Here below is a case where you have defined three frames and you want to
47 : calculate the progress along the path and the distance from it in p1
48 :
49 : \plumedfile
50 : p1: PATHMSD REFERENCE=file.pdb LAMBDA=500.0 NEIGH_STRIDE=4 NEIGH_SIZE=8
51 : PRINT ARG=p1.sss,p1.zzz STRIDE=1 FILE=colvar FMT=%8.4f
52 : \endplumedfile
53 :
54 : note that NEIGH_STRIDE=4 NEIGH_SIZE=8 control the neighborlist parameter (optional but
55 : recommended for performance) and states that the neighbor list will be calculated every 4
56 : timesteps and consider only the closest 8 member to the actual md snapshots.
57 :
58 : In the REFERENCE PDB file the frames must be separated either using END or ENDMDL.
59 :
60 : \note
61 : The implementation of this collective variable and of \ref PROPERTYMAP
62 : is shared, as well as most input options.
63 :
64 :
65 : */
66 : //+ENDPLUMEDOC
67 :
68 26 : class PathMSD : public PathMSDBase {
69 : public:
70 : explicit PathMSD(const ActionOptions&);
71 : static void registerKeywords(Keywords& keys);
72 : };
73 :
74 6465 : PLUMED_REGISTER_ACTION(PathMSD,"PATHMSD")
75 :
76 14 : void PathMSD::registerKeywords(Keywords& keys) {
77 14 : PathMSDBase::registerKeywords(keys);
78 14 : componentsAreNotOptional(keys);
79 56 : keys.addOutputComponent("sss","default","the position on the path");
80 56 : keys.addOutputComponent("zzz","default","the distance from the path");
81 14 : }
82 :
83 13 : PathMSD::PathMSD(const ActionOptions&ao):
84 13 : Action(ao),PathMSDBase(ao)
85 : {
86 13 : checkRead();
87 :
88 13 : log<<" Bibliography "
89 39 : <<plumed.cite("Branduardi, Gervasio, Parrinello J. Chem. Phys. 126, 054103 (2007)")
90 26 : <<"\n";
91 : // no need to read anything
92 39 : addComponentWithDerivatives("sss"); componentIsNotPeriodic("sss");
93 39 : addComponentWithDerivatives("zzz"); componentIsNotPeriodic("zzz");
94 26 : requestAtoms(pdbv[0].getAtomNumbers());
95 :
96 13 : double i=1.;
97 1153 : for(unsigned it=0 ; it<nframes ; ++it) {
98 570 : vector<double> v; v.push_back(i);
99 570 : indexvec.push_back(v); i+=1.;
100 : }
101 13 : }
102 :
103 : }
104 :
105 4839 : }
|