Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2015-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 "Bias.h"
23 : #include "ActionRegister.h"
24 : #include "tools/Random.h"
25 : #include "core/PlumedMain.h"
26 : #include "core/Atoms.h"
27 :
28 : #include <iostream>
29 :
30 :
31 : using namespace std;
32 :
33 :
34 : namespace PLMD {
35 : namespace bias {
36 :
37 : //+PLUMEDOC BIAS EXTENDED_LAGRANGIAN
38 : /*
39 : Add extended Lagrangian.
40 :
41 : This action can be used to create fictitious collective variables coupled to the real ones.
42 : Given \f$x_i\f$ the i-th argument of this bias potential, potential
43 : and kinetic contributions are added to the energy of the system as
44 : \f[
45 : V=\sum_i \frac{k_i}{2} (x_i-s_i)^2 + \sum_i \frac{\dot{s}_i^2}{2m_i}
46 : \f].
47 :
48 : The resulting potential is thus similar to a \ref RESTRAINT,
49 : but the restraint center moved with time following Hamiltonian
50 : dynamics with mass \f$m_i\f$.
51 :
52 : This bias potential accepts thus vectorial keywords (one element per argument)
53 : to define the coupling constant (KAPPA) and a relaxation time \f$tau\f$ (TAU).
54 : The mass is them computed as \f$m=k(\frac{\tau}{2\pi})^2\f$.
55 :
56 : Notice that this action creates several components.
57 : The ones named XX_fict are the fictitious coordinates. It is possible
58 : to add further forces on them by means of other bias potential,
59 : e.g. to obtain an indirect \ref METAD as in \cite continua .
60 : Also notice that the velocities of the fictitious coordinates
61 : are reported (XX_vfict). However, printed velocities are the ones
62 : at the previous step.
63 :
64 : It is also possible to provide a non-zero friction (one value per component).
65 : This is then used to implement a Langevin thermostat, so as to implement
66 : TAMD/dAFED method \cite Maragliano2006 \cite AbramsJ2008 . Notice that
67 : here a massive Langevin thermostat is used, whereas usually
68 : TAMD employs an overamped Langevin dynamics and dAFED
69 : a Gaussian thermostat.
70 :
71 : \warning
72 : The bias potential is reported in the component bias.
73 : Notice that this bias potential, although formally compatible with
74 : replica exchange framework, probably does not work as expected in that case.
75 : Indeed, since fictitious coordinates are not swapped upon exchange,
76 : acceptace can be expected to be extremely low unless (by chance) two neighboring
77 : replicas have the fictitious variables located properly in space.
78 :
79 : \warning
80 : \ref RESTART is not properly supported by this action. Indeed,
81 : at every start the postion of the fictitious variable is reset to the value
82 : of the real variable, and its velocity is set to zero.
83 : This is not expected to introduce big errors, but certainly is
84 : introducing a small inconsistency between a single long run
85 : and many shorter runs.
86 :
87 : \par Examples
88 :
89 : The following input tells plumed to perform a metadynamics
90 : with an extended Lagrangian on two torsional angles.
91 : \plumedfile
92 : phi: TORSION ATOMS=5,7,9,15
93 : psi: TORSION ATOMS=7,9,15,17
94 : ex: EXTENDED_LAGRANGIAN ARG=phi,psi KAPPA=20,20.0 TAU=0.1,0.1
95 : METAD ARG=ex.phi_fict,ex.psi_fict PACE=100 SIGMA=0.35,0.35 HEIGHT=0.1
96 : # monitor the two variables
97 : PRINT STRIDE=10 ARG=phi,psi,ex.phi_fict,ex.psi_fict FILE=COLVAR
98 : \endplumedfile
99 :
100 : The following input tells plumed to perform a TAMD (or dAFED)
101 : calculation on two torsional angles, keeping the two variables
102 : at a fictitious temperature of 3000K with a Langevin thermostat
103 : with friction 10
104 : \plumedfile
105 : phi: TORSION ATOMS=5,7,9,15
106 : psi: TORSION ATOMS=7,9,15,17
107 : ex: EXTENDED_LAGRANGIAN ARG=phi,psi KAPPA=20,20.0 TAU=0.1,0.1 FRICTION=10,10 TEMP=3000
108 : # monitor the two variables
109 : PRINT STRIDE=10 ARG=phi,psi,ex.phi_fict,ex.psi_fict FILE=COLVAR
110 : \endplumedfile
111 :
112 : */
113 : //+ENDPLUMEDOC
114 :
115 6 : class ExtendedLagrangian : public Bias {
116 : bool firsttime;
117 : std::vector<double> fict;
118 : std::vector<double> vfict;
119 : std::vector<double> vfict_laststep;
120 : std::vector<double> ffict;
121 : std::vector<double> kappa;
122 : std::vector<double> tau;
123 : std::vector<double> friction;
124 : std::vector<Value*> fictValue;
125 : std::vector<Value*> vfictValue;
126 : double kbt;
127 : Random rand;
128 : public:
129 : explicit ExtendedLagrangian(const ActionOptions&);
130 : void calculate();
131 : void update();
132 : static void registerKeywords(Keywords& keys);
133 : };
134 :
135 6454 : PLUMED_REGISTER_ACTION(ExtendedLagrangian,"EXTENDED_LAGRANGIAN")
136 :
137 3 : void ExtendedLagrangian::registerKeywords(Keywords& keys) {
138 3 : Bias::registerKeywords(keys);
139 6 : keys.use("ARG");
140 12 : keys.add("compulsory","KAPPA","specifies that the restraint is harmonic and what the values of the force constants on each of the variables are");
141 12 : keys.add("compulsory","TAU","specifies that the restraint is harmonic and what the values of the force constants on each of the variables are");
142 15 : keys.add("compulsory","FRICTION","0.0","add a friction to the variable");
143 12 : keys.add("optional","TEMP","the system temperature - needed when FRICTION is present. If not provided will be taken from MD code (if available)");
144 3 : componentsAreNotOptional(keys);
145 12 : keys.addOutputComponent("_fict","default","one or multiple instances of this quantity will be refereceable elsewhere in the input file. "
146 : "These quantities will named with the arguments of the bias followed by "
147 : "the character string _tilde. It is possible to add forces on these variable.");
148 12 : keys.addOutputComponent("_vfict","default","one or multiple instances of this quantity will be refereceable elsewhere in the input file. "
149 : "These quantities will named with the arguments of the bias followed by "
150 : "the character string _tilde. It is NOT possible to add forces on these variable.");
151 3 : }
152 :
153 2 : ExtendedLagrangian::ExtendedLagrangian(const ActionOptions&ao):
154 : PLUMED_BIAS_INIT(ao),
155 : firsttime(true),
156 : fict(getNumberOfArguments(),0.0),
157 : vfict(getNumberOfArguments(),0.0),
158 : vfict_laststep(getNumberOfArguments(),0.0),
159 : ffict(getNumberOfArguments(),0.0),
160 : kappa(getNumberOfArguments(),0.0),
161 : tau(getNumberOfArguments(),0.0),
162 : friction(getNumberOfArguments(),0.0),
163 : fictValue(getNumberOfArguments(),NULL),
164 : vfictValue(getNumberOfArguments(),NULL),
165 20 : kbt(0.0)
166 : {
167 4 : parseVector("TAU",tau);
168 4 : parseVector("FRICTION",friction);
169 4 : parseVector("KAPPA",kappa);
170 2 : double temp=-1.0;
171 4 : parse("TEMP",temp);
172 2 : if(temp>=0.0) kbt=plumed.getAtoms().getKBoltzmann()*temp;
173 4 : else kbt=plumed.getAtoms().getKbT();
174 2 : checkRead();
175 :
176 2 : log.printf(" with harmonic force constant");
177 16 : for(unsigned i=0; i<kappa.size(); i++) log.printf(" %f",kappa[i]);
178 2 : log.printf("\n");
179 :
180 2 : log.printf(" with relaxation time");
181 16 : for(unsigned i=0; i<tau.size(); i++) log.printf(" %f",tau[i]);
182 2 : log.printf("\n");
183 :
184 : bool hasFriction=false;
185 10 : for(unsigned i=0; i<getNumberOfArguments(); ++i) if(friction[i]>0.0) hasFriction=true;
186 :
187 2 : if(hasFriction) {
188 2 : log.printf(" with friction");
189 16 : for(unsigned i=0; i<friction.size(); i++) log.printf(" %f",friction[i]);
190 2 : log.printf("\n");
191 : }
192 :
193 2 : log.printf(" and kbt");
194 2 : log.printf(" %f",kbt);
195 2 : log.printf("\n");
196 :
197 10 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
198 4 : std::string comp=getPntrToArgument(i)->getName()+"_fict";
199 4 : addComponentWithDerivatives(comp);
200 4 : if(getPntrToArgument(i)->isPeriodic()) {
201 : std::string a,b;
202 4 : getPntrToArgument(i)->getDomain(a,b);
203 4 : componentIsPeriodic(comp,a,b);
204 0 : } else componentIsNotPeriodic(comp);
205 4 : fictValue[i]=getPntrToComponent(comp);
206 8 : comp=getPntrToArgument(i)->getName()+"_vfict";
207 4 : addComponent(comp);
208 4 : componentIsNotPeriodic(comp);
209 4 : vfictValue[i]=getPntrToComponent(comp);
210 : }
211 :
212 6 : log<<" Bibliography "<<plumed.cite("Iannuzzi, Laio, and Parrinello, Phys. Rev. Lett. 90, 238302 (2003)");
213 2 : if(hasFriction) {
214 6 : log<<plumed.cite("Maragliano and Vanden-Eijnden, Chem. Phys. Lett. 426, 168 (2006)");
215 6 : log<<plumed.cite("Abrams and Tuckerman, J. Phys. Chem. B 112, 15742 (2008)");
216 : }
217 2 : log<<"\n";
218 2 : }
219 :
220 :
221 24 : void ExtendedLagrangian::calculate() {
222 :
223 24 : if(firsttime) {
224 10 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
225 8 : fict[i]=getArgument(i);
226 : }
227 2 : firsttime=false;
228 : }
229 : double ene=0.0;
230 120 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
231 48 : const double cv=difference(i,fict[i],getArgument(i));
232 48 : const double k=kappa[i];
233 48 : const double f=-k*cv;
234 48 : ene+=0.5*k*cv*cv;
235 : setOutputForce(i,f);
236 48 : ffict[i]=-f;
237 : };
238 : setBias(ene);
239 120 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
240 192 : fict[i]=fictValue[i]->bringBackInPbc(fict[i]);
241 96 : fictValue[i]->set(fict[i]);
242 96 : vfictValue[i]->set(vfict_laststep[i]);
243 : }
244 24 : }
245 :
246 24 : void ExtendedLagrangian::update() {
247 24 : double dt=getTimeStep()*getStride();
248 120 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
249 144 : double mass=kappa[i]*tau[i]*tau[i]/(4*pi*pi); // should be k/omega**2
250 48 : double c1=exp(-0.5*friction[i]*dt);
251 48 : double c2=sqrt(kbt*(1.0-c1*c1)/mass);
252 : // consider additional forces on the fictitious particle
253 : // (e.g. MetaD stuff)
254 96 : ffict[i]+=fictValue[i]->getForce();
255 :
256 : // update velocity (half step)
257 96 : vfict[i]+=ffict[i]*0.5*dt/mass;
258 : // thermostat (half step)
259 96 : vfict[i]=c1*vfict[i]+c2*rand.Gaussian();
260 : // save full step velocity to be dumped at next step
261 48 : vfict_laststep[i]=vfict[i];
262 : // thermostat (half step)
263 96 : vfict[i]=c1*vfict[i]+c2*rand.Gaussian();
264 : // update velocity (half step)
265 96 : vfict[i]+=ffict[i]*0.5*dt/mass;
266 : // update position (full step)
267 96 : fict[i]+=vfict[i]*dt;
268 : }
269 24 : }
270 :
271 : }
272 :
273 4839 : }
|