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