LCOV - code coverage report
Current view: top level - ves - VesDeltaF.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 343 351 97.7 %
Date: 2025-03-25 09:33:27 Functions: 7 8 87.5 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2016-2021 The VES code team
       3             :    (see the PEOPLE-VES file at the root of this folder for a list of names)
       4             : 
       5             :    See http://www.ves-code.org for more information.
       6             : 
       7             :    This file is part of VES code module.
       8             : 
       9             :    The VES code module 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             :    The VES code module 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 the VES code module.  If not, see <http://www.gnu.org/licenses/>.
      21             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      22             : 
      23             : #include "bias/Bias.h"
      24             : #include "core/PlumedMain.h"
      25             : #include "core/ActionRegister.h"
      26             : #include "tools/Communicator.h"
      27             : #include "tools/Grid.h"
      28             : #include "tools/File.h"
      29             : //#include <algorithm> //std::fill
      30             : 
      31             : namespace PLMD {
      32             : namespace ves {
      33             : 
      34             : //+PLUMEDOC VES_BIAS VES_DELTA_F
      35             : /*
      36             : Implementation of VES Delta F method
      37             : 
      38             : Implementation of VES\f$\Delta F\f$ method \cite Invernizzi2019vesdeltaf (step two only).
      39             : 
      40             : \warning
      41             :   Notice that this is a stand-alone bias Action, it does not need any of the other VES module components
      42             : 
      43             : First you should create some estimate of the local free energy basins of your system,
      44             : using e.g. multiple \ref METAD short runs, and combining them with the \ref sum_hills utility.
      45             : Once you have them, you can use this bias Action to perform the VES optimization part of the method.
      46             : 
      47             : These \f$N+1\f$ local basins are used to model the global free energy.
      48             : In particular, given the conditional probabilities \f$P(\mathbf{s}|i)\propto e^{-\beta F_i(\mathbf{s})}\f$
      49             : and the probabilities of being in a given basin \f$P_i\f$, we can write:
      50             : \f[
      51             :   e^{-\beta F(\mathbf{s})}\propto P(\mathbf{s})=\sum_{i=0}^N P(\mathbf{s}|i)P_i \, .
      52             : \f]
      53             : We use this free energy model and the chosen bias factor \f$\gamma\f$ to build the bias potential:
      54             : \f$V(\mathbf{s})=-(1-1/\gamma)F(\mathbf{s})\f$.
      55             : Or, more explicitly:
      56             : \f[
      57             :   V(\mathbf{s})=(1-1/\gamma)\frac{1}{\beta}\log\left[e^{-\beta F_0(\mathbf{s})}
      58             :   +\sum_{i=1}^{N} e^{-\beta F_i(\mathbf{s})} e^{-\beta \alpha_i}\right] \, ,
      59             : \f]
      60             : where the parameters \f$\boldsymbol{\alpha}\f$ are the \f$N\f$ free energy differences (see below) from the \f$F_0\f$ basin.
      61             : 
      62             : By default the \f$F_i(\mathbf{s})\f$ are shifted so that \f$\min[F_i(\mathbf{s})]=0\f$ for all \f$i=\{0,...,N\}\f$.
      63             : In this case the optimization parameters \f$\alpha_i\f$ are the difference in height between the minima of the basins.
      64             : Using the keyword `NORMALIZE`, you can also decide to normalize the local free energies so that
      65             : \f$\int d\mathbf{s}\, e^{-\beta F_i(\mathbf{s})}=1\f$.
      66             : In this case the parameters will represent not the difference in height (which depends on the chosen CVs),
      67             : but the actual free energy difference, \f$\alpha_i=\Delta F_i\f$.
      68             : 
      69             : However, as discussed in Ref. \cite Invernizzi2019vesdeltaf, a better estimate of \f$\Delta F_i\f$ should be obtained through the reweighting procedure.
      70             : 
      71             : \par Examples
      72             : 
      73             : The following performs the optimization of the free energy difference between two metastable basins:
      74             : 
      75             : \plumedfile
      76             : cv: DISTANCE ATOMS=1,2
      77             : ves: VES_DELTA_F ...
      78             :   ARG=cv
      79             :   TEMP=300
      80             :   FILE_F0=fesA.data
      81             :   FILE_F1=fesB.data
      82             :   BIASFACTOR=10.0
      83             :   M_STEP=0.1
      84             :   AV_STRIDE=500
      85             :   PRINT_STRIDE=100
      86             : ...
      87             : PRINT FMT=%g STRIDE=500 FILE=Colvar.data ARG=cv,ves.bias,ves.rct
      88             : \endplumedfile
      89             : 
      90             : The local FES files can be obtained as described in Sec. 4.2 of Ref. \cite Invernizzi2019vesdeltaf, i.e. for example:
      91             : - run 4 independent metad runs, all starting from basin A, and kill them as soon as they make the first transition (see e.g. \ref COMMITTOR)
      92             : - \verbatim cat HILLS* > all_HILLS \endverbatim
      93             : - \verbatim plumed sum_hills --hills all_HILLS --outfile all_fesA.dat --mintozero --min 0 --max 1 --bin 100 \endverbatim
      94             : - \verbatim awk -v n_rep=4 '{if($1!="#!" && $1!="") {for(i=1+(NF-1)/2; i<=NF; i++) $i/=n_rep;} print $0}' all_fesA.dat > fesA.data \endverbatim
      95             : 
      96             : The header of both FES files must be identical, and should be similar to the following:
      97             : 
      98             : \auxfile{fesA.data}
      99             : #! FIELDS cv file.free der_cv
     100             : #! SET min_cv 0
     101             : #! SET max_cv 1
     102             : #! SET nbins_cv  100
     103             : #! SET periodic_cv false
     104             : 0 0 0
     105             : \endauxfile
     106             : \auxfile{fesB.data}
     107             : #! FIELDS cv file.free der_cv
     108             : #! SET min_cv 0
     109             : #! SET max_cv 1
     110             : #! SET nbins_cv  100
     111             : #! SET periodic_cv false
     112             : 0 0 0
     113             : \endauxfile
     114             : 
     115             : */
     116             : //+ENDPLUMEDOC
     117             : 
     118             : class VesDeltaF : public bias::Bias {
     119             : 
     120             : private:
     121             :   double beta_;
     122             :   unsigned NumParallel_;
     123             :   unsigned rank_;
     124             :   unsigned NumWalkers_;
     125             :   bool isFirstStep_;
     126             :   bool afterCalculate_;
     127             : 
     128             : //prob
     129             :   double tot_prob_;
     130             :   std::vector<double> prob_;
     131             :   std::vector< std::vector<double> > der_prob_;
     132             : 
     133             : //local basins
     134             :   std::vector< std::unique_ptr<Grid> > grid_p_; //pointers because of GridBase::create
     135             :   std::vector<double> norm_;
     136             : 
     137             : //optimizer-related stuff
     138             :   long long unsigned mean_counter_;
     139             :   unsigned mean_weight_tau_;
     140             :   unsigned alpha_size_;
     141             :   unsigned sym_alpha_size_;
     142             :   std::vector<double> mean_alpha_;
     143             :   std::vector<double> inst_alpha_;
     144             :   std::vector<double> past_increment2_;
     145             :   double minimization_step_;
     146             :   bool damping_off_;
     147             : //'tg' -> 'target distribution'
     148             :   double inv_gamma_;
     149             :   unsigned tg_counter_;
     150             :   unsigned tg_stride_;
     151             :   std::vector<double> tg_dV_dAlpha_;
     152             :   std::vector<double> tg_d2V_dAlpha2_;
     153             : //'av' -> 'ensemble average'
     154             :   unsigned av_counter_;
     155             :   unsigned av_stride_;
     156             :   std::vector<double> av_dV_dAlpha_;
     157             :   std::vector<double> av_dV_dAlpha_prod_;
     158             :   std::vector<double> av_d2V_dAlpha2_;
     159             : //printing
     160             :   unsigned print_stride_;
     161             :   OFile alphaOfile_;
     162             : //other
     163             :   std::vector<double> exp_alpha_;
     164             :   std::vector<double> prev_exp_alpha_;
     165             :   double work_;
     166             : 
     167             : //functions
     168             :   void update_alpha();
     169             :   void update_tg_and_rct();
     170             :   inline unsigned get_index(const unsigned, const unsigned) const;
     171             : 
     172             : public:
     173             :   explicit VesDeltaF(const ActionOptions&);
     174             :   void calculate() override;
     175             :   void update() override;
     176             :   static void registerKeywords(Keywords& keys);
     177             : };
     178             : 
     179             : PLUMED_REGISTER_ACTION(VesDeltaF,"VES_DELTA_F")
     180             : 
     181           6 : void VesDeltaF::registerKeywords(Keywords& keys) {
     182           6 :   Bias::registerKeywords(keys);
     183           6 :   keys.add("optional","TEMP","temperature is compulsory, but it can be sometimes fetched from the MD engine");
     184             : //local free energies
     185           6 :   keys.add("numbered","FILE_F","names of files containing local free energies and derivatives. "
     186             :            "The first one, FILE_F0, is used as reference for all the free energy differences.");
     187          12 :   keys.reset_style("FILE_F","compulsory");
     188           6 :   keys.addFlag("NORMALIZE",false,"normalize all local free energies so that alpha will be (approx) Delta F");
     189           6 :   keys.addFlag("NO_MINTOZERO",false,"leave local free energies as provided, without shifting them to zero min");
     190             : //target distribution
     191           6 :   keys.add("compulsory","BIASFACTOR","0","the gamma bias factor used for well-tempered target p(s)."
     192             :            " Set to 0 for non-tempered flat target");
     193           6 :   keys.add("optional","TG_STRIDE","( default=1 ) number of AV_STRIDE between updates"
     194             :            " of target p(s) and reweighing factor c(t)");
     195             : //optimization
     196           6 :   keys.add("compulsory","M_STEP","1.0","the mu step used for the Omega functional minimization");
     197           6 :   keys.add("compulsory","AV_STRIDE","500","number of simulation steps between alpha updates");
     198           6 :   keys.add("optional","TAU_MEAN","exponentially decaying average for alpha (rescaled using AV_STRIDE)."
     199             :            " Should be used only in very specific cases");
     200           6 :   keys.add("optional","INITIAL_ALPHA","( default=0 ) an initial guess for the bias potential parameter alpha");
     201           6 :   keys.addFlag("DAMPING_OFF",false,"do not use an AdaGrad-like term to rescale M_STEP");
     202             : //output parameters file
     203           6 :   keys.add("compulsory","ALPHA_FILE","ALPHA","file name for output minimization parameters");
     204           6 :   keys.add("optional","PRINT_STRIDE","( default=10 ) stride for printing to ALPHA_FILE");
     205           6 :   keys.add("optional","FMT","specify format for ALPHA_FILE");
     206             : //debug flags
     207           6 :   keys.addFlag("SERIAL",false,"perform the calculation in serial even if multiple tasks are available");
     208           6 :   keys.addFlag("MULTIPLE_WALKERS",false,"use multiple walkers connected via MPI for the optimization");
     209           6 :   keys.use("RESTART");
     210             : 
     211             : //output components
     212          12 :   keys.addOutputComponent("rct","default","scalar","the reweighting factor c(t)");
     213          12 :   keys.addOutputComponent("work","default","scalar","the work done by the bias in one AV_STRIDE");
     214           6 : }
     215             : 
     216           4 : VesDeltaF::VesDeltaF(const ActionOptions&ao)
     217             :   : PLUMED_BIAS_INIT(ao)
     218           4 :   , isFirstStep_(true)
     219           4 :   , afterCalculate_(false)
     220           4 :   , mean_counter_(0)
     221           4 :   , av_counter_(0)
     222           4 :   , work_(0) {
     223             : //set beta
     224           4 :   const double Kb=getKBoltzmann();
     225           4 :   double KbT=getkBT();
     226           4 :   plumed_massert(KbT>0,"your MD engine does not pass the temperature to plumed, you must specify it using TEMP");
     227           4 :   beta_=1.0/KbT;
     228             : 
     229             : //initialize probability grids using local free energies
     230             :   bool spline=true;
     231             :   bool sparsegrid=false;
     232           4 :   std::string funcl="file.free"; //typical name given by sum_hills
     233             : 
     234             :   std::vector<std::string> fes_names;
     235           8 :   for(unsigned n=0;; n++) { //NB: here we start from FILE_F0 not from FILE_F1
     236             :     std::string filename;
     237          24 :     if(!parseNumbered("FILE_F",n,filename)) {
     238             :       break;
     239             :     }
     240           8 :     fes_names.push_back(filename);
     241           8 :     IFile gridfile;
     242           8 :     gridfile.open(filename);
     243           8 :     auto g=GridBase::create(funcl,getArguments(),gridfile,sparsegrid,spline,true);
     244             : // we assume this cannot be sparse. in case we want it to be sparse, some of the methods
     245             : // that are available only in Grid should be ported to GridBase
     246           8 :     auto gg=dynamic_cast<Grid*>(g.get());
     247             : // if this throws, g is deleted
     248           8 :     plumed_assert(gg);
     249             : // release ownership in order to transfer it to emplaced pointer
     250             : // cppcheck-suppress ignoredReturnValue
     251             :     g.release();
     252           8 :     grid_p_.emplace_back(gg);
     253          16 :   }
     254           4 :   plumed_massert(grid_p_.size()>1,"at least 2 basins must be defined, starting from FILE_F0");
     255           4 :   alpha_size_=grid_p_.size()-1;
     256           4 :   sym_alpha_size_=alpha_size_*(alpha_size_+1)/2; //useful for symmetric matrix [alpha_size_]x[alpha_size_]
     257             :   //check for consistency with first local free energy
     258           8 :   for(unsigned n=1; n<grid_p_.size(); n++) {
     259           8 :     std::string error_tag="FILE_F"+std::to_string(n)+" '"+fes_names[n]+"' not compatible with reference one, FILE_F0";
     260           4 :     plumed_massert(grid_p_[n]->getSize()==grid_p_[0]->getSize(),error_tag);
     261           4 :     plumed_massert(grid_p_[n]->getMin()==grid_p_[0]->getMin(),error_tag);
     262           4 :     plumed_massert(grid_p_[n]->getMax()==grid_p_[0]->getMax(),error_tag);
     263           4 :     plumed_massert(grid_p_[n]->getBinVolume()==grid_p_[0]->getBinVolume(),error_tag);
     264             :   }
     265             : 
     266           4 :   bool no_mintozero=false;
     267           4 :   parseFlag("NO_MINTOZERO",no_mintozero);
     268           4 :   if(!no_mintozero) {
     269           6 :     for(unsigned n=0; n<grid_p_.size(); n++) {
     270           4 :       grid_p_[n]->setMinToZero();
     271             :     }
     272             :   }
     273           4 :   bool normalize=false;
     274           4 :   parseFlag("NORMALIZE",normalize);
     275           4 :   norm_.resize(grid_p_.size(),0);
     276           4 :   std::vector<double> c_norm(grid_p_.size());
     277             :   //convert the FESs to probability distributions
     278             :   //NB: the spline interpolation will be done on the probability distributions, not on the given FESs
     279             :   const unsigned ncv=getNumberOfArguments(); //just for ease
     280          12 :   for(unsigned n=0; n<grid_p_.size(); n++) {
     281         808 :     for(Grid::index_t t=0; t<grid_p_[n]->getSize(); t++) {
     282         800 :       std::vector<double> der(ncv);
     283         800 :       const double val=std::exp(-beta_*grid_p_[n]->getValueAndDerivatives(t,der));
     284        1600 :       for(unsigned s=0; s<ncv; s++) {
     285         800 :         der[s]*=-beta_*val;
     286             :       }
     287         800 :       grid_p_[n]->setValueAndDerivatives(t,val,der);
     288         800 :       norm_[n]+=val;
     289             :     }
     290           8 :     c_norm[n]=1./beta_*std::log(norm_[n]);
     291           8 :     if(normalize) {
     292           4 :       grid_p_[n]->scaleAllValuesAndDerivatives(1./norm_[n]);
     293           4 :       norm_[n]=1;
     294             :     }
     295             :   }
     296             : 
     297             : //get target
     298           4 :   double biasfactor=0;
     299           4 :   parse("BIASFACTOR",biasfactor);
     300           4 :   plumed_massert(biasfactor==0 || biasfactor>1,"BIASFACTOR must be zero (for uniform target) or greater than one");
     301           4 :   if(biasfactor==0) {
     302           2 :     inv_gamma_=0;
     303             :   } else {
     304           2 :     inv_gamma_=1./biasfactor;
     305             :   }
     306           4 :   tg_counter_=0;
     307           4 :   tg_stride_=1;
     308           4 :   parse("TG_STRIDE",tg_stride_);
     309           4 :   tg_dV_dAlpha_.resize(alpha_size_,0);
     310           4 :   tg_d2V_dAlpha2_.resize(sym_alpha_size_,0);
     311             : 
     312             : //setup optimization stuff
     313           4 :   minimization_step_=1;
     314           4 :   parse("M_STEP",minimization_step_);
     315             : 
     316           4 :   av_stride_=500;
     317           4 :   parse("AV_STRIDE",av_stride_);
     318           4 :   av_dV_dAlpha_.resize(alpha_size_,0);
     319           4 :   av_dV_dAlpha_prod_.resize(sym_alpha_size_,0);
     320           4 :   av_d2V_dAlpha2_.resize(sym_alpha_size_,0);
     321             : 
     322           4 :   mean_weight_tau_=0;
     323           4 :   parse("TAU_MEAN",mean_weight_tau_);
     324           4 :   if(mean_weight_tau_!=1) { //set it to 1 for basic SGD
     325           4 :     plumed_massert((mean_weight_tau_==0 || mean_weight_tau_>av_stride_),"TAU_MEAN is rescaled with AV_STRIDE, so it has to be greater");
     326           4 :     mean_weight_tau_/=av_stride_; //this way you can look at the number of simulation steps to choose TAU_MEAN
     327             :   }
     328             : 
     329           8 :   parseVector("INITIAL_ALPHA",mean_alpha_);
     330           4 :   if(mean_alpha_.size()>0) {
     331           2 :     plumed_massert(mean_alpha_.size()==alpha_size_,"provide one INITIAL_ALPHA for each basin beyond the first one");
     332             :   } else {
     333           2 :     mean_alpha_.resize(alpha_size_,0);
     334             :   }
     335           4 :   inst_alpha_=mean_alpha_;
     336           4 :   exp_alpha_.resize(alpha_size_);
     337           8 :   for(unsigned i=0; i<alpha_size_; i++) {
     338           4 :     exp_alpha_[i]=std::exp(-beta_*mean_alpha_[i]);
     339             :   }
     340           4 :   prev_exp_alpha_=exp_alpha_;
     341             : 
     342           4 :   damping_off_=false;
     343           4 :   parseFlag("DAMPING_OFF",damping_off_);
     344           4 :   if(damping_off_) {
     345           2 :     past_increment2_.resize(alpha_size_,1);
     346             :   } else {
     347           2 :     past_increment2_.resize(alpha_size_,0);
     348             :   }
     349             : 
     350             : //file printing options
     351           4 :   std::string alphaFileName("ALPHA");
     352           4 :   parse("ALPHA_FILE",alphaFileName);
     353           4 :   print_stride_=10;
     354           8 :   parse("PRINT_STRIDE",print_stride_);
     355             :   std::string fmt;
     356           4 :   parse("FMT",fmt);
     357             : 
     358             : //other flags, mainly for debugging
     359           4 :   NumParallel_=comm.Get_size();
     360           4 :   rank_=comm.Get_rank();
     361           4 :   bool serial=false;
     362           4 :   parseFlag("SERIAL",serial);
     363           4 :   if(serial) {
     364           2 :     log.printf(" -- SERIAL: running without loop parallelization\n");
     365           2 :     NumParallel_=1;
     366           2 :     rank_=0;
     367             :   }
     368             : 
     369           4 :   bool multiple_walkers=false;
     370           4 :   parseFlag("MULTIPLE_WALKERS",multiple_walkers);
     371           4 :   if(!multiple_walkers) {
     372           2 :     NumWalkers_=1;
     373             :   } else {
     374           2 :     if(comm.Get_rank()==0) { //multi_sim_comm works well on first rank only
     375           2 :       NumWalkers_=multi_sim_comm.Get_size();
     376             :     }
     377           2 :     if(comm.Get_size()>1) { //if each walker has more than one processor update them all
     378           0 :       comm.Bcast(NumWalkers_,0);
     379             :     }
     380             :   }
     381             : 
     382           4 :   checkRead();
     383             : 
     384             : //restart if needed
     385           4 :   if(getRestart()) {
     386           2 :     IFile ifile;
     387           2 :     ifile.link(*this);
     388           2 :     if(NumWalkers_>1) {
     389           4 :       ifile.enforceSuffix("");
     390             :     }
     391           2 :     if(ifile.FileExist(alphaFileName)) {
     392           2 :       log.printf("  Restarting from: %s\n",alphaFileName.c_str());
     393           2 :       log.printf("    all options (also PRINT_STRIDE) must be consistent!\n");
     394           2 :       log.printf("    any INITIAL_ALPHA will be overwritten\n");
     395           2 :       ifile.open(alphaFileName);
     396             :       double time;
     397           2 :       std::vector<double> damping(alpha_size_);
     398          20 :       while(ifile.scanField("time",time)) { //room for improvements: only last line is important
     399          16 :         for(unsigned i=0; i<alpha_size_; i++) {
     400           8 :           const std::string index(std::to_string(i+1));
     401           8 :           prev_exp_alpha_[i]=std::exp(-beta_*mean_alpha_[i]);
     402          16 :           ifile.scanField("alpha_"+index,mean_alpha_[i]);
     403          16 :           ifile.scanField("auxiliary_"+index,inst_alpha_[i]);
     404          16 :           ifile.scanField("damping_"+index,damping[i]);
     405             :         }
     406           8 :         ifile.scanField();
     407           8 :         mean_counter_+=print_stride_;
     408             :       }
     409           4 :       for(unsigned i=0; i<alpha_size_; i++) {
     410           2 :         exp_alpha_[i]=std::exp(-beta_*mean_alpha_[i]);
     411           2 :         past_increment2_[i]=damping[i]*damping[i];
     412             :       }
     413             :       //sync all walkers and treads. Not sure is mandatory but is no harm
     414           2 :       comm.Barrier();
     415           2 :       if(comm.Get_rank()==0) {
     416           2 :         multi_sim_comm.Barrier();
     417             :       }
     418             :     } else {
     419           0 :       log.printf("  -- WARNING: restart requested, but no '%s' file found!\n",alphaFileName.c_str());
     420             :     }
     421           2 :   }
     422             : 
     423             : //setup output file with Alpha values
     424           4 :   alphaOfile_.link(*this);
     425           4 :   if(NumWalkers_>1) {
     426           2 :     if(comm.Get_rank()==0 && multi_sim_comm.Get_rank()>0) {
     427             :       alphaFileName="/dev/null";  //only first walker writes on file
     428             :     }
     429           4 :     alphaOfile_.enforceSuffix("");
     430             :   }
     431           4 :   alphaOfile_.open(alphaFileName);
     432           4 :   if(fmt.length()>0) {
     433           8 :     alphaOfile_.fmtField(" "+fmt);
     434             :   }
     435             : 
     436             : //add other output components
     437           8 :   addComponent("rct");
     438           8 :   componentIsNotPeriodic("rct");
     439           8 :   addComponent("work");
     440           4 :   componentIsNotPeriodic("work");
     441             : 
     442             : //print some info
     443           4 :   log.printf("  Temperature T: %g\n",1./(Kb*beta_));
     444           4 :   log.printf("  Beta (1/Kb*T): %g\n",beta_);
     445           4 :   log.printf("  Local free energy basins files and normalization constants:\n");
     446          12 :   for(unsigned n=0; n<grid_p_.size(); n++) {
     447           8 :     log.printf("    F_%d filename: %s  c_%d=%g\n",n,fes_names[n].c_str(),n,c_norm[n]);
     448             :   }
     449           4 :   if(no_mintozero) {
     450           2 :     log.printf(" -- NO_MINTOZERO: local free energies are not shifted to be zero at minimum\n");
     451             :   }
     452           4 :   if(normalize) {
     453           2 :     log.printf(" -- NORMALIZE: F_n+=c_n, alpha=DeltaF\n");
     454             :   }
     455           4 :   log.printf("  Using target distribution with 1/gamma = %g\n",inv_gamma_);
     456           4 :   log.printf("    and updated with stride %d\n",tg_stride_);
     457           4 :   log.printf("  Step for the minimization algorithm: %g\n",minimization_step_);
     458           4 :   log.printf("  Stride for the ensemble average: %d\n",av_stride_);
     459           4 :   if(mean_weight_tau_>1) {
     460           2 :     log.printf("  Exponentially decaying average with weight=tau/av_stride=%d\n",mean_weight_tau_);
     461             :   }
     462           4 :   if(mean_weight_tau_==1) {
     463           0 :     log.printf(" +++ WARNING +++ setting TAU_MEAN=1 is equivalent to use simple SGD, without mean alpha nor hessian contribution\n");
     464             :   }
     465           4 :   log.printf("  Initial guess for alpha:\n");
     466           8 :   for(unsigned i=0; i<alpha_size_; i++) {
     467           4 :     log.printf("    alpha_%d = %g\n",i+1,mean_alpha_[i]);
     468             :   }
     469           4 :   if(damping_off_) {
     470           2 :     log.printf(" -- DAMPING_OFF: the minimization step will NOT become smaller as the simulation goes on\n");
     471             :   }
     472           4 :   log.printf("  Printing on file %s with stride %d\n",alphaFileName.c_str(),print_stride_);
     473           4 :   if(serial) {
     474           2 :     log.printf(" -- SERIAL: running without loop parallelization\n");
     475             :   }
     476           4 :   if(NumParallel_>1) {
     477           2 :     log.printf("  Using multiple threads per simulation: %d\n",NumParallel_);
     478             :   }
     479           4 :   if(multiple_walkers) {
     480           2 :     log.printf(" -- MULTIPLE_WALKERS: multiple simulations will combine statistics for the optimization\n");
     481           2 :     if(NumWalkers_>1) {
     482           2 :       log.printf("    number of walkers: %d\n",NumWalkers_);
     483           2 :       log.printf("    walker rank: %d\n",multi_sim_comm.Get_rank()); //only comm.Get_rank()=0 prints, so this is fine
     484             :     } else {
     485           0 :       log.printf(" +++ WARNING +++ only one replica found: are you sure you are running MPI-connected simulations?\n");
     486             :     }
     487             :   }
     488           4 :   log.printf(" Bibliography ");
     489           8 :   log<<plumed.cite("Invernizzi and Parrinello, J. Chem. Theory Comput. 15, 2187-2194 (2019)");
     490           8 :   log<<plumed.cite("Valsson and Parrinello, Phys. Rev. Lett. 113, 090601 (2014)");
     491           4 :   if(inv_gamma_>0) {
     492           4 :     log<<plumed.cite("Valsson and Parrinello, J. Chem. Theory Comput. 11, 1996-2002 (2015)");
     493             :   }
     494             : 
     495             : //last initializations
     496           4 :   prob_.resize(grid_p_.size());
     497           4 :   der_prob_.resize(grid_p_.size(),std::vector<double>(getNumberOfArguments()));
     498           4 :   update_tg_and_rct();
     499           8 : }
     500             : 
     501         804 : void VesDeltaF::calculate() {
     502             : //get CVs
     503         804 :   const unsigned ncv=getNumberOfArguments(); //just for ease
     504         804 :   std::vector<double> cv(ncv);
     505        1608 :   for(unsigned s=0; s<ncv; s++) {
     506         804 :     cv[s]=getArgument(s);
     507             :   }
     508             : //get probabilities for each basin, and total one
     509        2412 :   for(unsigned n=0; n<grid_p_.size(); n++) {
     510        1608 :     prob_[n]=grid_p_[n]->getValueAndDerivatives(cv,der_prob_[n]);
     511             :   }
     512         804 :   tot_prob_=prob_[0];
     513        1608 :   for(unsigned i=0; i<alpha_size_; i++) {
     514         804 :     tot_prob_+=prob_[i+1]*exp_alpha_[i];
     515             :   }
     516             : 
     517             : //update bias and forces: V=-(1-inv_gamma_)*fes
     518         804 :   setBias((1-inv_gamma_)/beta_*std::log(tot_prob_));
     519        1608 :   for(unsigned s=0; s<ncv; s++) {
     520         804 :     double dProb_dCV_s=der_prob_[0][s];
     521        1608 :     for(unsigned i=0; i<alpha_size_; i++) {
     522         804 :       dProb_dCV_s+=der_prob_[i+1][s]*exp_alpha_[i];
     523             :     }
     524         804 :     setOutputForce(s,-(1-inv_gamma_)/beta_/tot_prob_*dProb_dCV_s);
     525             :   }
     526         804 :   afterCalculate_=true;
     527         804 : }
     528             : 
     529         804 : void VesDeltaF::update() {
     530             : //skip first step to sync getTime() and av_counter_, as in METAD
     531         804 :   if(isFirstStep_) {
     532           4 :     isFirstStep_=false;
     533           4 :     return;
     534             :   }
     535         800 :   plumed_massert(afterCalculate_,"VesDeltaF::update() must be called after VesDeltaF::calculate() to work properly");
     536         800 :   afterCalculate_=false;
     537             : 
     538             : //calculate derivatives for ensemble averages
     539         800 :   std::vector<double> dV_dAlpha(alpha_size_);
     540         800 :   std::vector<double> d2V_dAlpha2(sym_alpha_size_);
     541        1600 :   for(unsigned i=0; i<alpha_size_; i++) {
     542         800 :     dV_dAlpha[i]=-(1-inv_gamma_)/tot_prob_*prob_[i+1]*exp_alpha_[i];
     543             :   }
     544        1600 :   for(unsigned i=0; i<alpha_size_; i++) {
     545         800 :     d2V_dAlpha2[get_index(i,i)]=-beta_*dV_dAlpha[i];
     546        1600 :     for(unsigned j=i; j<alpha_size_; j++) {
     547         800 :       d2V_dAlpha2[get_index(i,j)]-=beta_/(1-inv_gamma_)*dV_dAlpha[i]*dV_dAlpha[j];
     548             :     }
     549             :   }
     550             : //update ensemble averages
     551         800 :   av_counter_++;
     552        1600 :   for(unsigned i=0; i<alpha_size_; i++) {
     553         800 :     av_dV_dAlpha_[i]+=(dV_dAlpha[i]-av_dV_dAlpha_[i])/av_counter_;
     554        1600 :     for(unsigned j=i; j<alpha_size_; j++) {
     555         800 :       const unsigned ij=get_index(i,j);
     556         800 :       av_dV_dAlpha_prod_[ij]+=(dV_dAlpha[i]*dV_dAlpha[j]-av_dV_dAlpha_prod_[ij])/av_counter_;
     557         800 :       av_d2V_dAlpha2_[ij]+=(d2V_dAlpha2[ij]-av_d2V_dAlpha2_[ij])/av_counter_;
     558             :     }
     559             :   }
     560             : //update work
     561         800 :   double prev_tot_prob=prob_[0];
     562        1600 :   for(unsigned i=0; i<alpha_size_; i++) {
     563         800 :     prev_tot_prob+=prob_[i+1]*prev_exp_alpha_[i];
     564             :   }
     565         800 :   work_+=(1-inv_gamma_)/beta_*std::log(tot_prob_/prev_tot_prob);
     566             : 
     567             : //update coefficients
     568         800 :   if(av_counter_==av_stride_) {
     569          16 :     update_alpha();
     570          16 :     tg_counter_++;
     571          16 :     if(tg_counter_==tg_stride_) {
     572          12 :       update_tg_and_rct();
     573          12 :       tg_counter_=0;
     574             :     }
     575             :     //reset the ensemble averages
     576          16 :     av_counter_=0;
     577             :     std::fill(av_dV_dAlpha_.begin(),av_dV_dAlpha_.end(),0);
     578             :     std::fill(av_dV_dAlpha_prod_.begin(),av_dV_dAlpha_prod_.end(),0);
     579             :     std::fill(av_d2V_dAlpha2_.begin(),av_d2V_dAlpha2_.end(),0);
     580             :   }
     581             : }
     582             : 
     583          16 : void VesDeltaF::update_tg_and_rct() {
     584             : //calculate target averages
     585          16 :   double Z_0=norm_[0];
     586          32 :   for(unsigned i=0; i<alpha_size_; i++) {
     587          16 :     Z_0+=norm_[i+1]*exp_alpha_[i];
     588             :   }
     589          16 :   double Z_tg=0;
     590             :   std::fill(tg_dV_dAlpha_.begin(),tg_dV_dAlpha_.end(),0);
     591             :   std::fill(tg_d2V_dAlpha2_.begin(),tg_d2V_dAlpha2_.end(),0);
     592        1116 :   for(Grid::index_t t=rank_; t<grid_p_[0]->getSize(); t+=NumParallel_) {
     593             :     //TODO can we recycle some code?
     594        1100 :     std::vector<double> prob(grid_p_.size());
     595        3300 :     for(unsigned n=0; n<grid_p_.size(); n++) {
     596        2200 :       prob[n]=grid_p_[n]->getValue(t);
     597             :     }
     598        1100 :     double tot_prob=prob[0];
     599        2200 :     for(unsigned i=0; i<alpha_size_; i++) {
     600        1100 :       tot_prob+=prob[i+1]*exp_alpha_[i];
     601             :     }
     602        1100 :     std::vector<double> dV_dAlpha(alpha_size_);
     603        1100 :     std::vector<double> d2V_dAlpha2(sym_alpha_size_);
     604        2200 :     for(unsigned i=0; i<alpha_size_; i++) {
     605        1100 :       dV_dAlpha[i]=-(1-inv_gamma_)/tot_prob*prob[i+1]*exp_alpha_[i];
     606             :     }
     607        2200 :     for(unsigned i=0; i<alpha_size_; i++) {
     608        1100 :       d2V_dAlpha2[get_index(i,i)]=-beta_*dV_dAlpha[i];
     609        2200 :       for(unsigned j=i; j<alpha_size_; j++) {
     610        1100 :         d2V_dAlpha2[get_index(i,j)]-=beta_/(1-inv_gamma_)*dV_dAlpha[i]*dV_dAlpha[j];
     611             :       }
     612             :     }
     613        1100 :     const double unnorm_tg_p=std::pow(tot_prob,inv_gamma_);
     614        1100 :     Z_tg+=unnorm_tg_p;
     615        2200 :     for(unsigned i=0; i<alpha_size_; i++) {
     616        1100 :       tg_dV_dAlpha_[i]+=unnorm_tg_p*dV_dAlpha[i];
     617             :     }
     618        2200 :     for(unsigned ij=0; ij<sym_alpha_size_; ij++) {
     619        1100 :       tg_d2V_dAlpha2_[ij]+=unnorm_tg_p*d2V_dAlpha2[ij];
     620             :     }
     621             :   }
     622          16 :   if(NumParallel_>1) {
     623          10 :     comm.Sum(Z_tg);
     624          10 :     comm.Sum(tg_dV_dAlpha_);
     625          10 :     comm.Sum(tg_d2V_dAlpha2_);
     626             :   }
     627          32 :   for(unsigned i=0; i<alpha_size_; i++) {
     628          16 :     tg_dV_dAlpha_[i]/=Z_tg;
     629             :   }
     630          32 :   for(unsigned ij=0; ij<sym_alpha_size_; ij++) {
     631          16 :     tg_d2V_dAlpha2_[ij]/=Z_tg;
     632             :   }
     633          16 :   getPntrToComponent("rct")->set(-1./beta_*std::log(Z_tg/Z_0)); //Z_tg is the best available estimate of Z_V
     634          16 : }
     635             : 
     636          16 : void VesDeltaF::update_alpha() {
     637             : //combining the averages of multiple walkers
     638          16 :   if(NumWalkers_>1) {
     639           8 :     if(comm.Get_rank()==0) { //sum only once: in the first rank of each walker
     640           8 :       multi_sim_comm.Sum(av_dV_dAlpha_);
     641           8 :       multi_sim_comm.Sum(av_dV_dAlpha_prod_);
     642           8 :       multi_sim_comm.Sum(av_d2V_dAlpha2_);
     643          16 :       for(unsigned i=0; i<alpha_size_; i++) {
     644           8 :         av_dV_dAlpha_[i]/=NumWalkers_;
     645             :       }
     646          16 :       for(unsigned ij=0; ij<sym_alpha_size_; ij++) {
     647           8 :         av_dV_dAlpha_prod_[ij]/=NumWalkers_;
     648           8 :         av_d2V_dAlpha2_[ij]/=NumWalkers_;
     649             :       }
     650             :     }
     651           8 :     if(comm.Get_size()>1) { //if there are more ranks for each walker, everybody has to know
     652           0 :       comm.Bcast(av_dV_dAlpha_,0);
     653           0 :       comm.Bcast(av_dV_dAlpha_prod_,0);
     654           0 :       comm.Bcast(av_d2V_dAlpha2_,0);
     655             :     }
     656             :   }
     657             :   //set work and reset it
     658          16 :   getPntrToComponent("work")->set(work_);
     659          16 :   work_=0;
     660             : 
     661             : //build the gradient and the Hessian of the functional
     662          16 :   std::vector<double> grad_omega(alpha_size_);
     663          16 :   std::vector<double> hess_omega(sym_alpha_size_);
     664          32 :   for(unsigned i=0; i<alpha_size_; i++) {
     665          16 :     grad_omega[i]=tg_dV_dAlpha_[i]-av_dV_dAlpha_[i];
     666          32 :     for(unsigned j=i; j<alpha_size_; j++) {
     667          16 :       const unsigned ij=get_index(i,j);
     668          16 :       hess_omega[ij]=beta_*(av_dV_dAlpha_prod_[ij]-av_dV_dAlpha_[i]*av_dV_dAlpha_[j])+tg_d2V_dAlpha2_[ij]-av_d2V_dAlpha2_[ij];
     669             :     }
     670             :   }
     671             : //calculate the increment and update alpha
     672          16 :   mean_counter_++;
     673             :   long long unsigned mean_weight=mean_counter_;
     674          16 :   if(mean_weight_tau_>0 && mean_weight_tau_<mean_counter_) {
     675             :     mean_weight=mean_weight_tau_;
     676             :   }
     677          16 :   std::vector<double> damping(alpha_size_);
     678          32 :   for(unsigned i=0; i<alpha_size_; i++) {
     679          16 :     double increment_i=grad_omega[i];
     680          32 :     for(unsigned j=0; j<alpha_size_; j++) {
     681          16 :       increment_i+=hess_omega[get_index(i,j)]*(inst_alpha_[j]-mean_alpha_[j]);
     682             :     }
     683          16 :     if(!damping_off_) {
     684           8 :       past_increment2_[i]+=increment_i*increment_i;
     685             :     }
     686          16 :     damping[i]=std::sqrt(past_increment2_[i]);
     687          16 :     prev_exp_alpha_[i]=std::exp(-beta_*mean_alpha_[i]);
     688          16 :     inst_alpha_[i]-=minimization_step_/damping[i]*increment_i;
     689          16 :     mean_alpha_[i]+=(inst_alpha_[i]-mean_alpha_[i])/mean_weight;
     690          16 :     exp_alpha_[i]=std::exp(-beta_*mean_alpha_[i]);
     691             :   }
     692             : 
     693             : //update the Alpha file
     694          16 :   if(mean_counter_%print_stride_==0) {
     695          16 :     alphaOfile_.printField("time",getTime());
     696          32 :     for(unsigned i=0; i<alpha_size_; i++) {
     697          16 :       const std::string index(std::to_string(i+1));
     698          32 :       alphaOfile_.printField("alpha_"+index,mean_alpha_[i]);
     699          32 :       alphaOfile_.printField("auxiliary_"+index,inst_alpha_[i]);
     700          32 :       alphaOfile_.printField("damping_"+index,damping[i]);
     701             :     }
     702          16 :     alphaOfile_.printField();
     703             :   }
     704          16 : }
     705             : 
     706             : //mapping of a [alpha_size_]x[alpha_size_] symmetric matrix into a vector of size sym_alpha_size_, useful for the communicator
     707        4632 : inline unsigned VesDeltaF::get_index(const unsigned i, const unsigned j) const {
     708        4632 :   if(i<=j) {
     709        4632 :     return j+i*(alpha_size_-1)-i*(i-1)/2;
     710             :   } else {
     711           0 :     return get_index(j,i);
     712             :   }
     713             : }
     714             : 
     715             : }
     716             : }

Generated by: LCOV version 1.16