LCOV - code coverage report
Current view: top level - bias - ExtendedLagrangian.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 88 89 98.9 %
Date: 2024-10-18 13:59:31 Functions: 4 5 80.0 %

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

Generated by: LCOV version 1.16