LCOV - code coverage report
Current view: top level - eds - EDS.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 452 488 92.6 %
Date: 2026-03-30 13:16:06 Functions: 20 23 87.0 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             : Copyright (c) 2017 of Glen Hocky and Andrew White
       3             : 
       4             : The eds module is free software: you can redistribute it and/or modify
       5             : it under the terms of the GNU Lesser General Public License as published by
       6             : the Free Software Foundation, either version 3 of the License, or
       7             : (at your option) any later version.
       8             : 
       9             : The eds module is distributed in the hope that it will be useful,
      10             : but WITHOUT ANY WARRANTY; without even the implied warranty of
      11             : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
      12             : GNU Lesser General Public License for more details.
      13             : 
      14             : You should have received a copy of the GNU Lesser General Public License
      15             : along with plumed.  If not, see <http://www.gnu.org/licenses/>.
      16             : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
      17             : #include "bias/Bias.h"
      18             : #include "bias/ReweightBase.h"
      19             : #include "core/ActionAtomistic.h"
      20             : #include "core/ActionRegister.h"
      21             : #include "core/Atoms.h"
      22             : #include "core/PlumedMain.h"
      23             : #include "tools/File.h"
      24             : #include "tools/Matrix.h"
      25             : #include "tools/Random.h"
      26             : 
      27             : #include <iostream>
      28             : 
      29             : using namespace PLMD;
      30             : using namespace bias;
      31             : 
      32             : // namespace is lowercase to match
      33             : // module names being all lowercase
      34             : 
      35             : namespace PLMD {
      36             : namespace eds {
      37             : 
      38             : //+PLUMEDOC EDSMOD_BIAS EDS
      39             : /*
      40             : Add a linear bias on a set of observables.
      41             : 
      42             : This force is the same as the linear part of the bias in \ref
      43             : RESTRAINT, but this bias has the ability to compute the prefactors
      44             : adaptively using the scheme of White and Voth \cite white2014efficient
      45             : in order to match target observable values for a set of CVs.
      46             : Further updates to the algorithm are described in \cite hocky2017cgds
      47             : and you can read a review on the method and its applications here: \cite Amirkulova2019Recent.
      48             : 
      49             : You can
      50             : see a tutorial on EDS specifically for biasing coordination number at
      51             : <a
      52             : href="http://thewhitelab.org/blog/tutorial/2017/05/10/lammps-coordination-number-tutorial/">
      53             : Andrew White's webpage</a>.
      54             : 
      55             : The addition to the potential is of the form
      56             : \f[
      57             :   \sum_i \frac{\alpha_i}{s_i} x_i
      58             : \f]
      59             : 
      60             : where for CV \f$x_i\f$, a coupling constant \f${\alpha}_i\f$ is determined
      61             : adaptively or set by the user to match a target value for
      62             : \f$x_i\f$. \f$s_i\f$ is a scale parameter, which by default is set to
      63             : the target value. It may also be set separately.
      64             : 
      65             : \warning
      66             : It is not possible to set the target value of the observable
      67             : to zero with the default value of \f$s_i\f$ as this will cause a
      68             : divide-by-zero error. Instead, set \f$s_i=1\f$ or modify the CV so the
      69             : desired target value is no longer zero.
      70             : 
      71             : Notice that a similar method is available as \ref MAXENT, although with different features and using a different optimization algorithm.
      72             : 
      73             : \par Virial
      74             : 
      75             : The bias forces modify the virial and this can change your simulation density if the bias is used in an NPT simulation.
      76             : One way to avoid changing the virial contribution from the bias is to add the keyword VIRIAL=1.0, which changes how the bias
      77             : is computed to minimize its contribution to the virial. This can also lead to smaller magnitude biases that are more robust if
      78             : transferred to other systems.  VIRIAL=1.0 can be a reasonable starting point and increasing the value changes the balance between matching
      79             : the set-points and minimizing the virial. See \cite Amirkulova2019Recent for details on the equations. Since the coupling constants
      80             : are unique with a single CV, VIRIAL is not applicable with a single CV. When used with multiple CVs, the CVs should be correlated
      81             : which is almost always the case.
      82             : 
      83             : \par Weighting
      84             : 
      85             : EDS computes means and variances as part of its algorithm. If you are
      86             : also using a biasing method like metadynamics, you may wish to remove
      87             : the effect of this bias in your EDS computations so that EDS works on
      88             : the canonical values (reweighted to be unbiased).  For example, you may be using
      89             : metadynamics to bias a dihedral angle to enhance sampling and be using
      90             : EDS to set the average distance between two particular atoms. Specifically:
      91             : 
      92             : \plumedfile
      93             : # set-up metadynamics
      94             : t: TORSION ATOMS=1,2,3,4
      95             : md: METAD ARG=d SIGMA=0.2 HEIGHT=0.3 PACE=500 TEMP=300
      96             : # compute bias weights
      97             : bias: REWEIGHT_METAD TEMP=300
      98             : # now do EDS on distance while removing effect of metadynamics
      99             : d: DISTANCE ATOMS=4,7
     100             : eds: EDS ARG=d CENTER=3.0 PERIOD=100 TEMP=300 LOGWEIGHTS=bias
     101             : \endplumedfile
     102             : 
     103             : This is an approximation though because EDS uses a finite samples while running to get means/variances.
     104             : At the end of a run,
     105             : you should ensure this approach worked and indeed your reweighted CV matches the target value.
     106             : 
     107             : \par Examples
     108             : 
     109             : The following input for a harmonic oscillator of two beads will
     110             : adaptively find a linear bias to change the mean and variance to the
     111             : target values. The PRINT line shows how to access the value of the
     112             : coupling constants.
     113             : 
     114             : \plumedfile
     115             : dist: DISTANCE ATOMS=1,2
     116             : # this is the squared of the distance
     117             : dist2: COMBINE ARG=dist POWERS=2 PERIODIC=NO
     118             : 
     119             : # bias mean and variance
     120             : eds: EDS ARG=dist,dist2 CENTER=2.0,1.0 PERIOD=100 TEMP=1.0
     121             : PRINT ARG=dist,dist2,eds.dist_coupling,eds.dist2_coupling,eds.bias,eds.force2 FILE=colvars.dat STRIDE=100
     122             : \endplumedfile
     123             : 
     124             : <hr>
     125             : 
     126             : Rather than trying to find the coupling constants adaptively, one can ramp up to a constant value.
     127             : \plumedfile
     128             : dist: DISTANCE ATOMS=1,2
     129             : dist2: COMBINE ARG=dist POWERS=2 PERIODIC=NO
     130             : 
     131             : # ramp couplings from 0,0 to -1,1 over 50000 steps
     132             : eds: EDS ARG=dist,dist2 CENTER=2.0,1.0 FIXED=-1,1 RAMP PERIOD=50000 TEMP=1.0
     133             : 
     134             : # same as above, except starting at -0.5,0.5 rather than default of 0,0
     135             : eds2: EDS ARG=dist,dist2 CENTER=2.0,1.0 FIXED=-1,1 INIT=-0.5,0.5 RAMP PERIOD=50000 TEMP=1.0
     136             : \endplumedfile
     137             : 
     138             : <hr>
     139             : A restart file can be added to dump information needed to restart/continue simulation using these parameters every PERIOD.
     140             : \plumedfile
     141             : dist: DISTANCE ATOMS=1,2
     142             : dist2: COMBINE ARG=dist POWERS=2 PERIODIC=NO
     143             : 
     144             : # add the option to write to a restart file
     145             : eds: EDS ARG=dist,dist2 CENTER=2.0,1.0 PERIOD=100 TEMP=1.0 OUT_RESTART=checkpoint.eds
     146             : \endplumedfile
     147             : 
     148             : The first few lines of the restart file that is output if we run a calculation with one CV will look something like this:
     149             : 
     150             : \auxfile{restart.eds}
     151             : #! FIELDS time d1_center d1_set d1_target d1_coupling d1_maxrange d1_maxgrad d1_accum d1_mean d1_std
     152             : #! SET adaptive  1
     153             : #! SET update_period  1
     154             : #! SET seed  0
     155             : #! SET kbt    2.4943
     156             :    0.0000   1.0000   0.0000   0.0000   0.0000   7.4830   0.1497   0.0000   0.0000   0.0000
     157             :    1.0000   1.0000   0.0000   0.0000   0.0000   7.4830   0.1497   0.0000   0.0000   0.0000
     158             :    2.0000   1.0000  -7.4830   0.0000   0.0000   7.4830   0.1497   0.0224   0.0000   0.0000
     159             :    3.0000   1.0000  -7.4830   0.0000  -7.4830   7.4830   0.1497   0.0224   0.0000   0.0000
     160             :    4.0000   1.0000  -7.4830   0.0000  -7.4830   7.4830   0.1497   0.0224   0.0000   0.0000
     161             : \endauxfile
     162             : 
     163             : <hr>
     164             : 
     165             : Read in a previous restart file. Adding RESTART flag makes output append
     166             : \plumedfile
     167             : d1: DISTANCE ATOMS=1,2
     168             : 
     169             : eds: EDS ARG=d1 CENTER=2.0 PERIOD=100 TEMP=1.0 IN_RESTART=restart.eds RESTART=YES
     170             : \endplumedfile
     171             : 
     172             : <hr>
     173             : 
     174             : Read in a previous restart file and freeze the bias at the final level from the previous simulation
     175             : \plumedfile
     176             : d1: DISTANCE ATOMS=1,2
     177             : 
     178             : eds: EDS ARG=d1 CENTER=2.0 TEMP=1.0 IN_RESTART=restart.eds FREEZE
     179             : \endplumedfile
     180             : 
     181             : <hr>
     182             : 
     183             : Read in a previous restart file and freeze the bias at the mean from the previous simulation
     184             : \plumedfile
     185             : d1: DISTANCE ATOMS=1,2
     186             : 
     187             : eds: EDS ARG=d1 CENTER=2.0 TEMP=1.0 IN_RESTART=restart.eds FREEZE MEAN
     188             : \endplumedfile
     189             : 
     190             : 
     191             : */
     192             : //+ENDPLUMEDOC
     193             : 
     194             : class EDS : public Bias {
     195             : 
     196             : private:
     197             :   /*We will get this and store it once, since on-the-fly changing number of CVs will be fatal*/
     198             :   const unsigned int ncvs_;
     199             :   std::vector<double> center_;
     200             :   std::vector<Value *> center_values_;
     201             :   ReweightBase *logweights_; // weights to use if reweighting averages
     202             :   std::vector<double> scale_;
     203             :   std::vector<double> current_coupling_;   // actually current coupling
     204             :   std::vector<double> set_coupling_;       // what our coupling is ramping up to. Equal to current_coupling when gathering stats
     205             :   std::vector<double> target_coupling_;    // used when loaded to reach a value
     206             :   std::vector<double> max_coupling_range_; // used for adaptive range
     207             :   std::vector<double> max_coupling_grad_;  // maximum allowed gradient
     208             :   std::vector<double> coupling_rate_;
     209             :   std::vector<double> coupling_accum_;
     210             :   std::vector<double> means_;
     211             :   std::vector<double> differences_;
     212             :   std::vector<double> alpha_vector_;
     213             :   std::vector<double> alpha_vector_2_;
     214             :   std::vector<double> ssds_;
     215             :   std::vector<double> step_size_;
     216             :   std::vector<double> pseudo_virial_;
     217             :   std::vector<Value *> out_coupling_;
     218             :   Matrix<double> covar_;
     219             :   Matrix<double> covar2_;
     220             :   Matrix<double> lm_inv_;
     221             :   std::string in_restart_name_;
     222             :   std::string out_restart_name_;
     223             :   std::string fmt_;
     224             :   OFile out_restart_;
     225             :   IFile in_restart_;
     226             :   bool b_c_values_;
     227             :   bool b_adaptive_;
     228             :   bool b_freeze_;
     229             :   bool b_equil_;
     230             :   bool b_ramp_;
     231             :   bool b_covar_;
     232             :   bool b_restart_;
     233             :   bool b_write_restart_;
     234             :   bool b_lm_;
     235             :   bool b_virial_;
     236             :   bool b_update_virial_;
     237             :   bool b_weights_;
     238             :   int seed_;
     239             :   int update_period_;
     240             :   int avg_coupling_count_;
     241             :   int update_calls_;
     242             :   double kbt_;
     243             :   double multi_prop_;
     244             :   double lm_mixing_par_;
     245             :   double virial_scaling_;
     246             :   double pseudo_virial_sum_; // net virial for all cvs in current period
     247             :   double max_logweight_;     // maximum observed max logweight for period
     248             :   double wsum_;              // sum of weights thus far
     249             :   Random rand_;
     250             :   Value *value_force2_;
     251             :   Value *value_pressure_;
     252             : 
     253             :   /*read input restart. b_mean sets if we use mean or final value for freeze*/
     254             :   void readInRestart(const bool b_mean);
     255             :   /*setup output restart*/
     256             :   void setupOutRestart();
     257             :   /*write output restart*/
     258             :   void writeOutRestart();
     259             :   void update_statistics();
     260             :   void update_pseudo_virial();
     261             :   void calc_lm_step_size();
     262             :   void calc_covar_step_size();
     263             :   void calc_ssd_step_size();
     264             :   void reset_statistics();
     265             :   void update_bias();
     266             :   void apply_bias();
     267             : 
     268             : public:
     269             :   explicit EDS(const ActionOptions &);
     270             :   void calculate();
     271             :   void update();
     272             :   void turnOnDerivatives();
     273             :   static void registerKeywords(Keywords &keys);
     274             :   ~EDS();
     275             : };
     276             : 
     277       13801 : PLUMED_REGISTER_ACTION(EDS, "EDS")
     278             : 
     279          12 : void EDS::registerKeywords(Keywords &keys) {
     280          12 :   Bias::registerKeywords(keys);
     281          12 :   keys.use("ARG");
     282          24 :   keys.add("optional", "CENTER", "The desired centers (equilibrium values) which will be sought during the adaptive linear biasing. This is for fixed centers");
     283          24 :   keys.add("optional", "CENTER_ARG", "The desired centers (equilibrium values) which will be sought during the adaptive linear biasing. "
     284             :            "CENTER_ARG is for calculated centers, e.g. from a CV or analysis. ");
     285             : 
     286          24 :   keys.add("optional", "PERIOD", "Steps over which to adjust bias for adaptive or ramping");
     287          24 :   keys.add("compulsory", "RANGE", "25.0", "The (starting) maximum increase in coupling constant per PERIOD (in \\f$k_B T\\f$/[BIAS_SCALE unit]) for each CV biased");
     288          24 :   keys.add("compulsory", "SEED", "0", "Seed for random order of changing bias");
     289          24 :   keys.add("compulsory", "INIT", "0", "Starting value for coupling constant");
     290          24 :   keys.add("compulsory", "FIXED", "0", "Fixed target values for coupling constant. Non-adaptive.");
     291          24 :   keys.add("optional", "BIAS_SCALE", "A divisor to set the units of the bias. "
     292             :            "If not set, this will be the CENTER value by default (as is done in White and Voth 2014).");
     293          24 :   keys.add("optional", "TEMP", "The system temperature. If not provided will be taken from MD code (if available)");
     294          24 :   keys.add("optional", "MULTI_PROP", "What proportion of dimensions to update at each step. "
     295             :            "Must be in interval [1,0), where 1 indicates all and any other indicates a stochastic update. "
     296             :            "If not set, default is 1 / N, where N is the number of CVs. ");
     297          24 :   keys.add("optional", "VIRIAL", "Add an update penalty for having non-zero virial contributions. Only makes sense with multiple correlated CVs.");
     298          24 :   keys.add("optional", "LOGWEIGHTS", "Add weights to use for computing statistics. For example, if biasing with metadynamics.");
     299          24 :   keys.addFlag("LM", false, "Use Levenberg-Marquadt algorithm along with simultaneous keyword. Otherwise use gradient descent.");
     300          24 :   keys.add("compulsory", "LM_MIXING", "1", "Initial mixing parameter when using Levenberg-Marquadt minimization.");
     301          24 :   keys.add("optional", "RESTART_FMT", "the format that should be used to output real numbers in EDS restarts");
     302          24 :   keys.add("optional", "OUT_RESTART", "Output file for all information needed to continue EDS simulation. "
     303             :            "If you have the RESTART directive set (global or for EDS), this file will be appended to. "
     304             :            "Note that the header will be printed again if appending.");
     305          24 :   keys.add("optional", "IN_RESTART", "Read this file to continue an EDS simulation. "
     306             :            "If same as OUT_RESTART and you have not set the RESTART directive, the file will be backed-up and overwritten with new output. "
     307             :            "If you do have the RESTART flag set and it is the same name as OUT_RESTART, this file will be appended.");
     308             : 
     309          24 :   keys.addFlag("RAMP", false, "Slowly increase bias constant to a fixed value");
     310          24 :   keys.addFlag("COVAR", false, "Utilize the covariance matrix when updating the bias. Default Off, but may be enabled due to other options");
     311          24 :   keys.addFlag("FREEZE", false, "Fix bias at current level (only used for restarting).");
     312          24 :   keys.addFlag("MEAN", false, "Instead of using final bias level from restart, use average. Can only be used in conjunction with FREEZE");
     313             : 
     314          12 :   keys.use("RESTART");
     315             : 
     316          24 :   keys.addOutputComponent("force2", "default", "squared value of force from the bias");
     317          24 :   keys.addOutputComponent("pressure", "default", "If using virial keyword, this is the current sum of virial terms. It is in units of pressure (energy / vol^3)");
     318          24 :   keys.addOutputComponent("_coupling", "default", "For each named CV biased, there will be a corresponding output CV_coupling storing the current linear bias prefactor.");
     319          12 : }
     320             : 
     321           8 : EDS::EDS(const ActionOptions &ao) : PLUMED_BIAS_INIT(ao),
     322           8 :   ncvs_(getNumberOfArguments()),
     323           8 :   scale_(ncvs_, 0.0),
     324           8 :   current_coupling_(ncvs_, 0.0),
     325           8 :   set_coupling_(ncvs_, 0.0),
     326           8 :   target_coupling_(ncvs_, 0.0),
     327           8 :   max_coupling_range_(ncvs_, 25.0),
     328           8 :   max_coupling_grad_(ncvs_, 0.0),
     329           8 :   coupling_rate_(ncvs_, 1.0),
     330           8 :   coupling_accum_(ncvs_, 0.0),
     331           8 :   means_(ncvs_, 0.0),
     332           8 :   step_size_(ncvs_, 0.0),
     333           8 :   pseudo_virial_(ncvs_),
     334           8 :   out_coupling_(ncvs_, NULL),
     335           8 :   in_restart_name_(""),
     336           8 :   out_restart_name_(""),
     337           8 :   fmt_("%f"),
     338           8 :   b_adaptive_(true),
     339           8 :   b_freeze_(false),
     340           8 :   b_equil_(true),
     341           8 :   b_ramp_(false),
     342           8 :   b_covar_(false),
     343           8 :   b_restart_(false),
     344           8 :   b_write_restart_(false),
     345           8 :   b_lm_(false),
     346           8 :   b_virial_(false),
     347           8 :   b_weights_(false),
     348           8 :   seed_(0),
     349           8 :   update_period_(0),
     350           8 :   avg_coupling_count_(1),
     351           8 :   update_calls_(0),
     352           8 :   kbt_(0.0),
     353           8 :   multi_prop_(-1.0),
     354           8 :   lm_mixing_par_(0.1),
     355           8 :   virial_scaling_(0.),
     356           8 :   pseudo_virial_sum_(0.0),
     357           8 :   max_logweight_(0.0),
     358           8 :   wsum_(0.0),
     359          32 :   value_force2_(NULL) {
     360           8 :   double temp = -1.0;
     361           8 :   bool b_mean = false;
     362             :   std::vector<Value *> wvalues;
     363             : 
     364           8 :   addComponent("force2");
     365           8 :   componentIsNotPeriodic("force2");
     366           8 :   value_force2_ = getPntrToComponent("force2");
     367             : 
     368          20 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     369          12 :     std::string comp = getPntrToArgument(i)->getName() + "_coupling";
     370          12 :     addComponent(comp);
     371          12 :     componentIsNotPeriodic(comp);
     372          12 :     out_coupling_[i] = getPntrToComponent(comp);
     373             :   }
     374             : 
     375           8 :   parseVector("CENTER", center_);
     376           8 :   parseArgumentList("CENTER_ARG", center_values_);
     377           8 :   parseArgumentList("LOGWEIGHTS", wvalues);
     378           8 :   parseVector("BIAS_SCALE", scale_);
     379           8 :   parseVector("RANGE", max_coupling_range_);
     380           8 :   parseVector("FIXED", target_coupling_);
     381           8 :   parseVector("INIT", set_coupling_);
     382           8 :   parse("PERIOD", update_period_);
     383           8 :   parse("TEMP", temp);
     384           8 :   parse("SEED", seed_);
     385           8 :   parse("MULTI_PROP", multi_prop_);
     386           8 :   parse("LM_MIXING", lm_mixing_par_);
     387           8 :   parse("RESTART_FMT", fmt_);
     388           8 :   parse("VIRIAL", virial_scaling_);
     389           8 :   fmt_ = " " + fmt_; // add space since parse strips them
     390           8 :   parse("OUT_RESTART", out_restart_name_);
     391           8 :   parseFlag("LM", b_lm_);
     392           8 :   parseFlag("RAMP", b_ramp_);
     393           8 :   parseFlag("FREEZE", b_freeze_);
     394           8 :   parseFlag("MEAN", b_mean);
     395           8 :   parseFlag("COVAR", b_covar_);
     396           8 :   parse("IN_RESTART", in_restart_name_);
     397           8 :   checkRead();
     398             : 
     399             :   /*
     400             :    * Things that are different when using changing centers:
     401             :    * 1. Scale
     402             :    * 2. The log file
     403             :    * 3. Reading Restarts
     404             :    */
     405             : 
     406           8 :   if (center_.size() == 0) {
     407           1 :     if (center_values_.size() == 0) {
     408           0 :       error("Must set either CENTER or CENTER_ARG");
     409           1 :     } else if (center_values_.size() != ncvs_) {
     410           0 :       error("CENTER_ARG must contain the same number of variables as ARG");
     411             :     }
     412           1 :     b_c_values_ = true;
     413           1 :     center_.resize(ncvs_);
     414           1 :     log.printf("  EDS will use possibly varying centers\n");
     415             :   } else {
     416           7 :     if (center_.size() != ncvs_) {
     417           0 :       error("Must have same number of CENTER arguments as ARG arguments");
     418           7 :     } else if (center_values_.size() != 0) {
     419           0 :       error("You can only set CENTER or CENTER_ARG. Not both");
     420             :     }
     421           7 :     b_c_values_ = false;
     422           7 :     log.printf("  EDS will use fixed centers\n");
     423             :   }
     424             : 
     425             :   // check for weights
     426           8 :   if (wvalues.size() > 1) {
     427           0 :     error("LOGWEIGHTS can only support one weight set. Please only pass one action");
     428           8 :   } else if (wvalues.size() == 1) {
     429           1 :     logweights_ = dynamic_cast<ReweightBase *>(wvalues[0]->getPntrToAction());
     430           1 :     b_weights_ = true;
     431             :   }
     432             : 
     433           8 :   log.printf("  setting scaling:");
     434           8 :   if (scale_.size() > 0 && scale_.size() < ncvs_) {
     435           0 :     error("the number of BIAS_SCALE values be the same as number of CVs");
     436           8 :   } else if (scale_.size() == 0 && b_c_values_) {
     437           0 :     log.printf(" Setting SCALE to be 1 for all CVs\n");
     438           0 :     scale_.resize(ncvs_);
     439           0 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     440           0 :       scale_[i] = 1;
     441             :     }
     442           8 :   } else if (scale_.size() == 0 && !b_c_values_) {
     443           2 :     log.printf(" (default) ");
     444             : 
     445           2 :     scale_.resize(ncvs_);
     446           6 :     for (unsigned int i = 0; i < scale_.size(); ++i) {
     447           4 :       if (center_[i] == 0) {
     448           0 :         error("BIAS_SCALE parameter has been set to CENTER value of 0 (as is default). This will divide by 0, so giving up. See doc for EDS bias");
     449             :       }
     450           4 :       scale_[i] = center_[i];
     451             :     }
     452             :   } else {
     453          14 :     for (unsigned int i = 0; i < scale_.size(); ++i) {
     454           8 :       log.printf(" %f", scale_[i]);
     455             :     }
     456             :   }
     457           8 :   log.printf("\n");
     458             : 
     459           8 :   if (b_lm_) {
     460           1 :     log.printf("  EDS will perform Levenberg-Marquardt minimization with mixing parameter = %f\n", lm_mixing_par_);
     461           1 :     differences_.resize(ncvs_);
     462           1 :     alpha_vector_.resize(ncvs_);
     463           1 :     alpha_vector_2_.resize(ncvs_);
     464           1 :     covar_.resize(ncvs_, ncvs_);
     465           1 :     covar2_.resize(ncvs_, ncvs_);
     466           1 :     lm_inv_.resize(ncvs_, ncvs_);
     467           1 :     covar2_ *= 0;
     468           1 :     lm_inv_ *= 0;
     469           1 :     if (multi_prop_ != 1) {
     470           0 :       log.printf("     WARNING - doing LM minimization but MULTI_PROP!=1\n");
     471             :     }
     472           7 :   } else if (b_covar_) {
     473           1 :     log.printf("  EDS will utilize covariance matrix for update steps\n");
     474           1 :     covar_.resize(ncvs_, ncvs_);
     475             :   } else {
     476           6 :     log.printf("  EDS will utilize variance for update steps\n");
     477           6 :     ssds_.resize(ncvs_);
     478             :   }
     479             : 
     480           8 :   b_virial_ = virial_scaling_;
     481             : 
     482           8 :   if (b_virial_) {
     483           1 :     if (ncvs_ == 1) {
     484           0 :       error("Minimizing the virial is only valid with multiply correlated collective variables.");
     485             :     }
     486             :     // check that the CVs can be used to compute pseudo-virial
     487           1 :     log.printf("  EDS will compute virials of CVs and penalize with scale of %f. Checking CVs are valid...", virial_scaling_);
     488           4 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     489           3 :       auto a = dynamic_cast<ActionAtomistic *>(getPntrToArgument(i)->getPntrToAction());
     490           3 :       if (!a) {
     491           0 :         error("If using VIRIAL keyword, you must have normal CVs as arguments to EDS. Offending action: " + getPntrToArgument(i)->getPntrToAction()->getName());
     492             :       }
     493             :       // cppcheck-suppress nullPointerRedundantCheck
     494           3 :       if (!(a->getPbc().isOrthorombic())) {
     495           3 :         log.printf("  WARNING: EDS Virial should have a orthorombic cell\n");
     496             :       }
     497             :     }
     498           1 :     log.printf("done.\n");
     499           1 :     addComponent("pressure");
     500           1 :     componentIsNotPeriodic("pressure");
     501           1 :     value_pressure_ = getPntrToComponent("pressure");
     502             :   }
     503             : 
     504           8 :   if (b_mean && !b_freeze_) {
     505           0 :     error("EDS keyword MEAN can only be used along with keyword FREEZE");
     506             :   }
     507             : 
     508           8 :   if (in_restart_name_ != "") {
     509           2 :     b_restart_ = true;
     510           2 :     log.printf("  reading simulation information from file: %s\n", in_restart_name_.c_str());
     511           2 :     readInRestart(b_mean);
     512             :   } else {
     513             : 
     514           6 :     if (temp >= 0.0) {
     515           6 :       kbt_ = plumed.getAtoms().getKBoltzmann() * temp;
     516             :     } else {
     517           0 :       kbt_ = plumed.getAtoms().getKbT();
     518             :     }
     519             : 
     520             :     // in driver, this results in kbt of 0
     521           6 :     if (kbt_ == 0) {
     522           0 :       error("  Unable to determine valid kBT. "
     523             :             "Could be because you are runnning from driver or MD didn't give temperature.\n"
     524             :             "Consider setting temperature manually with the TEMP keyword.");
     525             :       kbt_ = 1;
     526             :     }
     527             : 
     528           6 :     log.printf("  kBT = %f\n", kbt_);
     529           6 :     log.printf("  Updating every %i steps\n", update_period_);
     530             : 
     531           6 :     if (!b_c_values_) {
     532           5 :       log.printf("  with centers:");
     533          14 :       for (unsigned int i = 0; i < ncvs_; ++i) {
     534           9 :         log.printf(" %f ", center_[i]);
     535             :       }
     536             :     } else {
     537           1 :       log.printf("  with actions centers:");
     538           2 :       for (unsigned int i = 0; i < ncvs_; ++i) {
     539           1 :         log.printf(" %s ", center_values_[i]->getName().c_str());
     540             :         // add dependency on these actions
     541           1 :         addDependency(center_values_[i]->getPntrToAction());
     542             :       }
     543             :     }
     544             : 
     545           6 :     log.printf("\n  with initial ranges / rates:\n");
     546          16 :     for (unsigned int i = 0; i < max_coupling_range_.size(); ++i) {
     547             :       // this is just an empirical guess. Bigger range, bigger grads. Less frequent updates, bigger changes
     548             :       //
     549             :       // using the current maxing out scheme, max_coupling_range is the biggest step that can be taken in any given interval
     550          10 :       max_coupling_range_[i] *= kbt_;
     551          10 :       max_coupling_grad_[i] = max_coupling_range_[i];
     552          10 :       log.printf("    %f / %f\n", max_coupling_range_[i], max_coupling_grad_[i]);
     553             :     }
     554             : 
     555           6 :     if (seed_ > 0) {
     556           2 :       log.printf("  setting random seed = %i", seed_);
     557           2 :       rand_.setSeed(seed_);
     558             :     }
     559             : 
     560          16 :     for (unsigned int i = 0; i < ncvs_; ++i)
     561          10 :       if (target_coupling_[i] != 0.0) {
     562           1 :         b_adaptive_ = false;
     563             :       }
     564             : 
     565           6 :     if (!b_adaptive_) {
     566           1 :       if (b_ramp_) {
     567           1 :         log.printf("  ramping up coupling constants over %i steps\n", update_period_);
     568             :       }
     569             : 
     570           1 :       log.printf("  with starting coupling constants");
     571           2 :       for (unsigned int i = 0; i < set_coupling_.size(); ++i) {
     572           1 :         log.printf(" %f", set_coupling_[i]);
     573             :       }
     574           1 :       log.printf("\n");
     575           1 :       log.printf("  and final coupling constants");
     576           2 :       for (unsigned int i = 0; i < target_coupling_.size(); ++i) {
     577           1 :         log.printf(" %f", target_coupling_[i]);
     578             :       }
     579           1 :       log.printf("\n");
     580             :     }
     581             : 
     582             :     // now do setup
     583           6 :     if (b_ramp_) {
     584           1 :       update_period_ *= -1;
     585             :     }
     586             : 
     587          16 :     for (unsigned int i = 0; i < set_coupling_.size(); ++i) {
     588          10 :       current_coupling_[i] = set_coupling_[i];
     589             :     }
     590             : 
     591             :     // if b_adaptive_, then first half will be used for equilibrating and second half for statistics
     592           6 :     if (update_period_ > 0) {
     593           5 :       update_period_ /= 2;
     594             :     }
     595             :   }
     596             : 
     597           8 :   if (b_freeze_) {
     598           1 :     b_adaptive_ = false;
     599           1 :     update_period_ = 0;
     600           1 :     if (b_mean) {
     601           1 :       log.printf("  freezing bias at the average level from the restart file\n");
     602             :     } else {
     603           0 :       log.printf("  freezing bias at current level\n");
     604             :     }
     605             :   }
     606             : 
     607           8 :   if (multi_prop_ == -1.0) {
     608           5 :     log.printf("  Will update each dimension stochastically with probability 1 / number of CVs\n");
     609           5 :     multi_prop_ = 1.0 / ncvs_;
     610           3 :   } else if (multi_prop_ > 0 && multi_prop_ <= 1.0) {
     611           3 :     log.printf("  Will update each dimension stochastically with probability %f\n", multi_prop_);
     612             :   } else {
     613           0 :     error("  MULTI_PROP must be between 0 and 1\n");
     614             :   }
     615             : 
     616           8 :   if (out_restart_name_.length() > 0) {
     617           8 :     log.printf("  writing restart information every %i steps to file %s with format %s\n", abs(update_period_), out_restart_name_.c_str(), fmt_.c_str());
     618           8 :     b_write_restart_ = true;
     619           8 :     setupOutRestart();
     620             :   }
     621             : 
     622          16 :   log << "  Bibliography " << plumed.cite("White and Voth, J. Chem. Theory Comput. 10 (8), 3023-3030 (2014)") << "\n";
     623          16 :   log << "  Bibliography " << plumed.cite("G. M. Hocky, T. Dannenhoffer-Lafage, G. A. Voth, J. Chem. Theory Comput. 13 (9), 4593-4603 (2017)") << "\n";
     624           8 : }
     625             : 
     626           2 : void EDS::readInRestart(const bool b_mean) {
     627           2 :   int adaptive_i = 0;
     628             : 
     629           2 :   in_restart_.open(in_restart_name_);
     630             : 
     631           4 :   if (in_restart_.FieldExist("kbt")) {
     632           2 :     in_restart_.scanField("kbt", kbt_);
     633             :   } else {
     634           0 :     error("No field 'kbt' in restart file");
     635             :   }
     636           2 :   log.printf("  with kBT = %f\n", kbt_);
     637             : 
     638           4 :   if (in_restart_.FieldExist("update_period")) {
     639           2 :     in_restart_.scanField("update_period", update_period_);
     640             :   } else {
     641           0 :     error("No field 'update_period' in restart file");
     642             :   }
     643           2 :   log.printf("  Updating every %i steps\n", update_period_);
     644             : 
     645           4 :   if (in_restart_.FieldExist("adaptive")) {
     646             :     // note, no version of scanField for boolean
     647           2 :     in_restart_.scanField("adaptive", adaptive_i);
     648             :   } else {
     649           0 :     error("No field 'adaptive' in restart file");
     650             :   }
     651           2 :   b_adaptive_ = bool(adaptive_i);
     652             : 
     653           4 :   if (in_restart_.FieldExist("seed")) {
     654           2 :     in_restart_.scanField("seed", seed_);
     655             :   } else {
     656           0 :     error("No field 'seed' in restart file");
     657             :   }
     658           2 :   if (seed_ > 0) {
     659           0 :     log.printf("  setting random seed = %i", seed_);
     660           0 :     rand_.setSeed(seed_);
     661             :   }
     662             : 
     663             :   double time, tmp;
     664           2 :   std::vector<double> avg_bias = std::vector<double>(center_.size());
     665             :   unsigned int N = 0;
     666             :   std::string cv_name;
     667             : 
     668          24 :   while (in_restart_.scanField("time", time)) {
     669             : 
     670          20 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     671             :       cv_name = getPntrToArgument(i)->getName();
     672          10 :       in_restart_.scanField(cv_name + "_center", set_coupling_[i]);
     673          20 :       in_restart_.scanField(cv_name + "_set", set_coupling_[i]);
     674          20 :       in_restart_.scanField(cv_name + "_target", target_coupling_[i]);
     675          20 :       in_restart_.scanField(cv_name + "_coupling", current_coupling_[i]);
     676          20 :       in_restart_.scanField(cv_name + "_maxrange", max_coupling_range_[i]);
     677          20 :       in_restart_.scanField(cv_name + "_maxgrad", max_coupling_grad_[i]);
     678          20 :       in_restart_.scanField(cv_name + "_accum", coupling_accum_[i]);
     679          10 :       in_restart_.scanField(cv_name + "_mean", means_[i]);
     680          20 :       if (in_restart_.FieldExist(cv_name + "_pseudovirial")) {
     681           0 :         if (b_virial_) {
     682           0 :           in_restart_.scanField(cv_name + "_pseudovirial", pseudo_virial_[i]);
     683             :         } else { // discard the field
     684           0 :           in_restart_.scanField(cv_name + "_pseudovirial", tmp);
     685             :         }
     686             :       }
     687             :       // unused due to difference between covar/nocovar
     688          20 :       in_restart_.scanField(cv_name + "_std", tmp);
     689             : 
     690          10 :       avg_bias[i] += current_coupling_[i];
     691             :     }
     692          10 :     N++;
     693             : 
     694          10 :     in_restart_.scanField();
     695             :   }
     696             : 
     697           2 :   log.printf("  with centers:");
     698           4 :   for (unsigned int i = 0; i < center_.size(); ++i) {
     699           2 :     log.printf(" %f", center_[i]);
     700             :   }
     701           2 :   log.printf("\n  and scaling:");
     702           4 :   for (unsigned int i = 0; i < scale_.size(); ++i) {
     703           2 :     log.printf(" %f", scale_[i]);
     704             :   }
     705             : 
     706           2 :   log.printf("\n  with initial ranges / rates:\n");
     707           4 :   for (unsigned int i = 0; i < max_coupling_range_.size(); ++i) {
     708           2 :     log.printf("    %f / %f\n", max_coupling_range_[i], max_coupling_grad_[i]);
     709             :   }
     710             : 
     711           2 :   if (!b_adaptive_ && update_period_ < 0) {
     712           0 :     log.printf("  ramping up coupling constants over %i steps\n", -update_period_);
     713             :   }
     714             : 
     715           2 :   if (b_mean) {
     716           1 :     log.printf("Loaded in averages for coupling constants...\n");
     717           2 :     for (unsigned int i = 0; i < current_coupling_.size(); ++i) {
     718           1 :       current_coupling_[i] = avg_bias[i] / N;
     719             :     }
     720           2 :     for (unsigned int i = 0; i < current_coupling_.size(); ++i) {
     721           1 :       set_coupling_[i] = avg_bias[i] / N;
     722             :     }
     723             :   }
     724             : 
     725           2 :   log.printf("  with current coupling constants:\n    ");
     726           4 :   for (unsigned int i = 0; i < current_coupling_.size(); ++i) {
     727           2 :     log.printf(" %f", current_coupling_[i]);
     728             :   }
     729           2 :   log.printf("\n");
     730           2 :   log.printf("  with initial coupling constants:\n    ");
     731           4 :   for (unsigned int i = 0; i < set_coupling_.size(); ++i) {
     732           2 :     log.printf(" %f", set_coupling_[i]);
     733             :   }
     734           2 :   log.printf("\n");
     735           2 :   log.printf("  and final coupling constants:\n    ");
     736           4 :   for (unsigned int i = 0; i < target_coupling_.size(); ++i) {
     737           2 :     log.printf(" %f", target_coupling_[i]);
     738             :   }
     739           2 :   log.printf("\n");
     740             : 
     741           2 :   in_restart_.close();
     742           2 : }
     743             : 
     744           8 : void EDS::setupOutRestart() {
     745           8 :   out_restart_.link(*this);
     746           8 :   out_restart_.fmtField(fmt_);
     747           8 :   out_restart_.open(out_restart_name_);
     748             :   out_restart_.setHeavyFlush();
     749             : 
     750          16 :   out_restart_.addConstantField("adaptive").printField("adaptive", b_adaptive_);
     751          16 :   out_restart_.addConstantField("update_period").printField("update_period", update_period_);
     752          16 :   out_restart_.addConstantField("seed").printField("seed", seed_);
     753          16 :   out_restart_.addConstantField("kbt").printField("kbt", kbt_);
     754           8 : }
     755             : 
     756          27 : void EDS::writeOutRestart() {
     757             :   std::string cv_name;
     758          27 :   out_restart_.printField("time", getTimeStep() * getStep());
     759             : 
     760          66 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     761             :     cv_name = getPntrToArgument(i)->getName();
     762          78 :     out_restart_.printField(cv_name + "_center", center_[i]);
     763          78 :     out_restart_.printField(cv_name + "_set", set_coupling_[i]);
     764          78 :     out_restart_.printField(cv_name + "_target", target_coupling_[i]);
     765          78 :     out_restart_.printField(cv_name + "_coupling", current_coupling_[i]);
     766          78 :     out_restart_.printField(cv_name + "_maxrange", max_coupling_range_[i]);
     767          78 :     out_restart_.printField(cv_name + "_maxgrad", max_coupling_grad_[i]);
     768          78 :     out_restart_.printField(cv_name + "_accum", coupling_accum_[i]);
     769          39 :     out_restart_.printField(cv_name + "_mean", means_[i]);
     770          39 :     if (b_virial_) {
     771          18 :       out_restart_.printField(cv_name + "_pseudovirial", pseudo_virial_[i]);
     772             :     }
     773          39 :     if (!b_covar_ && !b_lm_) {
     774          42 :       out_restart_.printField(cv_name + "_std", ssds_[i] / (fmax(1, update_calls_ - 1)));
     775             :     } else {
     776          36 :       out_restart_.printField(cv_name + "_std", covar_(i, i) / (fmax(1, update_calls_ - 1)));
     777             :     }
     778             :   }
     779          27 :   out_restart_.printField();
     780          27 : }
     781             : 
     782          40 : void EDS::calculate() {
     783             : 
     784             :   // get center values from action if necessary
     785          40 :   if (b_c_values_)
     786          10 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     787           5 :       center_[i] = center_values_[i]->get();
     788             :     }
     789             : 
     790          40 :   apply_bias();
     791          40 : }
     792             : 
     793          40 : void EDS::apply_bias() {
     794             :   // Compute linear force as in "restraint"
     795             :   double ene = 0, totf2 = 0, cv, m, f;
     796             : 
     797         100 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     798          60 :     cv = difference(i, center_[i], getArgument(i));
     799          60 :     m = current_coupling_[i];
     800          60 :     f = -m;
     801          60 :     ene += m * cv;
     802             :     setOutputForce(i, f);
     803          60 :     totf2 += f * f;
     804             :   }
     805             : 
     806             :   setBias(ene);
     807          40 :   value_force2_->set(totf2);
     808          40 : }
     809             : 
     810          12 : void EDS::update_statistics() {
     811             :   double s, N, w = 1.0;
     812          12 :   std::vector<double> deltas(ncvs_);
     813             : 
     814             :   // update weight max, if necessary
     815          12 :   if (b_weights_) {
     816           2 :     w = logweights_->getLogWeight();
     817           2 :     if (max_logweight_ < w) {
     818             :       // we have new max. Need to shift existing values
     819           0 :       wsum_ *= exp(max_logweight_ - w);
     820           0 :       max_logweight_ = w;
     821             :     }
     822             :     // convert to weight
     823           2 :     w = exp(w - max_logweight_);
     824           2 :     wsum_ += w;
     825             :     N = wsum_;
     826             :   } else {
     827          10 :     N = fmax(1, update_calls_);
     828             :   }
     829             : 
     830             :   // Welford, West, and Hanso online variance method
     831             :   // with weights (default =  1.0)
     832          32 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     833          20 :     deltas[i] = difference(i, means_[i], getArgument(i)) * w;
     834          20 :     means_[i] += deltas[i] / N;
     835          20 :     if (!b_covar_ && !b_lm_) {
     836           8 :       ssds_[i] += deltas[i] * difference(i, means_[i], getArgument(i));
     837             :     }
     838             :   }
     839          12 :   if (b_covar_ || b_lm_) {
     840          16 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     841          36 :       for (unsigned int j = i; j < ncvs_; ++j) {
     842          24 :         s = (N - 1) * deltas[i] * deltas[j] / N / N - covar_(i, j) / N;
     843          24 :         covar_(i, j) += s;
     844             :         // do this so we don't double count
     845          24 :         covar_(j, i) = covar_(i, j);
     846             :       }
     847             :     }
     848             :   }
     849          12 :   if (b_virial_) {
     850           2 :     update_pseudo_virial();
     851             :   }
     852          12 : }
     853             : 
     854           8 : void EDS::reset_statistics() {
     855          20 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     856          12 :     means_[i] = 0;
     857          12 :     if (!b_covar_ && !b_lm_) {
     858           6 :       ssds_[i] = 0;
     859             :     }
     860             :   }
     861           8 :   if (b_covar_ || b_lm_)
     862           8 :     for (unsigned int i = 0; i < ncvs_; ++i)
     863          24 :       for (unsigned int j = 0; j < ncvs_; ++j) {
     864          18 :         covar_(i, j) = 0;
     865             :       }
     866           8 :   if (b_virial_) {
     867           4 :     for (unsigned int i = 0; i < ncvs_; ++i) {
     868           3 :       pseudo_virial_[i] = 0;
     869             :     }
     870           1 :     pseudo_virial_sum_ = 0;
     871             :   }
     872           8 :   if (b_weights_) {
     873           2 :     wsum_ = 0;
     874           2 :     max_logweight_ = 0;
     875             :   }
     876           8 : }
     877             : 
     878           1 : void EDS::calc_lm_step_size() {
     879             :   // calulcate step size
     880             :   // uses scale here, which by default is center
     881             : 
     882           1 :   mult(covar_, covar_, covar2_);
     883           4 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     884           3 :     differences_[i] = difference(i, center_[i], means_[i]);
     885           3 :     covar2_[i][i] += lm_mixing_par_ * covar2_[i][i];
     886             :   }
     887             : 
     888             :   // "step_size_vec" = 2*inv(covar*covar+ lambda diag(covar*covar))*covar*(mean-center)
     889           1 :   mult(covar_, differences_, alpha_vector_);
     890           1 :   Invert(covar2_, lm_inv_);
     891           1 :   mult(lm_inv_, alpha_vector_, alpha_vector_2_);
     892             : 
     893           4 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     894           3 :     step_size_[i] = 2 * alpha_vector_2_[i] / kbt_ / scale_[i];
     895             :   }
     896           1 : }
     897             : 
     898           1 : void EDS::calc_covar_step_size() {
     899             :   // calulcate step size
     900             :   // uses scale here, which by default is center
     901             :   double tmp;
     902           4 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     903             :     tmp = 0;
     904          12 :     for (unsigned int j = 0; j < ncvs_; ++j) {
     905           9 :       tmp += difference(i, center_[i], means_[i]) * covar_(i, j);
     906             :     }
     907           3 :     step_size_[i] = 2 * tmp / kbt_ / scale_[i] * update_calls_ / fmax(1, update_calls_ - 1);
     908             :   }
     909           1 : }
     910             : 
     911           6 : void EDS::calc_ssd_step_size() {
     912             :   double tmp;
     913          12 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     914           6 :     tmp = 2. * difference(i, center_[i], means_[i]) * ssds_[i] / fmax(1, update_calls_ - 1);
     915           6 :     step_size_[i] = tmp / kbt_ / scale_[i];
     916             :   }
     917           6 : }
     918             : 
     919           2 : void EDS::update_pseudo_virial() {
     920             :   // We want to compute the bias force on each atom times the position
     921             :   //  of the atoms.
     922             :   double p, netp = 0, netpv = 0;
     923             :   double volume = 0;
     924           8 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     925             :     // checked in setup to ensure this cast is valid.
     926           6 :     ActionAtomistic *cv = dynamic_cast<ActionAtomistic *>(getPntrToArgument(i)->getPntrToAction());
     927             :     Tensor &v(cv->modifyVirial());
     928           6 :     Tensor box(cv->getBox());
     929             :     const unsigned int natoms = cv->getNumberOfAtoms();
     930           6 :     if (!volume) {
     931           2 :       volume = box.determinant();
     932             :     }
     933             : 
     934             :     // pressure contribution is -dBias / dV
     935             :     // dBias / dV = alpha / w * dCV / dV
     936             :     // to get partial of CV wrt to volume
     937             :     // dCV/dV = sum dCV/dvij * vij / V
     938             :     // where vij is box element
     939             :     // add diagonal of virial tensor to get net pressure
     940             :     // TODO: replace this with adjugate (Jacobi's Formula)   for non-orthorombic case(?)
     941           6 :     p = v(0, 0) * box(0, 0) + v(1, 1) * box(1, 1) + v(2, 2) * box(2, 2);
     942           6 :     p /= volume;
     943             : 
     944           6 :     netp += p;
     945             : 
     946             :     // now scale for correct units in EDS algorithm
     947           6 :     p *= (volume) / (kbt_ * natoms);
     948             : 
     949             :     // compute running mean of scaled
     950           6 :     if (set_coupling_[i] != 0) {
     951           0 :       pseudo_virial_[i] = (p - pseudo_virial_[i]) / (fmax(1, update_calls_));
     952             :     } else {
     953           6 :       pseudo_virial_[i] = 0;
     954             :     }
     955             :     // update net pressure
     956           6 :     netpv += pseudo_virial_[i];
     957             :   }
     958             :   // update pressure
     959           2 :   value_pressure_->set(netp);
     960           2 :   pseudo_virial_sum_ = netpv;
     961           2 : }
     962             : 
     963           8 : void EDS::update_bias() {
     964           8 :   log.flush();
     965           8 :   if (b_lm_) {
     966           1 :     calc_lm_step_size();
     967           7 :   } else if (b_covar_) {
     968           1 :     calc_covar_step_size();
     969             :   } else {
     970           6 :     calc_ssd_step_size();
     971             :   }
     972             : 
     973          20 :   for (unsigned int i = 0; i < ncvs_; ++i) {
     974             : 
     975             :     // multidimesional stochastic step
     976          12 :     if (ncvs_ == 1 || (rand_.RandU01() < (multi_prop_))) {
     977             : 
     978          12 :       if (b_virial_) {
     979             :         // apply virial regularization
     980             :         //  P * dP/dcoupling
     981             :         //  coupling is already included in virial term due to plumed propogating from bias to forces
     982             :         //  thus we need to divide by it to get the derivative (since force is linear in coupling)
     983           3 :         if (fabs(set_coupling_[i]) > 0.000000001) // my heuristic for if EDS has started to prevent / 0
     984             :           // scale^2 here is to align units
     985             :         {
     986           0 :           step_size_[i] -= 2 * scale_[i] * scale_[i] * virial_scaling_ * pseudo_virial_sum_ * pseudo_virial_sum_ / set_coupling_[i];
     987             :         }
     988             :       }
     989          12 :       if (step_size_[i] == 0) {
     990           4 :         continue;
     991             :       }
     992             : 
     993             :       // clip gradient
     994           8 :       step_size_[i] = copysign(fmin(fabs(step_size_[i]), max_coupling_grad_[i]), step_size_[i]);
     995           8 :       coupling_accum_[i] += step_size_[i] * step_size_[i];
     996             : 
     997             :       // equation 5 in White and Voth, JCTC 2014
     998             :       // no negative sign because it's in step_size
     999           8 :       set_coupling_[i] += step_size_[i] * max_coupling_range_[i] / sqrt(coupling_accum_[i]);
    1000           8 :       coupling_rate_[i] = (set_coupling_[i] - current_coupling_[i]) / update_period_;
    1001             :     } else {
    1002             :       // do not change the bias
    1003           0 :       coupling_rate_[i] = 0;
    1004             :     }
    1005             :   }
    1006             : 
    1007             :   // reset means/vars
    1008           8 :   reset_statistics();
    1009           8 : }
    1010             : 
    1011          40 : void EDS::update() {
    1012             :   // adjust parameters according to EDS recipe
    1013          40 :   update_calls_++;
    1014             : 
    1015             :   // if we aren't wating for the bias to equilibrate, set flag to collect data
    1016             :   // want statistics before writing restart
    1017          40 :   if (!b_equil_ && update_period_ > 0) {
    1018          12 :     update_statistics();
    1019             :   }
    1020             : 
    1021             :   // write restart with correct statistics before bias update
    1022             :   // check if we're ramping or doing normal updates and then restart if needed. The ramping check
    1023             :   // is complicated because we could be frozen, finished ramping or not ramping.
    1024             :   // The + 2 is so we have an extra line showing that the bias isn't changing (for my sanity and yours)
    1025          40 :   if (b_write_restart_) {
    1026          40 :     if (getStep() == 0 ||
    1027          32 :         ((update_period_ < 0 && !b_freeze_ && update_calls_ <= fabs(update_period_) + 2) ||
    1028          24 :          (update_period_ > 0 && update_calls_ % update_period_ == 0))) {
    1029          27 :       writeOutRestart();
    1030             :     }
    1031             :   }
    1032             : 
    1033             :   int b_finished_equil_flag = 1;
    1034             : 
    1035             :   // assume forces already applied and saved
    1036             :   // are we ramping to a constant value and not done equilibrating?
    1037          40 :   if (update_period_ < 0) {
    1038           5 :     if (update_calls_ <= fabs(update_period_) && !b_freeze_) {
    1039           4 :       for (unsigned int i = 0; i < ncvs_; ++i) {
    1040           2 :         current_coupling_[i] += (target_coupling_[i] - set_coupling_[i]) / fabs(update_period_);
    1041             :       }
    1042             :     }
    1043             :     // make sure we don't reset update calls
    1044             :     b_finished_equil_flag = 0;
    1045          35 :   } else if (update_period_ == 0) {
    1046             :     // do we have a no-update case?
    1047             :     // not updating
    1048             :     // pass
    1049          30 :   } else if (b_equil_) {
    1050             :     // equilibrating
    1051             :     // check if we've reached the setpoint
    1052          48 :     for (unsigned int i = 0; i < ncvs_; ++i) {
    1053          30 :       if (coupling_rate_[i] == 0 || pow(current_coupling_[i] - set_coupling_[i], 2) < pow(coupling_rate_[i], 2)) {
    1054          14 :         b_finished_equil_flag &= 1;
    1055             :       } else {
    1056          16 :         current_coupling_[i] += coupling_rate_[i];
    1057             :         b_finished_equil_flag = 0;
    1058             :       }
    1059             :     }
    1060             :   }
    1061             : 
    1062             :   // reduce all the flags
    1063          40 :   if (b_equil_ && b_finished_equil_flag) {
    1064          11 :     b_equil_ = false;
    1065          11 :     update_calls_ = 0;
    1066             :   }
    1067             : 
    1068             :   // Now we update coupling constant, if necessary
    1069          40 :   if (!b_equil_ && update_period_ > 0 && update_calls_ == update_period_ && !b_freeze_) {
    1070           8 :     update_bias();
    1071           8 :     update_calls_ = 0;
    1072           8 :     avg_coupling_count_++;
    1073           8 :     b_equil_ = true; // back to equilibration now
    1074             :   }                  // close update if
    1075             : 
    1076             :   // pass couplings out so they are accessible
    1077         100 :   for (unsigned int i = 0; i < ncvs_; ++i) {
    1078          60 :     out_coupling_[i]->set(current_coupling_[i]);
    1079             :   }
    1080          40 : }
    1081             : 
    1082          16 : EDS::~EDS() {
    1083           8 :   out_restart_.close();
    1084          24 : }
    1085             : 
    1086           0 : void EDS::turnOnDerivatives() {
    1087             :   // do nothing
    1088             :   // this is to avoid errors triggered when a bias is used as a CV
    1089             :   // (This is done in ExtendedLagrangian.cpp)
    1090           0 : }
    1091             : 
    1092             : }
    1093             : } // close the 2 namespaces

Generated by: LCOV version 1.16