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