LCOV - code coverage report
Current view: top level - opes - OPESexpanded.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 430 442 97.3 %
Date: 2024-10-18 14:00:25 Functions: 12 13 92.3 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2020-2021 of Michele Invernizzi.
       3             : 
       4             :    This file is part of the OPES plumed module.
       5             : 
       6             :    The OPES plumed module is free software: you can redistribute it and/or modify
       7             :    it under the terms of the GNU Lesser General Public License as published by
       8             :    the Free Software Foundation, either version 3 of the License, or
       9             :    (at your option) any later version.
      10             : 
      11             :    The OPES plumed module is distributed in the hope that it will be useful,
      12             :    but WITHOUT ANY WARRANTY; without even the implied warranty of
      13             :    MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      14             :    GNU Lesser General Public License for more details.
      15             : 
      16             :    You should have received a copy of the GNU Lesser General Public License
      17             :    along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      18             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      19             : #include "bias/Bias.h"
      20             : #include "core/PlumedMain.h"
      21             : #include "core/ActionRegister.h"
      22             : #include "core/ActionSet.h"
      23             : #include "tools/Communicator.h"
      24             : #include "tools/File.h"
      25             : #include "tools/OpenMP.h"
      26             : 
      27             : #include "ExpansionCVs.h"
      28             : 
      29             : namespace PLMD {
      30             : namespace opes {
      31             : 
      32             : //+PLUMEDOC OPES_BIAS OPES_EXPANDED
      33             : /*
      34             : On-the-fly probability enhanced sampling with expanded ensembles for the target distribution.
      35             : 
      36             : This method is similar to the OPES method (\ref OPES "OPES") with expanded ensembles target distribution (replica-exchange-like) \cite Invernizzi2020unified.
      37             : 
      38             : An expanded ensemble is obtained by summing a set of ensembles at slightly different termodynamic conditions, or with slightly different Hamiltonians.
      39             : Such ensembles can be sampled via methods like replica exchange, or this \ref OPES_EXPANDED bias action.
      40             : A typical example is a multicanonical simulation, in which a whole range of temperatures is sampled instead of a single one.
      41             : 
      42             : In oreder to define an expanded target ensemble we use \ref EXPANSION_CV "expansion collective variables" (ECVs), \f$\Delta u_\lambda(\mathbf{x})\f$.
      43             : The bias at step \f$n\f$ is
      44             : \f[
      45             :   V_n(\mathbf{x})=-\frac{1}{\beta}\log \left(\frac{1}{N_{\{\lambda\}}}\sum_\lambda e^{-\Delta u_\lambda(\mathbf{x})+\beta\Delta F_n(\lambda)}\right)\, .
      46             : \f]
      47             : See Ref.\cite Invernizzi2020unified for more details on the method.
      48             : 
      49             : Notice that the estimates in the DELTAFS file are expressed in energy units, and should be multiplied by \f$\beta\f$ to be dimensionless as in Ref.\cite Invernizzi2020unified.
      50             : The DELTAFS file also contains an estimate of \f$c(t)=\frac{1}{\beta} \log \langle e^{\beta V}\rangle\f$.
      51             : Similarly to \ref OPES_METAD, it is printed only for reference, since \f$c(t)\f$ reaches a fixed value as the bias converges, and should NOT be used for reweighting.
      52             : Its value is also needed for restarting a simulation.
      53             : 
      54             : You can store the instantaneous \f$\Delta F_n(\lambda)\f$ estimates also in a more readable format using STATE_WFILE and STATE_WSTRIDE.
      55             : Restart can be done either from a DELTAFS file or from a STATE_RFILE, it is equivalent.
      56             : 
      57             : Contrary to \ref OPES_METAD, \ref OPES_EXPANDED does not use kernel density estimation.
      58             : 
      59             : \par Examples
      60             : 
      61             : \plumedfile
      62             : # simulate multiple temperatures, as in parallel tempering
      63             : ene: ENERGY
      64             : ecv: ECV_MULTITHERMAL ARG=ene TEMP_MAX=1000
      65             : opes: OPES_EXPANDED ARG=ecv.* PACE=500
      66             : PRINT FILE=COLVAR STRIDE=500 ARG=ene,opes.bias
      67             : \endplumedfile
      68             : 
      69             : You can easily combine multiple ECVs.
      70             : The \ref OPES_EXPANDED bias will create a multidimensional target grid to sample all the combinations.
      71             : 
      72             : \plumedfile
      73             : # simulate multiple temperatures while biasing a CV
      74             : ene: ENERGY
      75             : dst: DISTANCE ATOMS=1,2
      76             : 
      77             : ecv1: ECV_MULTITHERMAL ARG=ene TEMP_SET_ALL=200,300,500,1000
      78             : ecv2: ECV_UMBRELLAS_LINE ARG=dst CV_MIN=1.2 CV_MAX=4.3 SIGMA=0.5
      79             : opes: OPES_EXPANDED ARG=ecv1.*,ecv2.* PACE=500 OBSERVATION_STEPS=1
      80             : 
      81             : PRINT FILE=COLVAR STRIDE=500 ARG=ene,dst,opes.bias
      82             : \endplumedfile
      83             : 
      84             : If an ECV is based on more than one CV you must provide all the output component, in the proper order.
      85             : You can use \ref Regex for that, like in the following example.
      86             : 
      87             : \plumedfile
      88             : # simulate multiple temperatures and pressures while biasing a two-CVs linear path
      89             : ene: ENERGY
      90             : vol: VOLUME
      91             : ecv_mtp: ECV_MULTITHERMAL_MULTIBARIC ...
      92             :   ARG=ene,vol
      93             :   TEMP=300
      94             :   TEMP_MIN=200
      95             :   TEMP_MAX=800
      96             :   PRESSURE=0.06022140857*1000 #1 kbar
      97             :   PRESSURE_MIN=0
      98             :   PRESSURE_MAX=0.06022140857*2000 #2 kbar
      99             : ...
     100             : 
     101             : cv1: DISTANCE ATOMS=1,2
     102             : cv2: DISTANCE ATOMS=3,4
     103             : ecv_umb: ECV_UMBRELLAS_LINE ARG=cv1,cv2 TEMP=300 CV_MIN=0.1,0.1 CV_MAX=1.5,1.5 SIGMA=0.2 BARRIER=70
     104             : 
     105             : opes: OPES_EXPANDED ARG=(ecv_.*) PACE=500 WALKERS_MPI PRINT_STRIDE=1000
     106             : 
     107             : PRINT FILE=COLVAR STRIDE=500 ARG=ene,vol,cv1,cv2,opes.bias
     108             : \endplumedfile
     109             : 
     110             : 
     111             : */
     112             : //+ENDPLUMEDOC
     113             : 
     114             : class OPESexpanded : public bias::Bias {
     115             : 
     116             : private:
     117             :   bool isFirstStep_;
     118             :   unsigned NumOMP_;
     119             :   unsigned NumParallel_;
     120             :   unsigned rank_;
     121             :   unsigned NumWalkers_;
     122             :   unsigned walker_rank_;
     123             :   unsigned long long counter_;
     124             :   std::size_t ncv_;
     125             : 
     126             :   std::vector<const double *> ECVs_;
     127             :   std::vector<const double *> derECVs_;
     128             :   std::vector<opes::ExpansionCVs*> pntrToECVsClass_;
     129             :   std::vector< std::vector<unsigned> > index_k_;
     130             : // A note on indexes usage:
     131             : //  j -> underlying CVs
     132             : //  i -> DeltaFs
     133             : //  k -> single ECVs, which might not be trivially numbered
     134             : //  l -> groups of ECVs, pntrToECVsClass
     135             : //  h -> subgroups of ECVs, arguments in ECVsClass
     136             : //  w -> walkers
     137             : 
     138             :   double kbt_;
     139             :   unsigned stride_;
     140             :   unsigned deltaF_size_; //different from deltaF_.size() if NumParallel_>1
     141             :   std::vector<double> deltaF_;
     142             :   std::vector<double> diff_;
     143             :   double rct_;
     144             : 
     145             :   std::vector<double> all_deltaF_;
     146             :   std::vector<int> all_size_;
     147             :   std::vector<int> disp_;
     148             : 
     149             :   unsigned obs_steps_;
     150             :   std::vector<double> obs_cvs_;
     151             : 
     152             :   bool calc_work_;
     153             :   double work_;
     154             : 
     155             :   unsigned print_stride_;
     156             :   OFile deltaFsOfile_;
     157             :   std::vector<std::string> deltaF_name_;
     158             : 
     159             :   OFile stateOfile_;
     160             :   int wStateStride_;
     161             :   bool storeOldStates_;
     162             : 
     163             :   void init_pntrToECVsClass();
     164             :   void init_linkECVs();
     165             :   void init_fromObs();
     166             : 
     167             :   void printDeltaF();
     168             :   void dumpStateToFile();
     169             :   void updateDeltaF(double);
     170             :   double getExpansion(const unsigned) const;
     171             : 
     172             : public:
     173             :   explicit OPESexpanded(const ActionOptions&);
     174             :   void calculate() override;
     175             :   void update() override;
     176             :   static void registerKeywords(Keywords& keys);
     177             : };
     178             : 
     179             : PLUMED_REGISTER_ACTION(OPESexpanded,"OPES_EXPANDED")
     180             : 
     181          32 : void OPESexpanded::registerKeywords(Keywords& keys)
     182             : {
     183          32 :   Bias::registerKeywords(keys);
     184          32 :   keys.remove("ARG");
     185          64 :   keys.add("compulsory","ARG","the label of the ECVs that define the expansion. You can use an * to make sure all the output components of the ECVs are used, as in the examples above");
     186          64 :   keys.add("compulsory","PACE","how often the bias is updated");
     187          64 :   keys.add("compulsory","OBSERVATION_STEPS","100","number of unbiased initial PACE steps to collect statistics for initialization");
     188             : //DeltaFs and state files
     189          64 :   keys.add("compulsory","FILE","DELTAFS","a file with the estimate of the relative Delta F for each component of the target and of the global c(t)");
     190          64 :   keys.add("compulsory","PRINT_STRIDE","100","stride for printing to DELTAFS file, in units of PACE");
     191          64 :   keys.add("optional","FMT","specify format for DELTAFS file");
     192          64 :   keys.add("optional","STATE_RFILE","read from this file the Delta F estimates and all the info needed to RESTART the simulation");
     193          64 :   keys.add("optional","STATE_WFILE","write to this file the Delta F estimates and all the info needed to RESTART the simulation");
     194          64 :   keys.add("optional","STATE_WSTRIDE","number of MD steps between writing the STATE_WFILE. Default is only on CPT events (but not all MD codes set them)");
     195          64 :   keys.addFlag("STORE_STATES",false,"append to STATE_WFILE instead of ovewriting it each time");
     196             : //miscellaneous
     197          64 :   keys.addFlag("CALC_WORK",false,"calculate the total accumulated work done by the bias since last restart");
     198          64 :   keys.addFlag("WALKERS_MPI",false,"switch on MPI version of multiple walkers");
     199          64 :   keys.addFlag("SERIAL",false,"perform calculations in serial");
     200          32 :   keys.use("RESTART");
     201          32 :   keys.use("UPDATE_FROM");
     202          32 :   keys.use("UPDATE_UNTIL");
     203             : 
     204             : //output components
     205          64 :   keys.addOutputComponent("work","CALC_WORK","total accumulated work done by the bias");
     206          32 : }
     207             : 
     208          30 : OPESexpanded::OPESexpanded(const ActionOptions&ao)
     209             :   : PLUMED_BIAS_INIT(ao)
     210          30 :   , isFirstStep_(true)
     211          30 :   , counter_(0)
     212          30 :   , ncv_(getNumberOfArguments())
     213          30 :   , deltaF_size_(0)
     214          30 :   , rct_(0)
     215          30 :   , work_(0)
     216             : {
     217             : //set pace
     218          30 :   parse("PACE",stride_);
     219          30 :   parse("OBSERVATION_STEPS",obs_steps_);
     220          30 :   plumed_massert(obs_steps_!=0,"minimum is OBSERVATION_STEPS=1");
     221          30 :   obs_cvs_.resize(obs_steps_*ncv_);
     222             : 
     223             : //deltaFs file
     224             :   std::string deltaFsFileName;
     225          30 :   parse("FILE",deltaFsFileName);
     226          60 :   parse("PRINT_STRIDE",print_stride_);
     227             :   std::string fmt;
     228          60 :   parse("FMT",fmt);
     229             : //output checkpoint of current state
     230             :   std::string restartFileName;
     231          60 :   parse("STATE_RFILE",restartFileName);
     232             :   std::string stateFileName;
     233          30 :   parse("STATE_WFILE",stateFileName);
     234          30 :   wStateStride_=0;
     235          30 :   parse("STATE_WSTRIDE",wStateStride_);
     236          30 :   storeOldStates_=false;
     237          30 :   parseFlag("STORE_STATES",storeOldStates_);
     238          30 :   if(wStateStride_!=0 || storeOldStates_)
     239           5 :     plumed_massert(stateFileName.length()>0,"filename for storing simulation status not specified, use STATE_WFILE");
     240          30 :   if(wStateStride_>0)
     241           5 :     plumed_massert(wStateStride_>=(int)stride_,"STATE_WSTRIDE is in units of MD steps, thus should be a multiple of PACE");
     242          30 :   if(stateFileName.length()>0 && wStateStride_==0)
     243           1 :     wStateStride_=-1;//will print only on CPT events (checkpoints set by some MD engines, like gromacs)
     244             : 
     245             : //work flag
     246          30 :   parseFlag("CALC_WORK",calc_work_);
     247             : 
     248             : //multiple walkers //external MW for cp2k not supported, but anyway cp2k cannot put bias on energy!
     249          30 :   bool walkers_mpi=false;
     250          30 :   parseFlag("WALKERS_MPI",walkers_mpi);
     251          30 :   if(walkers_mpi)
     252             :   {
     253             :     //If this Action is not compiled with MPI the user is informed and we exit gracefully
     254           4 :     plumed_massert(Communicator::plumedHasMPI(),"Invalid walkers configuration: WALKERS_MPI flag requires MPI compilation");
     255           4 :     plumed_massert(Communicator::initialized(),"Invalid walkers configuration: WALKERS_MPI needs the communicator correctly initialized.");
     256             : 
     257           4 :     if(comm.Get_rank()==0) //multi_sim_comm works on first rank only
     258             :     {
     259           4 :       NumWalkers_=multi_sim_comm.Get_size();
     260           4 :       walker_rank_=multi_sim_comm.Get_rank();
     261             :     }
     262           4 :     comm.Bcast(NumWalkers_,0); //if each walker has more than one processor update them all
     263           4 :     comm.Bcast(walker_rank_,0);
     264             :   }
     265             :   else
     266             :   {
     267          26 :     NumWalkers_=1;
     268          26 :     walker_rank_=0;
     269             :   }
     270             : 
     271             : //parallelization stuff
     272          30 :   NumOMP_=OpenMP::getNumThreads();
     273          30 :   NumParallel_=comm.Get_size();
     274          30 :   rank_=comm.Get_rank();
     275          30 :   bool serial=false;
     276          30 :   parseFlag("SERIAL",serial);
     277          30 :   if(serial)
     278             :   {
     279           5 :     NumOMP_=1;
     280           5 :     NumParallel_=1;
     281           5 :     rank_=0;
     282             :   }
     283             : 
     284          30 :   checkRead();
     285             : 
     286             : //check ECVs and link them
     287          30 :   init_pntrToECVsClass();
     288             : //set kbt_
     289          30 :   kbt_=pntrToECVsClass_[0]->getKbT();
     290          67 :   for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     291          37 :     plumed_massert(std::abs(kbt_-pntrToECVsClass_[l]->getKbT())<1e-4,"must set same TEMP for each ECV");
     292             : 
     293             : //restart if needed
     294          30 :   if(getRestart())
     295             :   {
     296             :     bool stateRestart=true;
     297          10 :     if(restartFileName.length()==0)
     298             :     {
     299             :       stateRestart=false;
     300             :       restartFileName=deltaFsFileName;
     301             :     }
     302          10 :     IFile ifile;
     303          10 :     ifile.link(*this);
     304          10 :     if(ifile.FileExist(restartFileName))
     305             :     {
     306          10 :       log.printf("  RESTART - make sure all ECVs used are the same as before\n");
     307          10 :       log.printf("    restarting from: %s\n",restartFileName.c_str());
     308          10 :       ifile.open(restartFileName);
     309          10 :       if(stateRestart) //get all info
     310             :       {
     311           2 :         log.printf("    it should be a STATE file (not a DELTAFS file)\n");
     312             :         double time; //not used
     313           2 :         ifile.scanField("time",time);
     314           2 :         ifile.scanField("counter",counter_);
     315           4 :         ifile.scanField("rct",rct_);
     316             :         std::string tmp_lambda;
     317          66 :         while(ifile.scanField(getPntrToArgument(0)->getName(),tmp_lambda))
     318             :         {
     319          64 :           std::string subs="DeltaF_"+tmp_lambda;
     320         128 :           for(unsigned jj=1; jj<ncv_; jj++)
     321             :           {
     322             :             tmp_lambda.clear();
     323          64 :             ifile.scanField(getPntrToArgument(jj)->getName(),tmp_lambda);
     324         128 :             subs+="_"+tmp_lambda;
     325             :           }
     326          64 :           deltaF_name_.push_back(subs);
     327             :           double tmp_deltaF;
     328          64 :           ifile.scanField("DeltaF",tmp_deltaF);
     329          64 :           deltaF_.push_back(tmp_deltaF);
     330          64 :           ifile.scanField();
     331             :           tmp_lambda.clear();
     332             :         }
     333           2 :         log.printf("  successfully read %lu DeltaF values\n",deltaF_name_.size());
     334           2 :         if(NumParallel_>1)
     335           2 :           all_deltaF_=deltaF_;
     336             :       }
     337             :       else //get just deltaFs names
     338             :       {
     339           8 :         ifile.scanFieldList(deltaF_name_);
     340           8 :         plumed_massert(deltaF_name_.size()>=4,"RESTART - fewer than expected FIELDS found in '"+deltaFsFileName+"' file");
     341           8 :         plumed_massert(deltaF_name_[deltaF_name_.size()-1]=="print_stride","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
     342           8 :         plumed_massert(deltaF_name_[0]=="time","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
     343           8 :         plumed_massert(deltaF_name_[1]=="rct","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
     344             :         deltaF_name_.pop_back();
     345             :         deltaF_name_.erase(deltaF_name_.begin(),deltaF_name_.begin()+2);
     346             :         std::size_t pos=5; //each name starts with "DeltaF"
     347          22 :         for(unsigned j=0; j<ncv_; j++)
     348          14 :           pos=deltaF_name_[0].find("_",pos+1); //checking only first one, hopefully is enough
     349           8 :         plumed_massert(pos<deltaF_name_[0].length(),"RESTART - fewer '_' than expected in DeltaF fields: did you remove any CV?");
     350           8 :         pos=deltaF_name_[0].find("_",pos+1);
     351           8 :         plumed_massert(pos>deltaF_name_[0].length(),"RESTART - more '_' than expected in DeltaF fields: did you add new CV?");
     352             :       }
     353             :       //get lambdas, init ECVs and Link them
     354          10 :       deltaF_size_=deltaF_name_.size();
     355         525 :       auto getLambdaName=[](const std::string& name,const unsigned start,const unsigned dim)
     356             :       {
     357             :         std::size_t pos_start=5; //each name starts with "DeltaF"
     358        1068 :         for(unsigned j=0; j<=start; j++)
     359         543 :           pos_start=name.find("_",pos_start+1);
     360             :         std::size_t pos_end=pos_start;
     361        1527 :         for(unsigned j=0; j<dim; j++)
     362        1002 :           pos_end=name.find("_",pos_end+1);
     363         525 :         pos_start++; //do not include heading "_"
     364         525 :         return name.substr(pos_start,pos_end-pos_start);
     365             :       };
     366          10 :       unsigned index_j=ncv_;
     367             :       unsigned sizeSkip=1;
     368          22 :       for(int l=pntrToECVsClass_.size()-1; l>=0; l--)
     369             :       {
     370          12 :         const unsigned dim_l=pntrToECVsClass_[l]->getNumberOfArguments();
     371          12 :         index_j-=dim_l;
     372          12 :         std::vector<std::string> lambdas_l(1);
     373          12 :         lambdas_l[0]=getLambdaName(deltaF_name_[0],index_j,dim_l);
     374         523 :         for(unsigned i=sizeSkip; i<deltaF_size_; i+=sizeSkip)
     375             :         {
     376         513 :           std::string tmp_lambda=getLambdaName(deltaF_name_[i],index_j,dim_l);
     377         513 :           if(tmp_lambda==lambdas_l[0])
     378             :             break;
     379         511 :           lambdas_l.push_back(tmp_lambda);
     380             :         }
     381          12 :         pntrToECVsClass_[l]->initECVs_restart(lambdas_l);
     382          12 :         sizeSkip*=lambdas_l.size();
     383          12 :       }
     384          10 :       plumed_massert(sizeSkip==deltaF_size_,"RESTART - this should not happen");
     385          10 :       init_linkECVs(); //link ECVs and initializes index_k_
     386          10 :       log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
     387          10 :       obs_steps_=0; //avoid initializing again
     388          10 :       if(stateRestart)
     389             :       {
     390           2 :         if(NumParallel_>1)
     391             :         {
     392           2 :           const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
     393             :           unsigned iter=0;
     394          34 :           for(unsigned i=start; i<start+deltaF_.size(); i++)
     395          32 :             deltaF_[iter++]=all_deltaF_[i];
     396             :         }
     397             :       }
     398             :       else //read each step
     399             :       {
     400           8 :         counter_=1;
     401             :         unsigned count_lines=0;
     402           8 :         ifile.allowIgnoredFields(); //this allows for multiple restart, but without checking for consistency between them!
     403             :         double time;
     404          48 :         while(ifile.scanField("time",time)) //only number of lines and last line is important
     405             :         {
     406             :           unsigned restart_stride;
     407          16 :           ifile.scanField("print_stride",restart_stride);
     408          16 :           ifile.scanField("rct",rct_);
     409          16 :           if(NumParallel_==1)
     410             :           {
     411        1014 :             for(unsigned i=0; i<deltaF_size_; i++)
     412         998 :               ifile.scanField(deltaF_name_[i],deltaF_[i]);
     413             :           }
     414             :           else
     415             :           {
     416           0 :             const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
     417             :             unsigned iter=0;
     418           0 :             for(unsigned i=start; i<start+deltaF_.size(); i++)
     419           0 :               ifile.scanField(deltaF_name_[i],deltaF_[iter++]);
     420             :           }
     421          16 :           ifile.scanField();
     422          16 :           if(count_lines>0)
     423           8 :             counter_+=restart_stride;
     424          16 :           count_lines++;
     425             :         }
     426           8 :         counter_*=NumWalkers_;
     427           8 :         log.printf("  successfully read %u lines, up to t=%g\n",count_lines,time);
     428             :       }
     429          10 :       ifile.reset(false);
     430          10 :       ifile.close();
     431             :     }
     432             :     else //same behaviour as METAD
     433           0 :       plumed_merror("RESTART requested, but file '"+restartFileName+"' was not found!\n  Set RESTART=NO or provide a restart file");
     434          10 :     if(NumWalkers_>1) //make sure that all walkers are doing the same thing
     435             :     {
     436           2 :       std::vector<unsigned long long> all_counter(NumWalkers_);
     437           2 :       if(comm.Get_rank()==0)
     438           2 :         multi_sim_comm.Allgather(counter_,all_counter);
     439           2 :       comm.Bcast(all_counter,0);
     440             :       bool same_number_of_steps=true;
     441           4 :       for(unsigned w=1; w<NumWalkers_; w++)
     442           2 :         if(all_counter[0]!=all_counter[w])
     443             :           same_number_of_steps=false;
     444           2 :       plumed_massert(same_number_of_steps,"RESTART - not all walkers are reading the same file!");
     445             :     }
     446          10 :   }
     447          20 :   else if(restartFileName.length()>0)
     448           0 :     log.printf(" +++ WARNING +++ the provided STATE_RFILE will be ignored, since RESTART was not requested\n");
     449             : 
     450             : //sync all walkers to avoid opening files before reding is over (see also METAD)
     451          30 :   comm.Barrier();
     452          30 :   if(comm.Get_rank()==0 && walkers_mpi)
     453           4 :     multi_sim_comm.Barrier();
     454             : 
     455             : //setup DeltaFs file
     456          30 :   deltaFsOfile_.link(*this);
     457          30 :   if(NumWalkers_>1)
     458             :   {
     459           4 :     if(walker_rank_>0)
     460             :       deltaFsFileName="/dev/null"; //only first walker writes on file
     461           8 :     deltaFsOfile_.enforceSuffix("");
     462             :   }
     463          30 :   deltaFsOfile_.open(deltaFsFileName);
     464          30 :   if(fmt.length()>0)
     465          60 :     deltaFsOfile_.fmtField(" "+fmt);
     466             :   deltaFsOfile_.setHeavyFlush(); //do I need it?
     467          30 :   deltaFsOfile_.addConstantField("print_stride");
     468          30 :   deltaFsOfile_.printField("print_stride",print_stride_);
     469             : 
     470             : //open file for storing state
     471          30 :   if(wStateStride_!=0)
     472             :   {
     473           6 :     stateOfile_.link(*this);
     474           6 :     if(NumWalkers_>1)
     475             :     {
     476           0 :       if(walker_rank_>0)
     477             :         stateFileName="/dev/null"; //only first walker writes on file
     478           0 :       stateOfile_.enforceSuffix("");
     479             :     }
     480           6 :     stateOfile_.open(stateFileName);
     481           6 :     if(fmt.length()>0)
     482          12 :       stateOfile_.fmtField(" "+fmt);
     483             :   }
     484             : 
     485             : //add output components
     486          30 :   if(calc_work_)
     487             :   {
     488          12 :     addComponent("work");
     489          12 :     componentIsNotPeriodic("work");
     490             :   }
     491             : 
     492             : //printing some info
     493          30 :   log.printf("  updating the bias with PACE = %u\n",stride_);
     494          30 :   log.printf("  initial unbiased OBSERVATION_STEPS = %u (in units of PACE)\n",obs_steps_);
     495          30 :   if(wStateStride_>0)
     496           5 :     log.printf("  state checkpoints are written on file %s every %d MD steps\n",stateFileName.c_str(),wStateStride_);
     497          30 :   if(wStateStride_==-1)
     498           1 :     log.printf("  state checkpoints are written on file '%s' only on CPT events (or never if MD code does define them!)\n",stateFileName.c_str());
     499          30 :   if(walkers_mpi)
     500           4 :     log.printf(" -- WALKERS_MPI: if multiple replicas are present, they will share the same bias via MPI\n");
     501          30 :   if(NumWalkers_>1)
     502             :   {
     503           4 :     log.printf("  using multiple walkers\n");
     504           4 :     log.printf("    number of walkers: %u\n",NumWalkers_);
     505           4 :     log.printf("    walker rank: %u\n",walker_rank_);
     506             :   }
     507          30 :   int mw_warning=0;
     508          30 :   if(!walkers_mpi && comm.Get_rank()==0 && multi_sim_comm.Get_size()>(int)NumWalkers_)
     509           0 :     mw_warning=1;
     510          30 :   comm.Bcast(mw_warning,0);
     511          30 :   if(mw_warning) //log.printf messes up with comm, so never use it without Bcast!
     512           0 :     log.printf(" +++ WARNING +++ multiple replicas will NOT communicate unless the flag WALKERS_MPI is used\n");
     513          30 :   if(NumParallel_>1)
     514           2 :     log.printf("  using multiple MPI processes per simulation: %u\n",NumParallel_);
     515          30 :   if(NumOMP_>1)
     516          25 :     log.printf("  using multiple OpenMP threads per simulation: %u\n",NumOMP_);
     517          30 :   if(serial)
     518           5 :     log.printf(" -- SERIAL: no loop parallelization, despite %d MPI processes and %u OpenMP threads available\n",comm.Get_size(),OpenMP::getNumThreads());
     519          30 :   log.printf("  Bibliography: ");
     520          60 :   log<<plumed.cite("M. Invernizzi, P.M. Piaggi, and M. Parrinello, Phys. Rev. X 10, 041034 (2020)");
     521          30 :   log.printf("\n");
     522          30 : }
     523             : 
     524        1490 : void OPESexpanded::calculate()
     525             : {
     526        1490 :   if(deltaF_size_==0) //no bias before initialization
     527         325 :     return;
     528             : 
     529             : //get diffMax, to avoid over/underflow
     530        1165 :   double diffMax=-std::numeric_limits<double>::max();
     531        1165 :   #pragma omp parallel num_threads(NumOMP_)
     532             :   {
     533             :     #pragma omp for reduction(max:diffMax)
     534             :     for(unsigned i=0; i<deltaF_.size(); i++)
     535             :     {
     536             :       diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
     537             :       if(diff_[i]>diffMax)
     538             :         diffMax=diff_[i];
     539             :     }
     540             :   }
     541        1165 :   if(NumParallel_>1)
     542         102 :     comm.Max(diffMax);
     543             : 
     544             : //calculate the bias and the forces
     545        1165 :   double sum=0;
     546        1165 :   std::vector<double> der_sum_cv(ncv_,0);
     547        1165 :   if(NumOMP_==1)
     548             :   {
     549        2730 :     for(unsigned i=0; i<deltaF_.size(); i++)
     550             :     {
     551        2520 :       double add_i=std::exp(diff_[i]-diffMax);
     552        2520 :       sum+=add_i;
     553             :       //set derivatives
     554        6960 :       for(unsigned j=0; j<ncv_; j++)
     555        4440 :         der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
     556             :     }
     557             :   }
     558             :   else
     559             :   {
     560         955 :     #pragma omp parallel num_threads(NumOMP_)
     561             :     {
     562             :       std::vector<double> omp_der_sum_cv(ncv_,0);
     563             :       #pragma omp for reduction(+:sum) nowait
     564             :       for(unsigned i=0; i<deltaF_.size(); i++)
     565             :       {
     566             :         double add_i=std::exp(diff_[i]-diffMax);
     567             :         sum+=add_i;
     568             :         //set derivatives
     569             :         for(unsigned j=0; j<ncv_; j++)
     570             :           omp_der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
     571             :       }
     572             :       #pragma omp critical
     573             :       for(unsigned j=0; j<ncv_; j++)
     574             :         der_sum_cv[j]+=omp_der_sum_cv[j];
     575             :     }
     576             :   }
     577        1165 :   if(NumParallel_>1)
     578             :   { //each MPI process has part of the full deltaF_ vector, so must Sum
     579         102 :     comm.Sum(sum);
     580         102 :     comm.Sum(der_sum_cv);
     581             :   }
     582             : 
     583             : //set bias and forces
     584        1165 :   const double bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
     585        1165 :   setBias(bias);
     586        3163 :   for(unsigned j=0; j<ncv_; j++)
     587        1998 :     setOutputForce(j,kbt_*der_sum_cv[j]/sum);
     588             : }
     589             : 
     590        1490 : void OPESexpanded::update()
     591             : {
     592        1490 :   if(isFirstStep_) //skip very first step, as in METAD
     593             :   {
     594          30 :     isFirstStep_=false;
     595          30 :     if(obs_steps_!=1) //if obs_steps_==1 go on with initialization
     596             :       return;
     597             :   }
     598        1464 :   if(getStep()%stride_==0)
     599             :   {
     600         739 :     if(obs_steps_>0)
     601             :     {
     602         463 :       for(unsigned j=0; j<ncv_; j++)
     603         304 :         obs_cvs_[counter_*ncv_+j]=getArgument(j);
     604         159 :       counter_++;
     605         159 :       if(counter_==obs_steps_)
     606             :       {
     607          20 :         log.printf("\nAction %s\n",getName().c_str());
     608          20 :         init_fromObs();
     609          20 :         log.printf("Finished initialization\n\n");
     610          20 :         counter_=NumWalkers_; //all preliminary observations count 1
     611          20 :         obs_steps_=0; //no more observation
     612             :       }
     613         159 :       return;
     614             :     }
     615             : 
     616             :     //update averages
     617         580 :     const double current_bias=getOutputQuantity(0); //the first value is always the bias
     618         580 :     if(NumWalkers_==1)
     619         500 :       updateDeltaF(current_bias);
     620             :     else
     621             :     {
     622          80 :       std::vector<double> cvs(ncv_);
     623         240 :       for(unsigned j=0; j<ncv_; j++)
     624         160 :         cvs[j]=getArgument(j);
     625          80 :       std::vector<double> all_bias(NumWalkers_);
     626          80 :       std::vector<double> all_cvs(NumWalkers_*ncv_);
     627          80 :       if(comm.Get_rank()==0)
     628             :       {
     629          80 :         multi_sim_comm.Allgather(current_bias,all_bias);
     630          80 :         multi_sim_comm.Allgather(cvs,all_cvs);
     631             :       }
     632          80 :       comm.Bcast(all_bias,0);
     633          80 :       comm.Bcast(all_cvs,0);
     634         240 :       for(unsigned w=0; w<NumWalkers_; w++)
     635             :       {
     636             :         //calculate ECVs
     637         160 :         unsigned index_wj=w*ncv_;
     638         380 :         for(unsigned k=0; k<pntrToECVsClass_.size(); k++)
     639             :         {
     640         220 :           pntrToECVsClass_[k]->calculateECVs(&all_cvs[index_wj]);
     641         220 :           index_wj+=pntrToECVsClass_[k]->getNumberOfArguments();
     642             :         }
     643         160 :         updateDeltaF(all_bias[w]);
     644             :       }
     645             :     }
     646             : 
     647             :     //write DeltaFs to file
     648         580 :     if((counter_/NumWalkers_-1)%print_stride_==0)
     649          44 :       printDeltaF();
     650             : 
     651             :     //calculate work if requested
     652         580 :     if(calc_work_)
     653             :     { //some copy and paste from calculate()
     654             :       //get diffMax, to avoid over/underflow
     655         110 :       double diffMax=-std::numeric_limits<double>::max();
     656         110 :       #pragma omp parallel num_threads(NumOMP_)
     657             :       {
     658             :         #pragma omp for reduction(max:diffMax)
     659             :         for(unsigned i=0; i<deltaF_.size(); i++)
     660             :         {
     661             :           diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
     662             :           if(diff_[i]>diffMax)
     663             :             diffMax=diff_[i];
     664             :         }
     665             :       }
     666         110 :       if(NumParallel_>1)
     667          50 :         comm.Max(diffMax);
     668             :       //calculate the bias
     669         110 :       double sum=0;
     670         110 :       #pragma omp parallel num_threads(NumOMP_)
     671             :       {
     672             :         #pragma omp for reduction(+:sum) nowait
     673             :         for(unsigned i=0; i<deltaF_.size(); i++)
     674             :           sum+=std::exp(diff_[i]-diffMax);
     675             :       }
     676         110 :       if(NumParallel_>1)
     677          50 :         comm.Sum(sum);
     678         110 :       const double new_bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
     679             :       //accumulate work
     680         110 :       work_+=new_bias-current_bias;
     681         220 :       getPntrToComponent("work")->set(work_);
     682             :     }
     683             :   }
     684             : 
     685             : //dump state if requested
     686        1305 :   if( (wStateStride_>0 && getStep()%wStateStride_==0) || (wStateStride_==-1 && getCPT()) )
     687          11 :     dumpStateToFile();
     688             : }
     689             : 
     690          30 : void OPESexpanded::init_pntrToECVsClass()
     691             : {
     692          30 :   std::vector<opes::ExpansionCVs*> all_pntrToECVsClass=plumed.getActionSet().select<opes::ExpansionCVs*>();
     693          30 :   plumed_massert(all_pntrToECVsClass.size()>0,"no Expansion CVs found");
     694          67 :   for(unsigned j=0; j<ncv_; j++)
     695             :   {
     696          74 :     std::string error_notECV("all the ARGs of "+getName()+" must be Expansion Collective Variables (ECV)");
     697          37 :     const unsigned dot_pos=getPntrToArgument(j)->getName().find(".");
     698          37 :     plumed_massert(dot_pos<getPntrToArgument(j)->getName().size(),error_notECV+", thus contain a dot in the name");
     699          37 :     unsigned foundECV_l=all_pntrToECVsClass.size();
     700          44 :     for(unsigned l=0; l<all_pntrToECVsClass.size(); l++)
     701             :     {
     702          44 :       if(getPntrToArgument(j)->getName().substr(0,dot_pos)==all_pntrToECVsClass[l]->getLabel())
     703             :       {
     704             :         foundECV_l=l;
     705          37 :         pntrToECVsClass_.push_back(all_pntrToECVsClass[l]);
     706          37 :         std::string missing_arg="some ECV component is missing from ARG";
     707          37 :         plumed_massert(j+all_pntrToECVsClass[l]->getNumberOfArguments()<=getNumberOfArguments(),missing_arg);
     708          90 :         for(unsigned h=0; h<all_pntrToECVsClass[l]->getNumberOfArguments(); h++)
     709             :         {
     710          53 :           std::string argName=getPntrToArgument(j+h)->getName();
     711          53 :           std::string expectedECVname=all_pntrToECVsClass[l]->getComponentsVector()[h];
     712          53 :           plumed_massert(argName==expectedECVname,missing_arg+", or is in wrong order: given ARG="+argName+" expected ARG="+expectedECVname);
     713             :         }
     714          37 :         j+=all_pntrToECVsClass[l]->getNumberOfArguments()-1;
     715             :         break;
     716             :       }
     717             :     }
     718          37 :     plumed_massert(foundECV_l<all_pntrToECVsClass.size(),error_notECV);
     719             :   }
     720          67 :   for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     721          44 :     for(unsigned ll=l+1; ll<pntrToECVsClass_.size(); ll++)
     722           7 :       plumed_massert(pntrToECVsClass_[l]->getLabel()!=pntrToECVsClass_[ll]->getLabel(),"cannot use same ECV twice");
     723          30 : }
     724             : 
     725          30 : void OPESexpanded::init_linkECVs()
     726             : {
     727             :   //TODO It should be possible to make all of this more straightforward (and probably also faster):
     728             :   //     - get rid of index_k_, making it trivial for each ECV
     729             :   //     - store the ECVs_ and derECVs_ vectors here as a contiguous vector, and use pointers in the ECV classes
     730             :   //     Some caveats:
     731             :   //     - ECVmultiThermalBaric has a nontrivial index_k_ to avoid duplicates. use duplicates instead
     732             :   //     - can the ECVs be MPI parallel or it's too complicated?
     733          30 :   plumed_massert(deltaF_size_>0,"must set deltaF_size_ before calling init_linkECVs()");
     734          30 :   if(NumParallel_==1)
     735          28 :     deltaF_.resize(deltaF_size_);
     736             :   else
     737             :   {
     738           2 :     const unsigned extra=(rank_<(deltaF_size_%NumParallel_)?1:0);
     739           2 :     deltaF_.resize(deltaF_size_/NumParallel_+extra);
     740             :     //these are used when printing deltaF_ to file
     741           2 :     all_deltaF_.resize(deltaF_size_);
     742           2 :     all_size_.resize(NumParallel_,deltaF_size_/NumParallel_);
     743           2 :     disp_.resize(NumParallel_);
     744           4 :     for(unsigned r=0; r<NumParallel_-1; r++)
     745             :     {
     746           2 :       if(r<deltaF_size_%NumParallel_)
     747           0 :         all_size_[r]++;
     748           2 :       disp_[r+1]=disp_[r]+all_size_[r];
     749             :     }
     750             :   }
     751          30 :   diff_.resize(deltaF_.size());
     752          30 :   ECVs_.resize(ncv_);
     753          30 :   derECVs_.resize(ncv_);
     754          30 :   index_k_.resize(deltaF_.size(),std::vector<unsigned>(ncv_));
     755             :   unsigned index_j=0;
     756          30 :   unsigned sizeSkip=deltaF_size_;
     757          67 :   for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     758             :   {
     759          37 :     std::vector< std::vector<unsigned> > l_index_k(pntrToECVsClass_[l]->getIndex_k());
     760          37 :     plumed_massert(deltaF_size_%l_index_k.size()==0,"buggy ECV: mismatch between getTotNumECVs() and getIndex_k().size()");
     761          37 :     plumed_massert(l_index_k[0].size()==pntrToECVsClass_[l]->getNumberOfArguments(),"buggy ECV: mismatch between number of ARG and underlying CVs");
     762          37 :     sizeSkip/=l_index_k.size();
     763          90 :     for(unsigned h=0; h<pntrToECVsClass_[l]->getNumberOfArguments(); h++)
     764             :     {
     765          53 :       ECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToECVs(h);
     766          53 :       derECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToDerECVs(h);
     767          53 :       if(NumParallel_==1)
     768             :       {
     769       45589 :         for(unsigned i=0; i<deltaF_size_; i++)
     770       45540 :           index_k_[i][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
     771             :       }
     772             :       else
     773             :       {
     774           4 :         const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
     775             :         unsigned iter=0;
     776          68 :         for(unsigned i=start; i<start+deltaF_.size(); i++)
     777          64 :           index_k_[iter++][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
     778             :       }
     779             :     }
     780          37 :     index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
     781          37 :   }
     782          30 :   plumed_massert(sizeSkip==1,"this should not happen!");
     783          30 : }
     784             : 
     785          20 : void OPESexpanded::init_fromObs() //This could probably be faster and/or require less memory...
     786             : {
     787             : //in case of multiple walkers gather all the statistics
     788          20 :   if(NumWalkers_>1)
     789             :   {
     790           2 :     std::vector<double> all_obs_cv(ncv_*obs_steps_*NumWalkers_);
     791           2 :     if(comm.Get_rank()==0)
     792           2 :       multi_sim_comm.Allgather(obs_cvs_,all_obs_cv);
     793           2 :     comm.Bcast(all_obs_cv,0);
     794           2 :     obs_cvs_=all_obs_cv; //could this lead to memory issues?
     795           2 :     obs_steps_*=NumWalkers_;
     796             :   }
     797             : 
     798             : //initialize ECVs from observations
     799             :   unsigned index_j=0;
     800          20 :   deltaF_size_=1;
     801          45 :   for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     802             :   {
     803          25 :     pntrToECVsClass_[l]->initECVs_observ(obs_cvs_,ncv_,index_j);
     804          25 :     deltaF_size_*=pntrToECVsClass_[l]->getTotNumECVs(); //ECVs from different exansions will be combined
     805          25 :     index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
     806             :   }
     807          20 :   plumed_massert(index_j==getNumberOfArguments(),"mismatch between number of linked CVs and number of ARG");
     808             : //link ECVs and initialize index_k_, mapping each deltaF to a single ECVs set
     809          20 :   init_linkECVs();
     810             : 
     811             : //initialize deltaF_ from obs
     812             : //for the first point, t=0, the ECVs are calculated by initECVs_observ, setting also any initial guess
     813             :   index_j=0;
     814       12379 :   for(unsigned i=0; i<deltaF_.size(); i++)
     815       56923 :     for(unsigned j=0; j<ncv_; j++)
     816       44564 :       deltaF_[i]+=kbt_*ECVs_[j][index_k_[i][j]];
     817         179 :   for(unsigned t=1; t<obs_steps_; t++) //starts from t=1
     818             :   {
     819             :     unsigned index_j=0;
     820         383 :     for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     821             :     {
     822         224 :       pntrToECVsClass_[l]->calculateECVs(&obs_cvs_[t*ncv_+index_j]);
     823         224 :       index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
     824             :     }
     825      102677 :     for(unsigned i=0; i<deltaF_.size(); i++)
     826             :     {
     827      102518 :       const double diff_i=(-getExpansion(i)+deltaF_[i]/kbt_-std::log(t));
     828      102518 :       if(diff_i>0) //save exp from overflow
     829       16071 :         deltaF_[i]-=kbt_*(diff_i+std::log1p(std::exp(-diff_i))+std::log1p(-1./(1.+t)));
     830             :       else
     831       86447 :         deltaF_[i]-=kbt_*(std::log1p(std::exp(diff_i))+std::log1p(-1./(1.+t)));
     832             :     }
     833             :   }
     834             :   obs_cvs_.clear();
     835             : 
     836             : //set deltaF_name_
     837          20 :   deltaF_name_.resize(deltaF_size_,"DeltaF");
     838          20 :   unsigned sizeSkip=deltaF_size_;
     839          45 :   for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
     840             :   {
     841          25 :     std::vector<std::string> lambdas_l=pntrToECVsClass_[l]->getLambdas();
     842          25 :     plumed_massert(lambdas_l.size()==pntrToECVsClass_[l]->getTotNumECVs(),"buggy ECV: mismatch between getTotNumECVs() and getLambdas().size()");
     843          25 :     sizeSkip/=lambdas_l.size();
     844       22457 :     for(unsigned i=0; i<deltaF_size_; i++)
     845       44864 :       deltaF_name_[i]+="_"+lambdas_l[(i/sizeSkip)%lambdas_l.size()];
     846          25 :   }
     847             : 
     848             : //print initialization to file
     849          20 :   log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
     850          20 :   printDeltaF();
     851          20 : }
     852             : 
     853          64 : void OPESexpanded::printDeltaF()
     854             : {
     855          64 :   deltaFsOfile_.printField("time",getTime());
     856          64 :   deltaFsOfile_.printField("rct",rct_);
     857          64 :   if(NumParallel_==1)
     858             :   {
     859       23988 :     for(unsigned i=0; i<deltaF_.size(); i++)
     860       23926 :       deltaFsOfile_.printField(deltaF_name_[i],deltaF_[i]);
     861             :   }
     862             :   else
     863             :   {
     864           2 :     comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
     865          66 :     for(unsigned i=0; i<deltaF_size_; i++)
     866          64 :       deltaFsOfile_.printField(deltaF_name_[i],all_deltaF_[i]);
     867             :   }
     868          64 :   deltaFsOfile_.printField();
     869          64 : }
     870             : 
     871          11 : void OPESexpanded::dumpStateToFile()
     872             : {
     873             : //rewrite header or rewind file
     874          11 :   if(storeOldStates_)
     875           3 :     stateOfile_.clearFields();
     876           8 :   else if(walker_rank_==0)
     877           8 :     stateOfile_.rewind();
     878             : //define fields
     879          11 :   stateOfile_.addConstantField("time");
     880          11 :   stateOfile_.addConstantField("counter");
     881          11 :   stateOfile_.addConstantField("rct");
     882             : //print
     883          11 :   stateOfile_.printField("time",getTime());
     884          11 :   stateOfile_.printField("counter",counter_);
     885          11 :   stateOfile_.printField("rct",rct_);
     886          11 :   if(NumParallel_>1)
     887           0 :     comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
     888         240 :   for(unsigned i=0; i<deltaF_size_; i++)
     889             :   {
     890             :     std::size_t pos_start=7; //skip "DeltaF_"
     891         687 :     for(unsigned j=0; j<ncv_; j++)
     892             :     {
     893             :       plumed_dbg_massert(pos_start>6,"not enought _ in deltaF_name_"+std::to_string(i-1)+" string?");
     894         458 :       const std::size_t pos_end=deltaF_name_[i].find("_",pos_start);
     895         916 :       stateOfile_.printField(getPntrToArgument(j)->getName(),"  "+deltaF_name_[i].substr(pos_start,pos_end-pos_start));
     896         458 :       pos_start=pos_end+1;
     897             :     }
     898         229 :     if(NumParallel_==1)
     899         458 :       stateOfile_.printField("DeltaF",deltaF_[i]);
     900             :     else
     901           0 :       stateOfile_.printField("DeltaF",all_deltaF_[i]);
     902         229 :     stateOfile_.printField();
     903             :   }
     904             : //make sure file is written even if small
     905          11 :   if(!storeOldStates_)
     906           8 :     stateOfile_.flush();
     907          11 : }
     908             : 
     909         660 : void OPESexpanded::updateDeltaF(double bias)
     910             : {
     911             :   plumed_dbg_massert(counter_>0,"deltaF_ must be initialized");
     912         660 :   counter_++;
     913         660 :   const double arg=(bias-rct_)/kbt_-std::log(counter_-1.);
     914             :   double increment;
     915         660 :   if(arg>0) //save exp from overflow
     916          28 :     increment=kbt_*(arg+std::log1p(std::exp(-arg)));
     917             :   else
     918         632 :     increment=kbt_*(std::log1p(std::exp(arg)));
     919         660 :   #pragma omp parallel num_threads(NumOMP_)
     920             :   {
     921             :     #pragma omp for
     922             :     for(unsigned i=0; i<deltaF_.size(); i++)
     923             :     {
     924             :       const double diff_i=(-getExpansion(i)+(bias-rct_+deltaF_[i])/kbt_-std::log(counter_-1.));
     925             :       if(diff_i>0) //save exp from overflow
     926             :         deltaF_[i]+=increment-kbt_*(diff_i+std::log1p(std::exp(-diff_i)));
     927             :       else
     928             :         deltaF_[i]+=increment-kbt_*std::log1p(std::exp(diff_i));
     929             :     }
     930             :   }
     931         660 :   rct_+=increment+kbt_*std::log1p(-1./counter_);
     932         660 : }
     933             : 
     934      644469 : double OPESexpanded::getExpansion(unsigned i) const
     935             : {
     936             :   double expansion=0;
     937     3003062 :   for(unsigned j=0; j<ncv_; j++)
     938     2358593 :     expansion+=ECVs_[j][index_k_[i][j]]; //the index_k could be trivially guessed for most ECVs, but unfourtunately not all
     939      644469 :   return expansion;
     940             : }
     941             : 
     942             : }
     943             : }

Generated by: LCOV version 1.16