LCOV - code coverage report
Current view: top level - bias - MovingRestraint.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 101 101 100.0 %
Date: 2024-10-11 08:09:47 Functions: 6 7 85.7 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2011-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             : 
      25             : namespace PLMD {
      26             : namespace bias {
      27             : 
      28             : //+PLUMEDOC BIAS MOVINGRESTRAINT
      29             : /*
      30             : Add a time-dependent, harmonic restraint on one or more variables.
      31             : 
      32             : This form of bias can be used to performed steered MD \cite Grubmuller3
      33             : and Jarzynski sampling \cite jarzynski.
      34             : 
      35             : The harmonic restraint on your system is given by:
      36             : 
      37             : \f[
      38             : V(\vec{s},t) = \frac{1}{2} \kappa(t) ( \vec{s} - \vec{s}_0(t) )^2
      39             : \f]
      40             : 
      41             : The time dependence of \f$\kappa\f$ and \f$\vec{s}_0\f$ are specified by a list of
      42             : STEP, KAPPA and AT keywords.  These keywords tell plumed what values \f$\kappa\f$ and \f$\vec{s}_0\f$
      43             : should have at the time specified by the corresponding STEP keyword.  In between these times
      44             : the values of \f$\kappa\f$ and \f$\vec{s}_0\f$ are linearly interpolated.
      45             : 
      46             : Additional material and examples can be also found in the tutorial \ref belfast-5
      47             : 
      48             : \par Examples
      49             : 
      50             : The following input is dragging the distance between atoms 2 and 4
      51             : from 1 to 2 in the first 1000 steps, then back in the next 1000 steps.
      52             : In the following 500 steps the restraint is progressively switched off.
      53             : \plumedfile
      54             : DISTANCE ATOMS=2,4 LABEL=d
      55             : MOVINGRESTRAINT ...
      56             :   ARG=d
      57             :   STEP0=0    AT0=1.0 KAPPA0=100.0
      58             :   STEP1=1000 AT1=2.0
      59             :   STEP2=2000 AT2=1.0
      60             :   STEP3=2500         KAPPA3=0.0
      61             : ... MOVINGRESTRAINT
      62             : \endplumedfile
      63             : The following input is progressively building restraints
      64             : distances between atoms 1 and 5 and between atoms 2 and 4
      65             : in the first 1000 steps. Afterwards, the restraint is kept
      66             : static.
      67             : \plumedfile
      68             : DISTANCE ATOMS=1,5 LABEL=d1
      69             : DISTANCE ATOMS=2,4 LABEL=d2
      70             : MOVINGRESTRAINT ...
      71             :   ARG=d1,d2
      72             :   STEP0=0    AT0=1.0,1.5 KAPPA0=0.0,0.0
      73             :   STEP1=1000 AT1=1.0,1.5 KAPPA1=1.0,1.0
      74             : ... MOVINGRESTRAINT
      75             : \endplumedfile
      76             : The following input is progressively bringing atoms 1 and 2
      77             : close to each other with an upper wall
      78             : \plumedfile
      79             : DISTANCE ATOMS=1,2 LABEL=d1
      80             : MOVINGRESTRAINT ...
      81             :   ARG=d1
      82             :   VERSE=U
      83             :   STEP0=0    AT0=1.0 KAPPA0=10.0
      84             :   STEP1=1000 AT1=0.0
      85             : ... MOVINGRESTRAINT
      86             : \endplumedfile
      87             : 
      88             : By default the Action is issuing some values which are
      89             : the work on each degree of freedom, the center of the harmonic potential,
      90             : the total bias deposited
      91             : 
      92             : (See also \ref DISTANCE).
      93             : 
      94             : \attention Work is not computed properly when KAPPA is time dependent.
      95             : 
      96             : */
      97             : //+ENDPLUMEDOC
      98             : 
      99             : 
     100             : class MovingRestraint : public Bias {
     101             :   std::vector<std::vector<double> > at;
     102             :   std::vector<std::vector<double> > kappa;
     103             :   std::vector<long int> step;
     104             :   std::vector<double> oldaa;
     105             :   std::vector<double> oldk;
     106             :   std::vector<double> olddpotdk;
     107             :   std::vector<double> oldf;
     108             :   std::vector<std::string> verse;
     109             :   std::vector<double> work;
     110             :   double tot_work;
     111             : public:
     112             :   explicit MovingRestraint(const ActionOptions&);
     113             :   void calculate() override;
     114             :   static void registerKeywords( Keywords& keys );
     115             : };
     116             : 
     117       10427 : PLUMED_REGISTER_ACTION(MovingRestraint,"MOVINGRESTRAINT")
     118             : 
     119           5 : void MovingRestraint::registerKeywords( Keywords& keys ) {
     120           5 :   Bias::registerKeywords(keys);
     121           5 :   keys.use("ARG");
     122          10 :   keys.add("compulsory","VERSE","B","Tells plumed whether the restraint is only acting for CV larger (U) or smaller (L) than "
     123             :            "the restraint or whether it is acting on both sides (B)");
     124          10 :   keys.add("numbered","STEP","This keyword appears multiple times as STEP\\f$x\\f$ with x=0,1,2,...,n. Each value given represents "
     125             :            "the MD step at which the restraint parameters take the values KAPPA\\f$x\\f$ and AT\\f$x\\f$.");
     126          10 :   keys.reset_style("STEP","compulsory");
     127          10 :   keys.add("numbered","AT","AT\\f$x\\f$ is equal to the position of the restraint at time STEP\\f$x\\f$. For intermediate times this parameter "
     128             :            "is linearly interpolated. If no AT\\f$x\\f$ is specified for STEP\\f$x\\f$ then the values of AT are kept constant "
     129             :            "during the interval of time between STEP\\f$x-1\\f$ and STEP\\f$x\\f$.");
     130          10 :   keys.reset_style("AT","compulsory");
     131          10 :   keys.add("numbered","KAPPA","KAPPA\\f$x\\f$ is equal to the value of the force constants at time STEP\\f$x\\f$. For intermediate times this "
     132             :            "parameter is linearly interpolated.  If no KAPPA\\f$x\\f$ is specified for STEP\\f$x\\f$ then the values of KAPPA\\f$x\\f$ "
     133             :            "are kept constant during the interval of time between STEP\\f$x-1\\f$ and STEP\\f$x\\f$.");
     134          10 :   keys.reset_style("KAPPA","compulsory");
     135          10 :   keys.addOutputComponent("work","default","the total work performed changing this restraint");
     136          10 :   keys.addOutputComponent("force2","default","the instantaneous value of the squared force due to this bias potential");
     137          10 :   keys.addOutputComponent("_cntr","default","one or multiple instances of this quantity can be referenced elsewhere in the input file. "
     138             :                           "these quantities will named with  the arguments of the bias followed by "
     139             :                           "the character string _cntr. These quantities give the instantaneous position "
     140             :                           "of the center of the harmonic potential.");
     141          10 :   keys.addOutputComponent("_work","default","one or multiple instances of this quantity can be referenced elsewhere in the input file. "
     142             :                           "These quantities will named with the arguments of the bias followed by "
     143             :                           "the character string _work. These quantities tell the user how much work has "
     144             :                           "been done by the potential in dragging the system along the various colvar axis.");
     145          10 :   keys.addOutputComponent("_kappa","default","one or multiple instances of this quantity can be referenced elsewhere in the input file. "
     146             :                           "These quantities will named with the arguments of the bias followed by "
     147             :                           "the character string _kappa. These quantities tell the user the time dependent value of kappa.");
     148           5 : }
     149             : 
     150           4 : MovingRestraint::MovingRestraint(const ActionOptions&ao):
     151             :   PLUMED_BIAS_INIT(ao),
     152           4 :   verse(getNumberOfArguments())
     153             : {
     154           4 :   parseVector("VERSE",verse);
     155           4 :   std::vector<long int> ss(1); ss[0]=-1;
     156           4 :   std::vector<double> kk( getNumberOfArguments() ), aa( getNumberOfArguments() );
     157          10 :   for(int i=0;; i++) {
     158             :     // Read in step
     159          28 :     if( !parseNumberedVector("STEP",i,ss) ) break;
     160          19 :     for(unsigned j=0; j<step.size(); j++) {
     161           9 :       if(ss[0]<step[j]) error("in moving restraint step number must always increase");
     162             :     }
     163          10 :     step.push_back(ss[0]);
     164             : 
     165             :     // Try to read kappa
     166          20 :     if( !parseNumberedVector("KAPPA",i,kk) ) kk=kappa[i-1];
     167          10 :     kappa.push_back(kk);
     168             : 
     169             :     // Now read AT
     170          20 :     if( !parseNumberedVector("AT",i,aa) ) aa=at[i-1];
     171          10 :     at.push_back(aa);
     172          10 :   }
     173           4 :   checkRead();
     174             : 
     175          14 :   for(unsigned i=0; i<step.size(); i++) {
     176          10 :     log.printf("  step%u %ld\n",i,step[i]);
     177          10 :     log.printf("  at");
     178          22 :     for(unsigned j=0; j<at[i].size(); j++) log.printf(" %f",at[i][j]);
     179          10 :     log.printf("\n");
     180          10 :     log.printf("  with force constant");
     181          22 :     for(unsigned j=0; j<kappa[i].size(); j++) log.printf(" %f",kappa[i][j]);
     182          10 :     log.printf("\n");
     183             :   };
     184             : 
     185          12 :   addComponent("force2"); componentIsNotPeriodic("force2");
     186             : 
     187             :   // add the centers of the restraint as additional components that can be retrieved (useful for debug)
     188             : 
     189             :   std::string comp;
     190           9 :   for(unsigned i=0; i< getNumberOfArguments() ; i++) {
     191           5 :     comp=getPntrToArgument(i)->getName()+"_cntr"; // each spring has its own center
     192           5 :     addComponent(comp); componentIsNotPeriodic(comp);
     193           5 :     comp=getPntrToArgument(i)->getName()+"_work"; // each spring has its own work
     194           5 :     addComponent(comp); componentIsNotPeriodic(comp);
     195           5 :     comp=getPntrToArgument(i)->getName()+"_kappa"; // each spring has its own kappa
     196           5 :     addComponent(comp); componentIsNotPeriodic(comp);
     197           5 :     work.push_back(0.); // initialize the work value
     198             :   }
     199           8 :   addComponent("work"); componentIsNotPeriodic("work");
     200           4 :   tot_work=0.0;
     201             : 
     202           4 :   log<<"  Bibliography ";
     203          12 :   log<<cite("Grubmuller, Heymann, and Tavan, Science 271, 997 (1996)")<<"\n";
     204             : 
     205           4 : }
     206             : 
     207             : 
     208         566 : void MovingRestraint::calculate() {
     209             :   double ene=0.0;
     210             :   double totf2=0.0;
     211         566 :   unsigned narg=getNumberOfArguments();
     212         566 :   long int now=getStep();
     213         566 :   std::vector<double> kk(narg),aa(narg),f(narg),dpotdk(narg);
     214         566 :   if(now<=step[0]) {
     215           3 :     kk=kappa[0];
     216           3 :     aa=at[0];
     217           3 :     oldaa=at[0];
     218           3 :     oldk=kappa[0];
     219           3 :     olddpotdk.resize(narg);
     220           3 :     oldf.resize(narg);
     221         563 :   } else if(now>=step[step.size()-1]) {
     222          47 :     kk=kappa[step.size()-1];
     223          47 :     aa=at[step.size()-1];
     224             :   } else {
     225             :     unsigned i=0;
     226         521 :     for(i=1; i<step.size(); i++) if(now<step[i]) break;
     227         516 :     double c2=(now-step[i-1])/double(step[i]-step[i-1]);
     228         516 :     double c1=1.0-c2;
     229        1040 :     for(unsigned j=0; j<narg; j++) kk[j]=(c1*kappa[i-1][j]+c2*kappa[i][j]);
     230        1040 :     for(unsigned j=0; j<narg; j++) {
     231         524 :       const double a1=at[i-1][j];
     232         524 :       const double a2=at[i][j];
     233         524 :       aa[j]=(c1*a1+c2*(a1+difference(j,a1,a2)));
     234             :     }
     235             :   }
     236         566 :   tot_work=0.0;
     237        1142 :   for(unsigned i=0; i<narg; ++i) {
     238         576 :     const double cv=difference(i,aa[i],getArgument(i)); // this gives: getArgument(i) - aa[i]
     239        1152 :     getPntrToComponent(getPntrToArgument(i)->getName()+"_cntr")->set(aa[i]);
     240         576 :     const double k=kk[i];
     241         576 :     f[i]=-k*cv;
     242         576 :     if(verse[i]=="U" && cv<0) continue;
     243         576 :     if(verse[i]=="L" && cv>0) continue;
     244        1728 :     plumed_assert(verse[i]=="U" || verse[i]=="L" || verse[i]=="B");
     245         576 :     dpotdk[i]=0.5*cv*cv;
     246         576 :     if(oldaa.size()==aa.size() && oldf.size()==f.size()) work[i]+=0.5*(oldf[i]+f[i])*(aa[i]-oldaa[i]) + 0.5*( dpotdk[i]+olddpotdk[i] )*(kk[i]-oldk[i]);
     247        1152 :     getPntrToComponent(getPntrToArgument(i)->getName()+"_work")->set(work[i]);
     248        1152 :     getPntrToComponent(getPntrToArgument(i)->getName()+"_kappa")->set(kk[i]);
     249         576 :     tot_work+=work[i];
     250         576 :     ene+=0.5*k*cv*cv;
     251         576 :     setOutputForce(i,f[i]);
     252         576 :     totf2+=f[i]*f[i];
     253             :   };
     254         566 :   getPntrToComponent("work")->set(tot_work);
     255         566 :   oldf=f;
     256         566 :   oldaa=aa;
     257         566 :   oldk=kk;
     258         566 :   olddpotdk=dpotdk;
     259             :   setBias(ene);
     260        1132 :   getPntrToComponent("force2")->set(totf2);
     261         566 : }
     262             : 
     263             : }
     264             : }
     265             : 
     266             : 

Generated by: LCOV version 1.15