Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2011-2023 The plumed team
3 : (see the PEOPLE file at the root of the distribution for a list of names)
4 :
5 : See http://www.plumed.org for more information.
6 :
7 : This file is part of plumed, version 2.
8 :
9 : plumed is free software: you can redistribute it and/or modify
10 : it under the terms of the GNU Lesser General Public License as published by
11 : the Free Software Foundation, either version 3 of the License, or
12 : (at your option) any later version.
13 :
14 : plumed is distributed in the hope that it will be useful,
15 : but WITHOUT ANY WARRANTY; without even the implied warranty of
16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 : GNU Lesser General Public License for more details.
18 :
19 : You should have received a copy of the GNU Lesser General Public License
20 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 : #include "Bias.h"
23 : #include "core/ActionRegister.h"
24 : #include "core/ActionSet.h"
25 : #include "core/PlumedMain.h"
26 : #include "core/FlexibleBin.h"
27 : #include "tools/Exception.h"
28 : #include "tools/Grid.h"
29 : #include "tools/Matrix.h"
30 : #include "tools/OpenMP.h"
31 : #include "tools/Random.h"
32 : #include "tools/File.h"
33 : #include "tools/Communicator.h"
34 : #include <ctime>
35 : #include <numeric>
36 :
37 : namespace PLMD {
38 : namespace bias {
39 :
40 : //+PLUMEDOC BIAS METAD
41 : /*
42 : Used to performed metadynamics on one or more collective variables.
43 :
44 : In a metadynamics simulations a history dependent bias composed of
45 : intermittently added Gaussian functions is added to the potential \cite metad.
46 :
47 : \f[
48 : V(\vec{s},t) = \sum_{ k \tau < t} W(k \tau)
49 : \exp\left(
50 : -\sum_{i=1}^{d} \frac{(s_i-s_i^{(0)}(k \tau))^2}{2\sigma_i^2}
51 : \right).
52 : \f]
53 :
54 : This potential forces the system away from the kinetic traps in the potential energy surface
55 : and out into the unexplored parts of the energy landscape. Information on the Gaussian
56 : functions from which this potential is composed is output to a file called HILLS, which
57 : is used both the restart the calculation and to reconstruct the free energy as a function of the CVs.
58 : The free energy can be reconstructed from a metadynamics calculation because the final bias is given
59 : by:
60 :
61 : \f[
62 : V(\vec{s}) = -F(\vec{s})
63 : \f]
64 :
65 : During post processing the free energy can be calculated in this way using the \ref sum_hills
66 : utility.
67 :
68 : In the simplest possible implementation of a metadynamics calculation the expense of a metadynamics
69 : calculation increases with the length of the simulation as one has to, at every step, evaluate
70 : the values of a larger and larger number of Gaussian kernels. To avoid this issue you can
71 : store the bias on a grid. This approach is similar to that proposed in \cite babi08jcp but has the
72 : advantage that the grid spacing is independent on the Gaussian width.
73 : Notice that you should provide the grid boundaries (GRID_MIN and GRID_MAX) and either the number of bins
74 : for every collective variable (GRID_BIN) or the desired grid spacing (GRID_SPACING).
75 : In case you provide both PLUMED will use the most conservative choice (highest number of bins) for each dimension.
76 : In case you do not provide any information about bin size (neither GRID_BIN nor GRID_SPACING)
77 : PLUMED will use 1/5 of the Gaussian width (SIGMA) as grid spacing if the width is fixed or 1/5 of the minimum
78 : Gaussian width (SIGMA_MIN) if the width is variable. This default choice should be reasonable for most applications.
79 :
80 : Alternatively to the use of grids, it is possible to use a neighbor list to decrease the cost of evaluating the bias,
81 : this can be enabled using NLIST. NLIST can be beneficial with more than 2 collective variables, where GRID becomes
82 : expensive and memory consuming. The neighbor list will be updated everytime the CVs go farther than a cut-off value
83 : from the position they were at last neighbor list update. Gaussians are added to the neigbhor list if their center
84 : is within 6.*DP2CUTOFF*sigma*sigma. While the list is updated if the CVs are farther from the center than 0.5 of the
85 : standard deviation of the Gaussian center distribution of the list. These parameters (6 and 0.5) can be modified using
86 : NLIST_PARAMETERS. Note that the use of neighbor list does not provide the exact bias.
87 :
88 : Metadynamics can be restarted either from a HILLS file as well as from a GRID, in this second
89 : case one can first save a GRID using GRID_WFILE (and GRID_WSTRIDE) and at a later stage read
90 : it using GRID_RFILE.
91 :
92 : The work performed by the METAD bias can be calculated using CALC_WORK, note that this is expensive when not using grids.
93 :
94 : Another option that is available in plumed is well-tempered metadynamics \cite Barducci:2008. In this
95 : variant of metadynamics the heights of the Gaussian hills are scaled at each step so the bias is now
96 : given by:
97 :
98 : \f[
99 : V({s},t)= \sum_{t'=0,\tau_G,2\tau_G,\dots}^{t'<t} W e^{-V({s}({q}(t'),t'))/\Delta T} \exp\left(
100 : -\sum_{i=1}^{d} \frac{(s_i({q})-s_i({q}(t')))^2}{2\sigma_i^2}
101 : \right),
102 : \f]
103 :
104 : This method ensures that the bias converges more smoothly. It should be noted that, in the case of well-tempered metadynamics, in
105 : the output printed the Gaussian height is re-scaled using the bias factor.
106 : Also notice that with well-tempered metadynamics the HILLS file does not contain the bias,
107 : but the negative of the free-energy estimate. This choice has the advantage that
108 : one can restart a simulation using a different value for the \f$\Delta T\f$. The applied bias will be scaled accordingly.
109 :
110 : Note that you can use here also the flexible Gaussian approach \cite Branduardi:2012dl
111 : in which you can adapt the Gaussian to the extent of Cartesian space covered by a variable or
112 : to the space in collective variable covered in a given time. In this case the width of the deposited
113 : Gaussian potential is denoted by one value only that is a Cartesian space (ADAPTIVE=GEOM) or a time
114 : (ADAPTIVE=DIFF). Note that a specific integration technique for the deposited Gaussian kernels
115 : should be used in this case. Check the documentation for utility sum_hills.
116 :
117 : With the keyword INTERVAL one changes the metadynamics algorithm setting the bias force equal to zero
118 : outside boundary \cite baftizadeh2012protein. If, for example, metadynamics is performed on a CV s and one is interested only
119 : to the free energy for s > boundary, the history dependent potential is still updated according to the above
120 : equations but the metadynamics force is set to zero for s < boundary. Notice that Gaussian kernels are added also
121 : if s < boundary, as the tails of these Gaussian kernels influence VG in the relevant region s > boundary. In this way, the
122 : force on the system in the region s > boundary comes from both metadynamics and the force field, in the region
123 : s < boundary only from the latter. This approach allows obtaining a history-dependent bias potential VG that
124 : fluctuates around a stable estimator, equal to the negative of the free energy far enough from the
125 : boundaries. Note that:
126 : - It works only for one-dimensional biases;
127 : - It works both with and without GRID;
128 : - The interval limit boundary in a region where the free energy derivative is not large;
129 : - If in the region outside the limit boundary the system has a free energy minimum, the INTERVAL keyword should
130 : be used together with a \ref UPPER_WALLS or \ref LOWER_WALLS at boundary.
131 :
132 : As a final note, since version 2.0.2 when the system is outside of the selected interval the force
133 : is set to zero and the bias value to the value at the corresponding boundary. This allows acceptances
134 : for replica exchange methods to be computed correctly.
135 :
136 : Multiple walkers \cite multiplewalkers can also be used. See below the examples.
137 :
138 :
139 : The \f$c(t)\f$ reweighting factor can also be calculated on the fly using the equations
140 : presented in \cite Tiwary_jp504920s.
141 : The expression used to calculate \f$c(t)\f$ follows directly from Eq. 3 in \cite Tiwary_jp504920s,
142 : where \f$F(\vec{s})=-\gamma/(\gamma-1) V(\vec{s})\f$.
143 : This gives smoother results than equivalent Eqs. 13 and Eqs. 14 in that paper.
144 : The \f$c(t)\f$ is given by the rct component while the bias
145 : normalized by \f$c(t)\f$ is given by the rbias component (rbias=bias-rct) which can be used
146 : to obtain a reweighted histogram.
147 : The calculation of \f$c(t)\f$ is enabled by using the keyword CALC_RCT.
148 : By default \f$c(t)\f$ is updated every time the bias changes, but if this slows down the simulation
149 : the keyword RCT_USTRIDE can be set to a value higher than 1.
150 : This option requires that a grid is used.
151 :
152 : Additional material and examples can be also found in the tutorials:
153 :
154 : - \ref lugano-3
155 :
156 : Concurrent metadynamics
157 : as done e.g. in Ref. \cite gil2015enhanced . This indeed can be obtained by using the METAD
158 : action multiple times in the same input file.
159 :
160 : \par Examples
161 :
162 : The following input is for a standard metadynamics calculation using as
163 : collective variables the distance between atoms 3 and 5
164 : and the distance between atoms 2 and 4. The value of the CVs and
165 : the metadynamics bias potential are written to the COLVAR file every 100 steps.
166 : \plumedfile
167 : DISTANCE ATOMS=3,5 LABEL=d1
168 : DISTANCE ATOMS=2,4 LABEL=d2
169 : METAD ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3 PACE=500 LABEL=restraint
170 : PRINT ARG=d1,d2,restraint.bias STRIDE=100 FILE=COLVAR
171 : \endplumedfile
172 : (See also \ref DISTANCE \ref PRINT).
173 :
174 : \par
175 : If you use adaptive Gaussian kernels, with diffusion scheme where you use
176 : a Gaussian that should cover the space of 20 time steps in collective variables.
177 : Note that in this case the histogram correction is needed when summing up hills.
178 : \plumedfile
179 : DISTANCE ATOMS=3,5 LABEL=d1
180 : DISTANCE ATOMS=2,4 LABEL=d2
181 : METAD ARG=d1,d2 SIGMA=20 HEIGHT=0.3 PACE=500 LABEL=restraint ADAPTIVE=DIFF
182 : PRINT ARG=d1,d2,restraint.bias STRIDE=100 FILE=COLVAR
183 : \endplumedfile
184 :
185 : \par
186 : If you use adaptive Gaussian kernels, with geometrical scheme where you use
187 : a Gaussian that should cover the space of 0.05 nm in Cartesian space.
188 : Note that in this case the histogram correction is needed when summing up hills.
189 : \plumedfile
190 : DISTANCE ATOMS=3,5 LABEL=d1
191 : DISTANCE ATOMS=2,4 LABEL=d2
192 : METAD ARG=d1,d2 SIGMA=0.05 HEIGHT=0.3 PACE=500 LABEL=restraint ADAPTIVE=GEOM
193 : PRINT ARG=d1,d2,restraint.bias STRIDE=100 FILE=COLVAR
194 : \endplumedfile
195 :
196 : \par
197 : When using adaptive Gaussian kernels you might want to limit how the hills width can change.
198 : You can use SIGMA_MIN and SIGMA_MAX keywords.
199 : The sigmas should specified in terms of CV so you should use the CV units.
200 : Note that if you use a negative number, this means that the limit is not set.
201 : Note also that in this case the histogram correction is needed when summing up hills.
202 : \plumedfile
203 : DISTANCE ATOMS=3,5 LABEL=d1
204 : DISTANCE ATOMS=2,4 LABEL=d2
205 : METAD ...
206 : ARG=d1,d2 SIGMA=0.05 HEIGHT=0.3 PACE=500 LABEL=restraint ADAPTIVE=GEOM
207 : SIGMA_MIN=0.2,0.1 SIGMA_MAX=0.5,1.0
208 : ... METAD
209 : PRINT ARG=d1,d2,restraint.bias STRIDE=100 FILE=COLVAR
210 : \endplumedfile
211 :
212 : \par
213 : Multiple walkers can be also use as in \cite multiplewalkers
214 : These are enabled by setting the number of walker used, the id of the
215 : current walker which interprets the input file, the directory where the
216 : hills containing files resides, and the frequency to read the other walkers.
217 : Here is an example
218 : \plumedfile
219 : DISTANCE ATOMS=3,5 LABEL=d1
220 : METAD ...
221 : ARG=d1 SIGMA=0.05 HEIGHT=0.3 PACE=500 LABEL=restraint
222 : WALKERS_N=10
223 : WALKERS_ID=3
224 : WALKERS_DIR=../
225 : WALKERS_RSTRIDE=100
226 : ... METAD
227 : \endplumedfile
228 : where WALKERS_N is the total number of walkers, WALKERS_ID is the
229 : id of the present walker (starting from 0 ) and the WALKERS_DIR is the directory
230 : where all the walkers are located. WALKERS_RSTRIDE is the number of step between
231 : one update and the other. Since version 2.2.5, hills files are automatically
232 : flushed every WALKERS_RSTRIDE steps.
233 :
234 : \par
235 : The \f$c(t)\f$ reweighting factor can be calculated on the fly using the equations
236 : presented in \cite Tiwary_jp504920s as described above.
237 : This is enabled by using the keyword CALC_RCT,
238 : and can be done only if the bias is defined on a grid.
239 : \plumedfile
240 : phi: TORSION ATOMS=1,2,3,4
241 : psi: TORSION ATOMS=5,6,7,8
242 :
243 : METAD ...
244 : LABEL=metad
245 : ARG=phi,psi SIGMA=0.20,0.20 HEIGHT=1.20 BIASFACTOR=5 TEMP=300.0 PACE=500
246 : GRID_MIN=-pi,-pi GRID_MAX=pi,pi GRID_BIN=150,150
247 : CALC_RCT
248 : RCT_USTRIDE=10
249 : ... METAD
250 : \endplumedfile
251 : Here we have asked that the calculation is performed every 10 hills deposition by using
252 : RCT_USTRIDE keyword. If this keyword is not given, the calculation will
253 : by default be performed every time the bias changes. The \f$c(t)\f$ reweighting factor will be given
254 : in the rct component while the instantaneous value of the bias potential
255 : normalized using the \f$c(t)\f$ reweighting factor is given in the rbias component
256 : [rbias=bias-rct] which can be used to obtain a reweighted histogram or
257 : free energy surface using the \ref HISTOGRAM analysis.
258 :
259 : \par
260 : The kinetics of the transitions between basins can also be analyzed on the fly as
261 : in \cite PRL230602. The flag ACCELERATION turn on accumulation of the acceleration
262 : factor that can then be used to determine the rate. This method can be used together
263 : with \ref COMMITTOR analysis to stop the simulation when the system get to the target basin.
264 : It must be used together with Well-Tempered Metadynamics. If restarting from a previous
265 : metadynamics you need to use the ACCELERATION_RFILE keyword to give the name of the
266 : data file from which the previous value of the acceleration factor should be read, otherwise the
267 : calculation of the acceleration factor will be wrong.
268 :
269 : \par
270 : By using the flag FREQUENCY_ADAPTIVE the frequency adaptive scheme introduced in \cite Wang-JCP-2018
271 : is turned on. The frequency for hill addition then changes dynamically based on the acceleration factor
272 : according to the following equation
273 : \f[
274 : \tau_{\mathrm{dep}}(t) =
275 : \min\left[
276 : \tau_0 \cdot
277 : \max\left[\frac{\alpha(t)}{\theta},1\right]
278 : ,\tau_{c}
279 : \right]
280 : \f]
281 : where \f$\tau_0\f$ is the initial hill addition frequency given by the PACE keyword,
282 : \f$\tau_{c}\f$ is the maximum allowed frequency given by the FA_MAX_PACE keyword,
283 : \f$\alpha(t)\f$ is the instantaneous acceleration factor at time \f$t\f$,
284 : and \f$\theta\f$ is a threshold value that acceleration factor has to reach before
285 : triggering a change in the hill addition frequency given by the FA_MIN_ACCELERATION keyword.
286 : The frequency for updating the hill addition frequency according to this equation is
287 : given by the FA_UPDATE_FREQUENCY keyword, by default it is the same as the value given
288 : in PACE. The hill hill addition frequency increase monotonously such that if the
289 : instantaneous acceleration factor is lower than in the previous updating step the
290 : previous \f$\tau_{\mathrm{dep}}\f$ is kept rather than updating it to a lower value.
291 : The instantaneous hill addition frequency \f$\tau_{\mathrm{dep}}(t)\f$ is outputted
292 : to pace component. Note that if restarting from a previous metadynamics run you need to
293 : use the ACCELERATION_RFILE keyword to read in the acceleration factors from the
294 : previous run, otherwise the hill addition frequency will start from the initial
295 : frequency.
296 :
297 :
298 : \par
299 : You can also provide a target distribution using the keyword TARGET
300 : \cite white2015designing
301 : \cite marinelli2015ensemble
302 : \cite gil2016empirical
303 : The TARGET should be a grid containing a free-energy (i.e. the -\f$k_B\f$T*log of the desired target distribution).
304 : Gaussian kernels will then be scaled by a factor
305 : \f[
306 : e^{\beta(\tilde{F}(s)-\tilde{F}_{max})}
307 : \f]
308 : Here \f$\tilde{F}(s)\f$ is the free energy defined on the grid and \f$\tilde{F}_{max}\f$ its maximum value.
309 : Notice that we here used the maximum value as in ref \cite gil2016empirical
310 : This choice allows to avoid exceedingly large Gaussian kernels to be added. However,
311 : it could make the Gaussian too small. You should always choose carefully the HEIGHT parameter
312 : in this case.
313 : The grid file should be similar to other PLUMED grid files in that it should contain
314 : both the target free-energy and its derivatives.
315 :
316 : Notice that if you wish your simulation to converge to the target free energy you should use
317 : the DAMPFACTOR command to provide a global tempering \cite dama2014well
318 : Alternatively, if you use a BIASFACTOR your simulation will converge to a free
319 : energy that is a linear combination of the target free energy and of the intrinsic free energy
320 : determined by the original force field.
321 :
322 : \plumedfile
323 : DISTANCE ATOMS=3,5 LABEL=d1
324 : METAD ...
325 : LABEL=t1
326 : ARG=d1 SIGMA=0.05 TAU=200 DAMPFACTOR=100 PACE=250
327 : GRID_MIN=1.14 GRID_MAX=1.32 GRID_BIN=6
328 : TARGET=dist.grid
329 : ... METAD
330 :
331 : PRINT ARG=d1,t1.bias STRIDE=100 FILE=COLVAR
332 : \endplumedfile
333 :
334 : The file dist.dat for this calculation would read:
335 :
336 : \auxfile{dist.grid}
337 : #! FIELDS d1 t1.target der_d1
338 : #! SET min_d1 1.14
339 : #! SET max_d1 1.32
340 : #! SET nbins_d1 6
341 : #! SET periodic_d1 false
342 : 1.1400 0.0031 0.1101
343 : 1.1700 0.0086 0.2842
344 : 1.2000 0.0222 0.6648
345 : 1.2300 0.0521 1.4068
346 : 1.2600 0.1120 2.6873
347 : 1.2900 0.2199 4.6183
348 : 1.3200 0.3948 7.1055
349 : \endauxfile
350 :
351 : Notice that BIASFACTOR can also be chosen as equal to 1. In this case one will perform
352 : unbiased sampling. Instead of using HEIGHT, one should provide the TAU parameter.
353 : \plumedfile
354 : d: DISTANCE ATOMS=3,5
355 : METAD ARG=d SIGMA=0.1 TAU=4.0 TEMP=300 PACE=100 BIASFACTOR=1.0
356 : \endplumedfile
357 : The HILLS file obtained will still work with `plumed sum_hills` so as to plot a free-energy.
358 : The case where this makes sense is probably that of RECT simulations.
359 :
360 : Regarding RECT simulations, you can also use the RECT keyword so as to avoid using multiple input files.
361 : For instance, a single input file will be
362 : \plumedfile
363 : d: DISTANCE ATOMS=3,5
364 : METAD ARG=d SIGMA=0.1 TAU=4.0 TEMP=300 PACE=100 RECT=1.0,1.5,2.0,3.0
365 : \endplumedfile
366 : The number of elements in the RECT array should be equal to the number of replicas.
367 :
368 : */
369 : //+ENDPLUMEDOC
370 :
371 : class MetaD : public Bias {
372 :
373 : private:
374 : struct Gaussian {
375 : bool multivariate; // this is required to discriminate the one dimensional case
376 : double height;
377 : std::vector<double> center;
378 : std::vector<double> sigma;
379 : std::vector<double> invsigma;
380 5325 : Gaussian(const bool m, const double h, const std::vector<double>& c, const std::vector<double>& s):
381 5325 : multivariate(m),height(h),center(c),sigma(s),invsigma(s) {
382 : // to avoid troubles from zero element in flexible hills
383 15512 : for(unsigned i=0; i<invsigma.size(); ++i) {
384 10187 : if(std::abs(invsigma[i])>1.e-20) {
385 10187 : invsigma[i]=1.0/invsigma[i] ;
386 : } else {
387 0 : invsigma[i]=0.0;
388 : }
389 : }
390 5325 : }
391 : };
392 8 : struct TemperingSpecs {
393 : bool is_active;
394 : std::string name_stem;
395 : std::string name;
396 : double biasf;
397 : double threshold;
398 : double alpha;
399 156 : inline TemperingSpecs(bool is_active, const std::string &name_stem, const std::string &name, double biasf, double threshold, double alpha) :
400 156 : is_active(is_active), name_stem(name_stem), name(name), biasf(biasf), threshold(threshold), alpha(alpha)
401 156 : {}
402 : };
403 : // general setup
404 : double kbt_;
405 : int stride_;
406 : bool calc_work_;
407 : // well-tempered MetaD
408 : bool welltemp_;
409 : double biasf_;
410 : // output files format
411 : std::string fmt_;
412 : // first step?
413 : bool isFirstStep_;
414 : // Gaussian starting parameters
415 : double height0_;
416 : std::vector<double> sigma0_;
417 : std::vector<double> sigma0min_;
418 : std::vector<double> sigma0max_;
419 : // Gaussians
420 : std::vector<Gaussian> hills_;
421 : std::unique_ptr<FlexibleBin> flexbin_;
422 : int adaptive_;
423 : OFile hillsOfile_;
424 : std::vector<std::unique_ptr<IFile>> ifiles_;
425 : std::vector<std::string> ifilesnames_;
426 : // Grids
427 : bool grid_;
428 : std::unique_ptr<GridBase> BiasGrid_;
429 : OFile gridfile_;
430 : bool storeOldGrids_;
431 : int wgridstride_;
432 : // multiple walkers
433 : int mw_n_;
434 : std::string mw_dir_;
435 : int mw_id_;
436 : int mw_rstride_;
437 : bool walkers_mpi_;
438 : unsigned mpi_nw_;
439 : // flying gaussians
440 : bool flying_;
441 : // kinetics from metadynamics
442 : bool acceleration_;
443 : double acc_;
444 : double acc_restart_mean_;
445 : // transition-tempering metadynamics
446 : bool calc_max_bias_;
447 : double max_bias_;
448 : bool calc_transition_bias_;
449 : double transition_bias_;
450 : std::vector<std::vector<double> > transitionwells_;
451 : static const size_t n_tempering_options_ = 1;
452 : static const std::string tempering_names_[1][2];
453 : double dampfactor_;
454 : struct TemperingSpecs tt_specs_;
455 : std::string targetfilename_;
456 : std::unique_ptr<GridBase> TargetGrid_;
457 : // frequency adaptive metadynamics
458 : int current_stride_;
459 : bool freq_adaptive_;
460 : int fa_update_frequency_;
461 : int fa_max_stride_;
462 : double fa_min_acceleration_;
463 : // intervals
464 : double uppI_;
465 : double lowI_;
466 : bool doInt_;
467 : // reweighting
468 : bool calc_rct_;
469 : double reweight_factor_;
470 : unsigned rct_ustride_;
471 : // work
472 : double work_;
473 : // neighbour list stuff
474 : bool nlist_;
475 : bool nlist_update_;
476 : unsigned nlist_steps_;
477 : std::array<double,2> nlist_param_;
478 : std::vector<Gaussian> nlist_hills_;
479 : std::vector<double> nlist_center_;
480 : std::vector<double> nlist_dev2_;
481 :
482 : double stretchA=1.0;
483 : double stretchB=0.0;
484 :
485 : bool noStretchWarningDone=false;
486 :
487 12 : void noStretchWarning() {
488 12 : if(!noStretchWarningDone) {
489 3 : log<<"\nWARNING: you are using a HILLS file with Gaussian kernels, PLUMED 2.8 uses stretched Gaussians by default\n";
490 : }
491 12 : noStretchWarningDone=true;
492 12 : }
493 :
494 : static void registerTemperingKeywords(const std::string &name_stem, const std::string &name, Keywords &keys);
495 : void readTemperingSpecs(TemperingSpecs &t_specs);
496 : void logTemperingSpecs(const TemperingSpecs &t_specs);
497 : void readGaussians(IFile*);
498 : void writeGaussian(const Gaussian&,OFile&);
499 : void addGaussian(const Gaussian&);
500 : double getHeight(const std::vector<double>&);
501 : void temperHeight(double &height, const TemperingSpecs &t_specs, const double tempering_bias);
502 : double getBias(const std::vector<double>&);
503 : double getBiasAndDerivatives(const std::vector<double>&, std::vector<double>&);
504 : double evaluateGaussian(const std::vector<double>&, const Gaussian&);
505 : double evaluateGaussianAndDerivatives(const std::vector<double>&, const Gaussian&,std::vector<double>&,std::vector<double>&);
506 : double getGaussianNormalization(const Gaussian&);
507 : std::vector<unsigned> getGaussianSupport(const Gaussian&);
508 : bool scanOneHill(IFile* ifile, std::vector<Value>& v, std::vector<double>& center, std::vector<double>& sigma, double& height, bool& multivariate);
509 : void computeReweightingFactor();
510 : double getTransitionBarrierBias();
511 : void updateFrequencyAdaptiveStride();
512 : void updateNlist();
513 :
514 : public:
515 : explicit MetaD(const ActionOptions&);
516 : void calculate() override;
517 : void update() override;
518 : static void registerKeywords(Keywords& keys);
519 : bool checkNeedsGradients()const override;
520 : };
521 :
522 : PLUMED_REGISTER_ACTION(MetaD,"METAD")
523 :
524 159 : void MetaD::registerKeywords(Keywords& keys) {
525 159 : Bias::registerKeywords(keys);
526 318 : keys.addOutputComponent("rbias","CALC_RCT","the instantaneous value of the bias normalized using the c(t) reweighting factor [rbias=bias-rct]."
527 : "This component can be used to obtain a reweighted histogram.");
528 318 : keys.addOutputComponent("rct","CALC_RCT","the reweighting factor c(t).");
529 318 : keys.addOutputComponent("work","CALC_WORK","accumulator for work");
530 318 : keys.addOutputComponent("acc","ACCELERATION","the metadynamics acceleration factor");
531 318 : keys.addOutputComponent("maxbias", "CALC_MAX_BIAS", "the maximum of the metadynamics V(s, t)");
532 318 : keys.addOutputComponent("transbias", "CALC_TRANSITION_BIAS", "the metadynamics transition bias V*(t)");
533 318 : keys.addOutputComponent("pace","FREQUENCY_ADAPTIVE","the hill addition frequency when employing frequency adaptive metadynamics");
534 318 : keys.addOutputComponent("nlker","NLIST","number of hills in the neighbor list");
535 318 : keys.addOutputComponent("nlsteps","NLIST","number of steps from last neighbor list update");
536 159 : keys.use("ARG");
537 318 : keys.add("compulsory","SIGMA","the widths of the Gaussian hills");
538 318 : keys.add("compulsory","PACE","the frequency for hill addition");
539 318 : keys.add("compulsory","FILE","HILLS","a file in which the list of added hills is stored");
540 318 : keys.add("optional","HEIGHT","the heights of the Gaussian hills. Compulsory unless TAU and either BIASFACTOR or DAMPFACTOR are given");
541 318 : keys.add("optional","FMT","specify format for HILLS files (useful for decrease the number of digits in regtests)");
542 318 : keys.add("optional","BIASFACTOR","use well tempered metadynamics and use this bias factor. Please note you must also specify temp");
543 318 : keys.addFlag("CALC_WORK",false,"calculate the total accumulated work done by the bias since last restart");
544 318 : keys.add("optional","RECT","list of bias factors for all the replicas");
545 318 : keys.add("optional","DAMPFACTOR","damp hills with exp(-max(V)/(kT*DAMPFACTOR)");
546 318 : for (size_t i = 0; i < n_tempering_options_; i++) {
547 159 : registerTemperingKeywords(tempering_names_[i][0], tempering_names_[i][1], keys);
548 : }
549 318 : keys.add("optional","TARGET","target to a predefined distribution");
550 318 : keys.add("optional","TEMP","the system temperature - this is only needed if you are doing well-tempered metadynamics");
551 318 : keys.add("optional","TAU","in well tempered metadynamics, sets height to (k_B Delta T*pace*timestep)/tau");
552 318 : keys.addFlag("CALC_RCT",false,"calculate the c(t) reweighting factor and use that to obtain the normalized bias [rbias=bias-rct]."
553 : "This method is not compatible with metadynamics not on a grid.");
554 318 : keys.add("optional","RCT_USTRIDE","the update stride for calculating the c(t) reweighting factor."
555 : "The default 1, so c(t) is updated every time the bias is updated.");
556 318 : keys.add("optional","GRID_MIN","the lower bounds for the grid");
557 318 : keys.add("optional","GRID_MAX","the upper bounds for the grid");
558 318 : keys.add("optional","GRID_BIN","the number of bins for the grid");
559 318 : keys.add("optional","GRID_SPACING","the approximate grid spacing (to be used as an alternative or together with GRID_BIN)");
560 318 : keys.addFlag("GRID_SPARSE",false,"use a sparse grid to store hills");
561 318 : keys.addFlag("GRID_NOSPLINE",false,"don't use spline interpolation with grids");
562 318 : keys.add("optional","GRID_WSTRIDE","write the grid to a file every N steps");
563 318 : keys.add("optional","GRID_WFILE","the file on which to write the grid");
564 318 : keys.add("optional","GRID_RFILE","a grid file from which the bias should be read at the initial step of the simulation");
565 318 : keys.addFlag("STORE_GRIDS",false,"store all the grid files the calculation generates. They will be deleted if this keyword is not present");
566 318 : keys.addFlag("NLIST",false,"Use neighbor list for kernels summation, faster but experimental");
567 318 : keys.add("optional", "NLIST_PARAMETERS","(default=6.,0.5) the two cutoff parameters for the Gaussians neighbor list");
568 318 : keys.add("optional","ADAPTIVE","use a geometric (=GEOM) or diffusion (=DIFF) based hills width scheme. Sigma is one number that has distance units or time step dimensions");
569 318 : keys.add("optional","SIGMA_MAX","the upper bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
570 318 : keys.add("optional","SIGMA_MIN","the lower bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
571 318 : keys.add("optional","WALKERS_ID", "walker id");
572 318 : keys.add("optional","WALKERS_N", "number of walkers");
573 318 : keys.add("optional","WALKERS_DIR", "shared directory with the hills files from all the walkers");
574 318 : keys.add("optional","WALKERS_RSTRIDE","stride for reading hills files");
575 318 : keys.addFlag("WALKERS_MPI",false,"Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR");
576 318 : keys.add("optional","INTERVAL","one dimensional lower and upper limits, outside the limits the system will not feel the biasing force.");
577 318 : keys.addFlag("FLYING_GAUSSIAN",false,"Switch on flying Gaussian method, must be used with WALKERS_MPI");
578 318 : keys.addFlag("ACCELERATION",false,"Set to TRUE if you want to compute the metadynamics acceleration factor.");
579 318 : keys.add("optional","ACCELERATION_RFILE","a data file from which the acceleration should be read at the initial step of the simulation");
580 318 : keys.addFlag("CALC_MAX_BIAS", false, "Set to TRUE if you want to compute the maximum of the metadynamics V(s, t)");
581 318 : keys.addFlag("CALC_TRANSITION_BIAS", false, "Set to TRUE if you want to compute a metadynamics transition bias V*(t)");
582 318 : keys.add("numbered", "TRANSITIONWELL", "This keyword appears multiple times as TRANSITIONWELL followed by an integer. Each specifies the coordinates for one well as in transition-tempered metadynamics. At least one must be provided.");
583 318 : keys.addFlag("FREQUENCY_ADAPTIVE",false,"Set to TRUE if you want to enable frequency adaptive metadynamics such that the frequency for hill addition to change dynamically based on the acceleration factor.");
584 318 : keys.add("optional","FA_UPDATE_FREQUENCY","the frequency for updating the hill addition pace in frequency adaptive metadynamics, by default this is equal to the value given in PACE");
585 318 : keys.add("optional","FA_MAX_PACE","the maximum hill addition frequency allowed in frequency adaptive metadynamics. By default there is no maximum value.");
586 318 : keys.add("optional","FA_MIN_ACCELERATION","only update the hill addition pace in frequency adaptive metadynamics after reaching the minimum acceleration factor given here. By default it is 1.0.");
587 159 : keys.use("RESTART");
588 159 : keys.use("UPDATE_FROM");
589 159 : keys.use("UPDATE_UNTIL");
590 159 : }
591 :
592 : const std::string MetaD::tempering_names_[1][2] = {{"TT", "transition tempered"}};
593 :
594 159 : void MetaD::registerTemperingKeywords(const std::string &name_stem, const std::string &name, Keywords &keys) {
595 318 : keys.add("optional", name_stem + "BIASFACTOR", "use " + name + " metadynamics with this bias factor. Please note you must also specify temp");
596 318 : keys.add("optional", name_stem + "BIASTHRESHOLD", "use " + name + " metadynamics with this bias threshold. Please note you must also specify " + name_stem + "BIASFACTOR");
597 318 : keys.add("optional", name_stem + "ALPHA", "use " + name + " metadynamics with this hill size decay exponent parameter. Please note you must also specify " + name_stem + "BIASFACTOR");
598 159 : }
599 :
600 157 : MetaD::MetaD(const ActionOptions& ao):
601 : PLUMED_BIAS_INIT(ao),
602 156 : kbt_(0.0),
603 156 : stride_(0),
604 156 : calc_work_(false),
605 156 : welltemp_(false),
606 156 : biasf_(-1.0),
607 156 : isFirstStep_(true),
608 156 : height0_(std::numeric_limits<double>::max()),
609 156 : adaptive_(FlexibleBin::none),
610 156 : grid_(false),
611 156 : wgridstride_(0),
612 156 : mw_n_(1), mw_dir_(""), mw_id_(0), mw_rstride_(1),
613 156 : walkers_mpi_(false), mpi_nw_(0),
614 156 : flying_(false),
615 156 : acceleration_(false), acc_(0.0), acc_restart_mean_(0.0),
616 156 : calc_max_bias_(false), max_bias_(0.0),
617 156 : calc_transition_bias_(false), transition_bias_(0.0),
618 156 : dampfactor_(0.0),
619 312 : tt_specs_(false, "TT", "Transition Tempered", -1.0, 0.0, 1.0),
620 156 : current_stride_(0),
621 156 : freq_adaptive_(false),
622 156 : fa_update_frequency_(0),
623 156 : fa_max_stride_(0),
624 156 : fa_min_acceleration_(1.0),
625 156 : uppI_(-1), lowI_(-1), doInt_(false),
626 156 : calc_rct_(false),
627 156 : reweight_factor_(0.0),
628 156 : rct_ustride_(1),
629 156 : work_(0),
630 156 : nlist_(false),
631 156 : nlist_update_(false),
632 469 : nlist_steps_(0) {
633 156 : if(!dp2cutoffNoStretch()) {
634 156 : stretchA=dp2cutoffA;
635 156 : stretchB=dp2cutoffB;
636 : }
637 : // parse the flexible hills
638 : std::string adaptiveoption;
639 : adaptiveoption="NONE";
640 312 : parse("ADAPTIVE",adaptiveoption);
641 156 : if(adaptiveoption=="GEOM") {
642 22 : log.printf(" Uses Geometry-based hills width: sigma must be in distance units and only one sigma is needed\n");
643 22 : adaptive_=FlexibleBin::geometry;
644 134 : } else if(adaptiveoption=="DIFF") {
645 3 : log.printf(" Uses Diffusion-based hills width: sigma must be in time steps and only one sigma is needed\n");
646 3 : adaptive_=FlexibleBin::diffusion;
647 131 : } else if(adaptiveoption=="NONE") {
648 130 : adaptive_=FlexibleBin::none;
649 : } else {
650 1 : error("I do not know this type of adaptive scheme");
651 : }
652 :
653 155 : parse("FMT",fmt_);
654 :
655 : // parse the sigma
656 155 : parseVector("SIGMA",sigma0_);
657 155 : if(adaptive_==FlexibleBin::none) {
658 : // if you use normal sigma you need one sigma per argument
659 130 : if( sigma0_.size()!=getNumberOfArguments() ) {
660 0 : error("number of arguments does not match number of SIGMA parameters");
661 : }
662 : } else {
663 : // if you use flexible hills you need one sigma
664 25 : if(sigma0_.size()!=1) {
665 1 : error("If you choose ADAPTIVE you need only one sigma according to your choice of type (GEOM/DIFF)");
666 : }
667 : // if adaptive then the number must be an integer
668 24 : if(adaptive_==FlexibleBin::diffusion) {
669 3 : if(int(sigma0_[0])-sigma0_[0]>1.e-9 || int(sigma0_[0])-sigma0_[0] <-1.e-9 || int(sigma0_[0])<1 ) {
670 0 : error("In case of adaptive hills with diffusion, the sigma must be an integer which is the number of time steps\n");
671 : }
672 : }
673 : // here evtl parse the sigma min and max values
674 48 : parseVector("SIGMA_MIN",sigma0min_);
675 24 : if(sigma0min_.size()>0 && sigma0min_.size()!=getNumberOfArguments()) {
676 1 : error("the number of SIGMA_MIN values be the same of the number of the arguments");
677 23 : } else if(sigma0min_.size()==0) {
678 23 : sigma0min_.resize(getNumberOfArguments());
679 67 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
680 44 : sigma0min_[i]=-1.;
681 : }
682 : }
683 :
684 46 : parseVector("SIGMA_MAX",sigma0max_);
685 23 : if(sigma0max_.size()>0 && sigma0max_.size()!=getNumberOfArguments()) {
686 1 : error("the number of SIGMA_MAX values be the same of the number of the arguments");
687 22 : } else if(sigma0max_.size()==0) {
688 22 : sigma0max_.resize(getNumberOfArguments());
689 64 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
690 42 : sigma0max_[i]=-1.;
691 : }
692 : }
693 :
694 44 : flexbin_=Tools::make_unique<FlexibleBin>(adaptive_,this,sigma0_[0],sigma0min_,sigma0max_);
695 : }
696 :
697 : // note: HEIGHT is not compulsory, since one could use the TAU keyword, see below
698 152 : parse("HEIGHT",height0_);
699 152 : parse("PACE",stride_);
700 151 : if(stride_<=0) {
701 0 : error("frequency for hill addition is nonsensical");
702 : }
703 151 : current_stride_ = stride_;
704 159 : std::string hillsfname="HILLS";
705 151 : parse("FILE",hillsfname);
706 :
707 : // Manually set to calculate special bias quantities
708 : // throughout the course of simulation. (These are chosen due to
709 : // relevance for tempering and event-driven logic as well.)
710 151 : parseFlag("CALC_MAX_BIAS", calc_max_bias_);
711 305 : parseFlag("CALC_TRANSITION_BIAS", calc_transition_bias_);
712 :
713 : std::vector<double> rect_biasf_;
714 302 : parseVector("RECT",rect_biasf_);
715 151 : if(rect_biasf_.size()>0) {
716 18 : int r=0;
717 18 : if(comm.Get_rank()==0) {
718 9 : r=multi_sim_comm.Get_rank();
719 : }
720 18 : comm.Bcast(r,0);
721 18 : biasf_=rect_biasf_[r];
722 18 : log<<" You are using RECT\n";
723 : } else {
724 266 : parse("BIASFACTOR",biasf_);
725 : }
726 151 : if( biasf_<1.0 && biasf_!=-1.0) {
727 0 : error("well tempered bias factor is nonsensical");
728 : }
729 151 : parse("DAMPFACTOR",dampfactor_);
730 151 : kbt_=getkBT();
731 151 : if(biasf_>=1.0) {
732 38 : if(kbt_==0.0) {
733 0 : error("Unless the MD engine passes the temperature to plumed, with well-tempered metad you must specify it using TEMP");
734 : }
735 38 : welltemp_=true;
736 : }
737 151 : if(dampfactor_>0.0) {
738 2 : if(kbt_==0.0) {
739 0 : error("Unless the MD engine passes the temperature to plumed, with damped metad you must specify it using TEMP");
740 : }
741 : }
742 :
743 151 : parseFlag("CALC_WORK",calc_work_);
744 :
745 : // Set transition tempering parameters.
746 : // Transition wells are read later via calc_transition_bias_.
747 151 : readTemperingSpecs(tt_specs_);
748 151 : if (tt_specs_.is_active) {
749 3 : calc_transition_bias_ = true;
750 : }
751 :
752 : // If any previous option specified to calculate a transition bias,
753 : // now read the transition wells for that quantity.
754 151 : if (calc_transition_bias_) {
755 13 : std::vector<double> tempcoords(getNumberOfArguments());
756 26 : for (unsigned i = 0; ; i++) {
757 78 : if (!parseNumberedVector("TRANSITIONWELL", i, tempcoords) ) {
758 : break;
759 : }
760 26 : if (tempcoords.size() != getNumberOfArguments()) {
761 0 : error("incorrect number of coordinates for transition tempering well");
762 : }
763 26 : transitionwells_.push_back(tempcoords);
764 : }
765 : }
766 :
767 302 : parse("TARGET",targetfilename_);
768 151 : if(targetfilename_.length()>0 && kbt_==0.0) {
769 0 : error("with TARGET temperature must be specified");
770 : }
771 151 : double tau=0.0;
772 151 : parse("TAU",tau);
773 151 : if(tau==0.0) {
774 129 : if(height0_==std::numeric_limits<double>::max()) {
775 0 : error("At least one between HEIGHT and TAU should be specified");
776 : }
777 : // if tau is not set, we compute it here from the other input parameters
778 129 : if(welltemp_) {
779 19 : tau=(kbt_*(biasf_-1.0))/height0_*getTimeStep()*stride_;
780 110 : } else if(dampfactor_>0.0) {
781 0 : tau=(kbt_*dampfactor_)/height0_*getTimeStep()*stride_;
782 : }
783 : } else {
784 22 : if(height0_!=std::numeric_limits<double>::max()) {
785 0 : error("At most one between HEIGHT and TAU should be specified");
786 : }
787 22 : if(welltemp_) {
788 19 : if(biasf_!=1.0) {
789 15 : height0_=(kbt_*(biasf_-1.0))/tau*getTimeStep()*stride_;
790 : } else {
791 4 : height0_=kbt_/tau*getTimeStep()*stride_; // special case for gamma=1
792 : }
793 3 : } else if(dampfactor_>0.0) {
794 2 : height0_=(kbt_*dampfactor_)/tau*getTimeStep()*stride_;
795 : } else {
796 1 : error("TAU only makes sense in well-tempered or damped metadynamics");
797 : }
798 : }
799 :
800 : // Grid Stuff
801 153 : std::vector<std::string> gmin(getNumberOfArguments());
802 300 : parseVector("GRID_MIN",gmin);
803 150 : if(gmin.size()!=getNumberOfArguments() && gmin.size()!=0) {
804 0 : error("not enough values for GRID_MIN");
805 : }
806 150 : std::vector<std::string> gmax(getNumberOfArguments());
807 300 : parseVector("GRID_MAX",gmax);
808 150 : if(gmax.size()!=getNumberOfArguments() && gmax.size()!=0) {
809 0 : error("not enough values for GRID_MAX");
810 : }
811 150 : std::vector<unsigned> gbin(getNumberOfArguments());
812 : std::vector<double> gspacing;
813 300 : parseVector("GRID_BIN",gbin);
814 150 : if(gbin.size()!=getNumberOfArguments() && gbin.size()!=0) {
815 0 : error("not enough values for GRID_BIN");
816 : }
817 300 : parseVector("GRID_SPACING",gspacing);
818 150 : if(gspacing.size()!=getNumberOfArguments() && gspacing.size()!=0) {
819 0 : error("not enough values for GRID_SPACING");
820 : }
821 150 : if(gmin.size()!=gmax.size()) {
822 0 : error("GRID_MAX and GRID_MIN should be either present or absent");
823 : }
824 150 : if(gspacing.size()!=0 && gmin.size()==0) {
825 0 : error("If GRID_SPACING is present also GRID_MIN and GRID_MAX should be present");
826 : }
827 150 : if(gbin.size()!=0 && gmin.size()==0) {
828 0 : error("If GRID_BIN is present also GRID_MIN and GRID_MAX should be present");
829 : }
830 150 : if(gmin.size()!=0) {
831 61 : if(gbin.size()==0 && gspacing.size()==0) {
832 6 : if(adaptive_==FlexibleBin::none) {
833 6 : log<<" Binsize not specified, 1/5 of sigma will be be used\n";
834 6 : plumed_assert(sigma0_.size()==getNumberOfArguments());
835 6 : gspacing.resize(getNumberOfArguments());
836 13 : for(unsigned i=0; i<gspacing.size(); i++) {
837 7 : gspacing[i]=0.2*sigma0_[i];
838 : }
839 : } else {
840 : // with adaptive hills and grid a sigma min must be specified
841 0 : for(unsigned i=0; i<sigma0min_.size(); i++)
842 0 : if(sigma0min_[i]<=0) {
843 0 : error("When using ADAPTIVE Gaussians on a grid SIGMA_MIN must be specified");
844 : }
845 0 : log<<" Binsize not specified, 1/5 of sigma_min will be be used\n";
846 0 : gspacing.resize(getNumberOfArguments());
847 0 : for(unsigned i=0; i<gspacing.size(); i++) {
848 0 : gspacing[i]=0.2*sigma0min_[i];
849 : }
850 : }
851 55 : } else if(gspacing.size()!=0 && gbin.size()==0) {
852 2 : log<<" The number of bins will be estimated from GRID_SPACING\n";
853 53 : } else if(gspacing.size()!=0 && gbin.size()!=0) {
854 1 : log<<" You specified both GRID_BIN and GRID_SPACING\n";
855 1 : log<<" The more conservative (highest) number of bins will be used for each variable\n";
856 : }
857 61 : if(gbin.size()==0) {
858 8 : gbin.assign(getNumberOfArguments(),1);
859 : }
860 61 : if(gspacing.size()!=0)
861 21 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
862 : double a,b;
863 13 : Tools::convert(gmin[i],a);
864 12 : Tools::convert(gmax[i],b);
865 12 : unsigned n=std::ceil(((b-a)/gspacing[i]));
866 12 : if(gbin[i]<n) {
867 11 : gbin[i]=n;
868 : }
869 : }
870 : }
871 149 : if(gbin.size()>0) {
872 60 : grid_=true;
873 : }
874 :
875 149 : bool sparsegrid=false;
876 149 : parseFlag("GRID_SPARSE",sparsegrid);
877 149 : bool nospline=false;
878 149 : parseFlag("GRID_NOSPLINE",nospline);
879 149 : bool spline=!nospline;
880 300 : parse("GRID_WSTRIDE",wgridstride_);
881 : std::string gridfilename_;
882 149 : parse("GRID_WFILE",gridfilename_);
883 149 : parseFlag("STORE_GRIDS",storeOldGrids_);
884 149 : if(grid_ && gridfilename_.length()>0) {
885 19 : if(wgridstride_==0 ) {
886 0 : error("frequency with which to output grid not specified use GRID_WSTRIDE");
887 : }
888 : }
889 149 : if(grid_ && wgridstride_>0) {
890 19 : if(gridfilename_.length()==0) {
891 1 : error("grid filename not specified use GRID_WFILE");
892 : }
893 : }
894 :
895 : std::string gridreadfilename_;
896 149 : parse("GRID_RFILE",gridreadfilename_);
897 :
898 149 : if(!grid_&&gridfilename_.length()> 0) {
899 0 : error("To write a grid you need first to define it!");
900 : }
901 149 : if(!grid_&&gridreadfilename_.length()>0) {
902 0 : error("To read a grid you need first to define it!");
903 : }
904 :
905 : /*setup neighbor list stuff*/
906 298 : parseFlag("NLIST", nlist_);
907 149 : nlist_center_.resize(getNumberOfArguments());
908 149 : nlist_dev2_.resize(getNumberOfArguments());
909 149 : if(nlist_&&grid_) {
910 1 : error("NLIST and GRID cannot be combined!");
911 : }
912 : std::vector<double> nlist_param;
913 298 : parseVector("NLIST_PARAMETERS",nlist_param);
914 149 : if(nlist_param.size()==0) {
915 149 : nlist_param_[0]=6.0;//*DP2CUTOFF -> max distance of neighbors
916 149 : nlist_param_[1]=0.5;//*nlist_dev2_[i] -> condition for rebuilding
917 : } else {
918 0 : plumed_massert(nlist_param.size()==2,"two cutoff parameters are needed for the neighbor list");
919 0 : plumed_massert(nlist_param[0]>1.0,"the first of NLIST_PARAMETERS must be greater than 1. The smaller the first, the smaller should be the second as well");
920 0 : const double min_PARAM_1=(1.-1./std::sqrt(nlist_param[0]/2))+0.16;
921 0 : plumed_massert(nlist_param[1]>0,"the second of NLIST_PARAMETERS must be greater than 0");
922 0 : plumed_massert(nlist_param[1]<=min_PARAM_1,"the second of NLIST_PARAMETERS must be smaller to avoid systematic errors. Largest suggested value is: 1.16-1/sqrt(PARAM_0/2) = "+std::to_string(min_PARAM_1));
923 0 : nlist_param_[0]=nlist_param[0];
924 0 : nlist_param_[1]=nlist_param[1];
925 : }
926 :
927 : // Reweighting factor rct
928 149 : parseFlag("CALC_RCT",calc_rct_);
929 149 : if (calc_rct_) {
930 6 : plumed_massert(grid_,"CALC_RCT is supported only if bias is on a grid");
931 : }
932 149 : parse("RCT_USTRIDE",rct_ustride_);
933 :
934 149 : if(dampfactor_>0.0) {
935 2 : if(!grid_) {
936 0 : error("With DAMPFACTOR you should use grids");
937 : }
938 : }
939 :
940 : // Multiple walkers
941 149 : parse("WALKERS_N",mw_n_);
942 149 : parse("WALKERS_ID",mw_id_);
943 149 : if(mw_n_<=mw_id_) {
944 0 : error("walker ID should be a numerical value less than the total number of walkers");
945 : }
946 149 : parse("WALKERS_DIR",mw_dir_);
947 149 : parse("WALKERS_RSTRIDE",mw_rstride_);
948 :
949 : // MPI version
950 149 : parseFlag("WALKERS_MPI",walkers_mpi_);
951 :
952 : //If this Action is not compiled with MPI the user is informed and we exit gracefully
953 149 : if(walkers_mpi_) {
954 39 : plumed_assert(Communicator::plumedHasMPI()) << "Invalid walkers configuration: WALKERS_MPI flag requires MPI compilation";
955 40 : plumed_assert(Communicator::initialized()) << "Invalid walkers configuration: WALKERS_MPI needs the communicator correctly initialized.";
956 : }
957 :
958 : // Flying Gaussian
959 148 : parseFlag("FLYING_GAUSSIAN", flying_);
960 :
961 : // Inteval keyword
962 149 : std::vector<double> tmpI(2);
963 296 : parseVector("INTERVAL",tmpI);
964 148 : if(tmpI.size()!=2&&tmpI.size()!=0) {
965 0 : error("both a lower and an upper limits must be provided with INTERVAL");
966 148 : } else if(tmpI.size()==2) {
967 2 : lowI_=tmpI.at(0);
968 2 : uppI_=tmpI.at(1);
969 2 : if(getNumberOfArguments()!=1) {
970 0 : error("INTERVAL limits correction works only for monodimensional metadynamics!");
971 : }
972 2 : if(uppI_<lowI_) {
973 0 : error("The Upper limit must be greater than the Lower limit!");
974 : }
975 2 : if(getPntrToArgument(0)->isPeriodic()) {
976 0 : error("INTERVAL cannot be used with periodic variables!");
977 : }
978 2 : doInt_=true;
979 : }
980 :
981 296 : parseFlag("ACCELERATION",acceleration_);
982 : // Check for a restart acceleration if acceleration is active.
983 : std::string acc_rfilename;
984 148 : if(acceleration_) {
985 8 : parse("ACCELERATION_RFILE", acc_rfilename);
986 : }
987 :
988 148 : freq_adaptive_=false;
989 148 : parseFlag("FREQUENCY_ADAPTIVE",freq_adaptive_);
990 : //
991 148 : fa_update_frequency_=0;
992 148 : parse("FA_UPDATE_FREQUENCY",fa_update_frequency_);
993 148 : if(fa_update_frequency_!=0 && !freq_adaptive_) {
994 0 : plumed_merror("It doesn't make sense to use the FA_MAX_PACE keyword if frequency adaptive METAD hasn't been activated by using the FREQUENCY_ADAPTIVE flag");
995 : }
996 148 : if(fa_update_frequency_==0 && freq_adaptive_) {
997 0 : fa_update_frequency_=stride_;
998 : }
999 : //
1000 148 : fa_max_stride_=0;
1001 148 : parse("FA_MAX_PACE",fa_max_stride_);
1002 148 : if(fa_max_stride_!=0 && !freq_adaptive_) {
1003 0 : plumed_merror("It doesn't make sense to use the FA_MAX_PACE keyword if frequency adaptive METAD hasn't been activated by using the FREQUENCY_ADAPTIVE flag");
1004 : }
1005 : //
1006 148 : fa_min_acceleration_=1.0;
1007 148 : parse("FA_MIN_ACCELERATION",fa_min_acceleration_);
1008 148 : if(fa_min_acceleration_!=1.0 && !freq_adaptive_) {
1009 0 : plumed_merror("It doesn't make sense to use the FA_MIN_ACCELERATION keyword if frequency adaptive METAD hasn't been activated by using the FREQUENCY_ADAPTIVE flag");
1010 : }
1011 :
1012 148 : checkRead();
1013 :
1014 148 : log.printf(" Gaussian width ");
1015 148 : if (adaptive_==FlexibleBin::diffusion) {
1016 2 : log.printf(" (Note: The units of sigma are in timesteps) ");
1017 : }
1018 148 : if (adaptive_==FlexibleBin::geometry) {
1019 19 : log.printf(" (Note: The units of sigma are in dist units) ");
1020 : }
1021 396 : for(unsigned i=0; i<sigma0_.size(); ++i) {
1022 248 : log.printf(" %f",sigma0_[i]);
1023 : }
1024 148 : log.printf(" Gaussian height %f\n",height0_);
1025 148 : log.printf(" Gaussian deposition pace %d\n",stride_);
1026 148 : log.printf(" Gaussian file %s\n",hillsfname.c_str());
1027 148 : if(welltemp_) {
1028 38 : log.printf(" Well-Tempered Bias Factor %f\n",biasf_);
1029 38 : log.printf(" Hills relaxation time (tau) %f\n",tau);
1030 38 : log.printf(" KbT %f\n",kbt_);
1031 : }
1032 :
1033 : // Transition tempered metadynamics options
1034 148 : if (tt_specs_.is_active) {
1035 3 : logTemperingSpecs(tt_specs_);
1036 : // Check that the appropriate transition bias quantity is calculated.
1037 : // (Should never trip, given that the flag is automatically set.)
1038 3 : if (!calc_transition_bias_) {
1039 0 : error(" transition tempering requires calculation of a transition bias");
1040 : }
1041 : }
1042 :
1043 : // Overall tempering sanity check (this gets tricky when multiple are active).
1044 : // When multiple temperings are active, it's fine to have one tempering attempt
1045 : // to increase hill size with increasing bias, so long as the others can shrink
1046 : // the hills faster than it increases their size in the long-time limit.
1047 : // This set of checks ensures that the hill sizes eventually decay to zero as c(t)
1048 : // diverges to infinity.
1049 : // The alpha parameter allows hills to decay as 1/t^alpha instead of 1/t,
1050 : // a slower decay, so as t -> infinity, only the temperings with the largest
1051 : // alphas govern the final asymptotic decay. (Alpha helps prevent false convergence.)
1052 148 : if (welltemp_ || dampfactor_ > 0.0 || tt_specs_.is_active) {
1053 : // Determine the number of active temperings.
1054 : int n_active = 0;
1055 43 : if (welltemp_) {
1056 : n_active++;
1057 : }
1058 43 : if (dampfactor_ > 0.0) {
1059 2 : n_active++;
1060 : }
1061 43 : if (tt_specs_.is_active) {
1062 3 : n_active++;
1063 : }
1064 : // Find the greatest alpha.
1065 43 : double greatest_alpha = 0.0;
1066 43 : if (welltemp_) {
1067 38 : greatest_alpha = std::max(greatest_alpha, 1.0);
1068 : }
1069 43 : if (dampfactor_ > 0.0) {
1070 4 : greatest_alpha = std::max(greatest_alpha, 1.0);
1071 : }
1072 43 : if (tt_specs_.is_active) {
1073 6 : greatest_alpha = std::max(greatest_alpha, tt_specs_.alpha);
1074 : }
1075 : // Find the least alpha.
1076 43 : double least_alpha = 1.0;
1077 : if (welltemp_) {
1078 : least_alpha = std::min(least_alpha, 1.0);
1079 : }
1080 43 : if (dampfactor_ > 0.0) {
1081 2 : least_alpha = std::min(least_alpha, 1.0);
1082 : }
1083 43 : if (tt_specs_.is_active) {
1084 4 : least_alpha = std::min(least_alpha, tt_specs_.alpha);
1085 : }
1086 : // Find the inverse harmonic average of the delta T parameters for all
1087 : // of the temperings with the greatest alpha values.
1088 : double total_governing_deltaT_inv = 0.0;
1089 43 : if (welltemp_ && 1.0 == greatest_alpha && biasf_ != 1.0) {
1090 34 : total_governing_deltaT_inv += 1.0 / (biasf_ - 1.0);
1091 : }
1092 43 : if (dampfactor_ > 0.0 && 1.0 == greatest_alpha) {
1093 2 : total_governing_deltaT_inv += 1.0 / (dampfactor_);
1094 : }
1095 43 : if (tt_specs_.is_active && tt_specs_.alpha == greatest_alpha) {
1096 3 : total_governing_deltaT_inv += 1.0 / (tt_specs_.biasf - 1.0);
1097 : }
1098 : // Give a newbie-friendly error message for people using one tempering if
1099 : // only one is active.
1100 43 : if (n_active == 1 && total_governing_deltaT_inv < 0.0) {
1101 0 : error("for stable tempering, the bias factor must be greater than one");
1102 : // Give a slightly more complex error message to users stacking multiple
1103 : // tempering options at a time, but all with uniform alpha values.
1104 43 : } else if (total_governing_deltaT_inv < 0.0 && greatest_alpha == least_alpha) {
1105 0 : error("for stable tempering, the sum of the inverse Delta T parameters must be greater than zero!");
1106 : // Give the most technical error message to users stacking multiple tempering
1107 : // options with different alpha parameters.
1108 43 : } else if (total_governing_deltaT_inv < 0.0 && greatest_alpha != least_alpha) {
1109 0 : error("for stable tempering, the sum of the inverse Delta T parameters for the greatest asymptotic hill decay exponents must be greater than zero!");
1110 : }
1111 : }
1112 :
1113 148 : if(doInt_) {
1114 2 : log.printf(" Upper and Lower limits boundaries for the bias are activated at %f - %f\n", lowI_, uppI_);
1115 : }
1116 :
1117 148 : if(grid_) {
1118 60 : log.printf(" Grid min");
1119 161 : for(unsigned i=0; i<gmin.size(); ++i) {
1120 101 : log.printf(" %s",gmin[i].c_str() );
1121 : }
1122 60 : log.printf("\n");
1123 60 : log.printf(" Grid max");
1124 161 : for(unsigned i=0; i<gmax.size(); ++i) {
1125 101 : log.printf(" %s",gmax[i].c_str() );
1126 : }
1127 60 : log.printf("\n");
1128 60 : log.printf(" Grid bin");
1129 161 : for(unsigned i=0; i<gbin.size(); ++i) {
1130 101 : log.printf(" %u",gbin[i]);
1131 : }
1132 60 : log.printf("\n");
1133 60 : if(spline) {
1134 60 : log.printf(" Grid uses spline interpolation\n");
1135 : }
1136 60 : if(sparsegrid) {
1137 6 : log.printf(" Grid uses sparse grid\n");
1138 : }
1139 60 : if(wgridstride_>0) {
1140 19 : log.printf(" Grid is written on file %s with stride %d\n",gridfilename_.c_str(),wgridstride_);
1141 : }
1142 : }
1143 :
1144 148 : if(mw_n_>1) {
1145 6 : if(walkers_mpi_) {
1146 0 : error("MPI version of multiple walkers is not compatible with filesystem version of multiple walkers");
1147 : }
1148 6 : log.printf(" %d multiple walkers active\n",mw_n_);
1149 6 : log.printf(" walker id %d\n",mw_id_);
1150 6 : log.printf(" reading stride %d\n",mw_rstride_);
1151 6 : if(mw_dir_!="") {
1152 3 : log.printf(" directory with hills files %s\n",mw_dir_.c_str());
1153 : }
1154 : } else {
1155 142 : if(walkers_mpi_) {
1156 38 : log.printf(" Multiple walkers active using MPI communnication\n");
1157 38 : if(mw_dir_!="") {
1158 0 : log.printf(" directory with hills files %s\n",mw_dir_.c_str());
1159 : }
1160 38 : if(comm.Get_rank()==0) {
1161 : // Only root of group can communicate with other walkers
1162 23 : mpi_nw_=multi_sim_comm.Get_size();
1163 : }
1164 : // Communicate to the other members of the same group
1165 : // info abount number of walkers and walker index
1166 38 : comm.Bcast(mpi_nw_,0);
1167 : }
1168 : }
1169 :
1170 148 : if(flying_) {
1171 6 : if(!walkers_mpi_) {
1172 0 : error("Flying Gaussian method must be used with MPI version of multiple walkers");
1173 : }
1174 6 : log.printf(" Flying Gaussian method with %d walkers active\n",mpi_nw_);
1175 : }
1176 :
1177 148 : if(nlist_) {
1178 2 : addComponent("nlker");
1179 2 : componentIsNotPeriodic("nlker");
1180 2 : addComponent("nlsteps");
1181 2 : componentIsNotPeriodic("nlsteps");
1182 : }
1183 :
1184 148 : if(calc_rct_) {
1185 12 : addComponent("rbias");
1186 12 : componentIsNotPeriodic("rbias");
1187 12 : addComponent("rct");
1188 6 : componentIsNotPeriodic("rct");
1189 6 : log.printf(" The c(t) reweighting factor will be calculated every %u hills\n",rct_ustride_);
1190 12 : getPntrToComponent("rct")->set(reweight_factor_);
1191 : }
1192 :
1193 148 : if(calc_work_) {
1194 2 : addComponent("work");
1195 2 : componentIsNotPeriodic("work");
1196 : }
1197 :
1198 148 : if(acceleration_) {
1199 4 : if (kbt_ == 0.0) {
1200 0 : error("The calculation of the acceleration works only if simulation temperature has been defined");
1201 : }
1202 4 : log.printf(" calculation on the fly of the acceleration factor\n");
1203 8 : addComponent("acc");
1204 8 : componentIsNotPeriodic("acc");
1205 : // Set the initial value of the the acceleration.
1206 : // If this is not a restart, set to 1.0.
1207 4 : if (acc_rfilename.length() == 0) {
1208 4 : getPntrToComponent("acc")->set(1.0);
1209 2 : if(getRestart()) {
1210 1 : log.printf(" WARNING: calculating the acceleration factor in a restarted run without reading in the previous value will most likely lead to incorrect results.\n");
1211 1 : log.printf(" You should use the ACCELERATION_RFILE keyword.\n");
1212 : }
1213 : // Otherwise, read and set the restart value.
1214 : } else {
1215 : // Restart of acceleration does not make sense if the restart timestep is zero.
1216 : //if (getStep() == 0) {
1217 : // error("Restarting calculation of acceleration factors works only if simulation timestep is restarted correctly");
1218 : //}
1219 : // Open the ACCELERATION_RFILE.
1220 2 : IFile acc_rfile;
1221 2 : acc_rfile.link(*this);
1222 2 : if(acc_rfile.FileExist(acc_rfilename)) {
1223 2 : acc_rfile.open(acc_rfilename);
1224 : } else {
1225 0 : error("The ACCELERATION_RFILE file you want to read: " + acc_rfilename + ", cannot be found!");
1226 : }
1227 : // Read the file to find the restart acceleration.
1228 2 : double acc_rmean=0.0;
1229 2 : double acc_rtime=0.0;
1230 : bool found=false;
1231 2 : std::string acclabel = getLabel() + ".acc";
1232 2 : acc_rfile.allowIgnoredFields();
1233 248 : while(acc_rfile.scanField("time", acc_rtime)) {
1234 122 : acc_rfile.scanField(acclabel, acc_rmean);
1235 122 : acc_rfile.scanField();
1236 : found=true;
1237 : }
1238 2 : if(!found) {
1239 0 : error("The ACCELERATION_RFILE file you want to read: " + acc_rfilename + ", does not contain a time field!");
1240 : }
1241 2 : acc_restart_mean_ = acc_rmean;
1242 : // Set component based on the read values.
1243 2 : getPntrToComponent("acc")->set(acc_rmean);
1244 2 : log.printf(" initial acceleration factor read from file %s: value of %f at time %f\n",acc_rfilename.c_str(),acc_rmean,acc_rtime);
1245 2 : }
1246 : }
1247 :
1248 148 : if (calc_max_bias_) {
1249 0 : if (!grid_) {
1250 0 : error("Calculating the maximum bias on the fly works only with a grid");
1251 : }
1252 0 : log.printf(" calculation on the fly of the maximum bias max(V(s,t)) \n");
1253 0 : addComponent("maxbias");
1254 0 : componentIsNotPeriodic("maxbias");
1255 : }
1256 :
1257 148 : if (calc_transition_bias_) {
1258 13 : if (!grid_) {
1259 0 : error("Calculating the transition bias on the fly works only with a grid");
1260 : }
1261 13 : log.printf(" calculation on the fly of the transition bias V*(t)\n");
1262 26 : addComponent("transbias");
1263 13 : componentIsNotPeriodic("transbias");
1264 13 : log<<" Number of transition wells "<<transitionwells_.size()<<"\n";
1265 13 : if (transitionwells_.size() == 0) {
1266 0 : error("Calculating the transition bias on the fly requires definition of at least one transition well");
1267 : }
1268 : // Check that a grid is in use.
1269 13 : if (!grid_) {
1270 0 : error(" transition barrier finding requires a grid for the bias");
1271 : }
1272 : // Log the wells and check that they are in the grid.
1273 39 : for (unsigned i = 0; i < transitionwells_.size(); i++) {
1274 : // Log the coordinate.
1275 26 : log.printf(" Transition well %d at coordinate ", i);
1276 64 : for (unsigned j = 0; j < getNumberOfArguments(); j++) {
1277 38 : log.printf("%f ", transitionwells_[i][j]);
1278 : }
1279 26 : log.printf("\n");
1280 : // Check that the coordinate is in the grid.
1281 64 : for (unsigned j = 0; j < getNumberOfArguments(); j++) {
1282 : double max, min;
1283 38 : Tools::convert(gmin[j], min);
1284 38 : Tools::convert(gmax[j], max);
1285 38 : if (transitionwells_[i][j] < min || transitionwells_[i][j] > max) {
1286 0 : error(" transition well is not in grid");
1287 : }
1288 : }
1289 : }
1290 : }
1291 :
1292 148 : if(freq_adaptive_) {
1293 2 : if(!acceleration_) {
1294 0 : plumed_merror("Frequency adaptive metadynamics only works if the calculation of the acceleration factor is enabled with the ACCELERATION keyword\n");
1295 : }
1296 2 : if(walkers_mpi_) {
1297 0 : plumed_merror("Combining frequency adaptive metadynamics with MPI multiple walkers is not allowed");
1298 : }
1299 :
1300 2 : log.printf(" Frequency adaptive metadynamics enabled\n");
1301 2 : if(getRestart() && acc_rfilename.length() == 0) {
1302 0 : log.printf(" WARNING: using the frequency adaptive scheme in a restarted run without reading in the previous value of the acceleration factor will most likely lead to incorrect results.\n");
1303 0 : log.printf(" You should use the ACCELERATION_RFILE keyword.\n");
1304 : }
1305 2 : log.printf(" The frequency for hill addition will change dynamically based on the metadynamics acceleration factor\n");
1306 2 : log.printf(" The hill addition frequency will be updated every %d steps\n",fa_update_frequency_);
1307 2 : if(fa_min_acceleration_>1.0) {
1308 2 : log.printf(" The hill addition frequency will only be updated once the metadynamics acceleration factor becomes larger than %.1f \n",fa_min_acceleration_);
1309 : }
1310 2 : if(fa_max_stride_!=0) {
1311 2 : log.printf(" The hill addition frequency will not become larger than %d steps\n",fa_max_stride_);
1312 : }
1313 4 : addComponent("pace");
1314 2 : componentIsNotPeriodic("pace");
1315 2 : updateFrequencyAdaptiveStride();
1316 : }
1317 :
1318 : // initializing and checking grid
1319 : bool restartedFromGrid=false; // restart from grid file
1320 148 : if(grid_) {
1321 60 : if(!(gridreadfilename_.length()>0)) {
1322 : // check for mesh and sigma size
1323 116 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
1324 : double a,b;
1325 74 : Tools::convert(gmin[i],a);
1326 74 : Tools::convert(gmax[i],b);
1327 74 : double mesh=(b-a)/((double)gbin[i]);
1328 74 : if(adaptive_==FlexibleBin::none) {
1329 74 : if(mesh>0.5*sigma0_[i]) {
1330 38 : log<<" WARNING: Using a METAD with a Grid Spacing larger than half of the Gaussians width (SIGMA) can produce artifacts\n";
1331 : }
1332 : } else {
1333 0 : if(sigma0min_[i]<0.) {
1334 0 : error("When using ADAPTIVE Gaussians on a grid SIGMA_MIN must be specified");
1335 : }
1336 0 : if(mesh>0.5*sigma0min_[i]) {
1337 0 : log<<" WARNING: to use a METAD with a GRID and ADAPTIVE you need to set a Grid Spacing lower than half of the Gaussians (SIGMA_MIN) \n";
1338 : }
1339 : }
1340 : }
1341 42 : std::string funcl=getLabel() + ".bias";
1342 42 : if(!sparsegrid) {
1343 36 : BiasGrid_=Tools::make_unique<Grid>(funcl,getArguments(),gmin,gmax,gbin,spline,true);
1344 : } else {
1345 6 : BiasGrid_=Tools::make_unique<SparseGrid>(funcl,getArguments(),gmin,gmax,gbin,spline,true);
1346 : }
1347 42 : std::vector<std::string> actualmin=BiasGrid_->getMin();
1348 42 : std::vector<std::string> actualmax=BiasGrid_->getMax();
1349 116 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
1350 : std::string is;
1351 74 : Tools::convert(i,is);
1352 74 : if(gmin[i]!=actualmin[i]) {
1353 0 : error("GRID_MIN["+is+"] must be adjusted to "+actualmin[i]+" to fit periodicity");
1354 : }
1355 74 : if(gmax[i]!=actualmax[i]) {
1356 0 : error("GRID_MAX["+is+"] must be adjusted to "+actualmax[i]+" to fit periodicity");
1357 : }
1358 : }
1359 42 : } else {
1360 : // read the grid in input, find the keys
1361 18 : if(walkers_mpi_&&gridreadfilename_.at(0)!='/') {
1362 : //if possible the root replica will share its current folder so that all walkers will read the same file
1363 0 : const std::string ret = std::filesystem::current_path();
1364 0 : gridreadfilename_ = "/" + gridreadfilename_;
1365 0 : gridreadfilename_ = ret + gridreadfilename_;
1366 0 : if(comm.Get_rank()==0) {
1367 0 : multi_sim_comm.Bcast(gridreadfilename_,0);
1368 : }
1369 0 : comm.Bcast(gridreadfilename_,0);
1370 : }
1371 18 : IFile gridfile;
1372 18 : gridfile.link(*this);
1373 18 : if(gridfile.FileExist(gridreadfilename_)) {
1374 18 : gridfile.open(gridreadfilename_);
1375 : } else {
1376 0 : error("The GRID file you want to read: " + gridreadfilename_ + ", cannot be found!");
1377 : }
1378 18 : std::string funcl=getLabel() + ".bias";
1379 36 : BiasGrid_=GridBase::create(funcl, getArguments(), gridfile, gmin, gmax, gbin, sparsegrid, spline, true);
1380 18 : if(BiasGrid_->getDimension()!=getNumberOfArguments()) {
1381 0 : error("mismatch between dimensionality of input grid and number of arguments");
1382 : }
1383 45 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1384 54 : if( getPntrToArgument(i)->isPeriodic()!=BiasGrid_->getIsPeriodic()[i] ) {
1385 0 : error("periodicity mismatch between arguments and input bias");
1386 : }
1387 : double a, b;
1388 27 : Tools::convert(gmin[i],a);
1389 27 : Tools::convert(gmax[i],b);
1390 27 : double mesh=(b-a)/((double)gbin[i]);
1391 27 : if(mesh>0.5*sigma0_[i]) {
1392 27 : log<<" WARNING: Using a METAD with a Grid Spacing larger than half of the Gaussians width can produce artifacts\n";
1393 : }
1394 : }
1395 18 : log.printf(" Restarting from %s\n",gridreadfilename_.c_str());
1396 18 : if(getRestart()) {
1397 : restartedFromGrid=true;
1398 : }
1399 18 : }
1400 : }
1401 :
1402 : // if we are restarting from GRID and using WALKERS_MPI we can check that all walkers have actually read the grid
1403 148 : if(getRestart()&&walkers_mpi_) {
1404 9 : std::vector<int> restarted(mpi_nw_,0);
1405 9 : if(comm.Get_rank()==0) {
1406 6 : multi_sim_comm.Allgather(int(restartedFromGrid), restarted);
1407 : }
1408 9 : comm.Bcast(restarted,0);
1409 : int result = std::accumulate(restarted.begin(),restarted.end(),0);
1410 9 : if(result!=0&&result!=mpi_nw_) {
1411 0 : error("in this WALKERS_MPI run some replica have restarted from GRID while other do not!");
1412 : }
1413 : }
1414 :
1415 186 : if(walkers_mpi_&&mw_dir_==""&&hillsfname.at(0)!='/') {
1416 : //if possible the root replica will share its current folder so that all walkers will read the same file
1417 76 : const std::string ret = std::filesystem::current_path();
1418 : mw_dir_ = ret;
1419 38 : mw_dir_ = mw_dir_ + "/";
1420 38 : if(comm.Get_rank()==0) {
1421 23 : multi_sim_comm.Bcast(mw_dir_,0);
1422 : }
1423 38 : comm.Bcast(mw_dir_,0);
1424 : }
1425 :
1426 : // creating std::vector of ifile* for hills reading
1427 : // open all files at the beginning and read Gaussians if restarting
1428 : bool restartedFromHills=false; // restart from hills files
1429 308 : for(int i=0; i<mw_n_; ++i) {
1430 : std::string fname;
1431 160 : if(mw_dir_!="") {
1432 47 : if(mw_n_>1) {
1433 9 : std::stringstream out;
1434 9 : out << i;
1435 18 : fname = mw_dir_+"/"+hillsfname+"."+out.str();
1436 47 : } else if(walkers_mpi_) {
1437 76 : fname = mw_dir_+"/"+hillsfname;
1438 : } else {
1439 : fname = hillsfname;
1440 : }
1441 : } else {
1442 113 : if(mw_n_>1) {
1443 9 : std::stringstream out;
1444 9 : out << i;
1445 18 : fname = hillsfname+"."+out.str();
1446 9 : } else {
1447 : fname = hillsfname;
1448 : }
1449 : }
1450 160 : ifiles_.emplace_back(Tools::make_unique<IFile>());
1451 : // this is just a shortcut pointer to the last element:
1452 : IFile *ifile = ifiles_.back().get();
1453 160 : ifilesnames_.push_back(fname);
1454 160 : ifile->link(*this);
1455 160 : if(ifile->FileExist(fname)) {
1456 36 : ifile->open(fname);
1457 36 : if(getRestart()&&!restartedFromGrid) {
1458 19 : log.printf(" Restarting from %s:",ifilesnames_[i].c_str());
1459 19 : readGaussians(ifiles_[i].get());
1460 : restartedFromHills=true;
1461 : }
1462 36 : ifiles_[i]->reset(false);
1463 : // close only the walker own hills file for later writing
1464 36 : if(i==mw_id_) {
1465 30 : ifiles_[i]->close();
1466 : }
1467 : } else {
1468 : // in case a file does not exist and we are restarting, complain that the file was not found
1469 124 : if(getRestart()&&!restartedFromGrid) {
1470 0 : error("restart file "+fname+" not found");
1471 : }
1472 : }
1473 : }
1474 :
1475 : // if we are restarting from FILE and using WALKERS_MPI we can check that all walkers have actually read the FILE
1476 148 : if(getRestart()&&walkers_mpi_) {
1477 9 : std::vector<int> restarted(mpi_nw_,0);
1478 9 : if(comm.Get_rank()==0) {
1479 6 : multi_sim_comm.Allgather(int(restartedFromHills), restarted);
1480 : }
1481 9 : comm.Bcast(restarted,0);
1482 : int result = std::accumulate(restarted.begin(),restarted.end(),0);
1483 9 : if(result!=0&&result!=mpi_nw_) {
1484 0 : error("in this WALKERS_MPI run some replica have restarted from FILE while other do not!");
1485 : }
1486 : }
1487 :
1488 148 : comm.Barrier();
1489 : // this barrier is needed when using walkers_mpi
1490 : // to be sure that all files have been read before
1491 : // backing them up
1492 : // it should not be used when walkers_mpi is false otherwise
1493 : // it would introduce troubles when using replicas without METAD
1494 : // (e.g. in bias exchange with a neutral replica)
1495 : // see issue #168 on github
1496 148 : if(comm.Get_rank()==0 && walkers_mpi_) {
1497 23 : multi_sim_comm.Barrier();
1498 : }
1499 :
1500 148 : if(targetfilename_.length()>0) {
1501 2 : IFile gridfile;
1502 2 : gridfile.open(targetfilename_);
1503 2 : std::string funcl=getLabel() + ".target";
1504 4 : TargetGrid_=GridBase::create(funcl,getArguments(),gridfile,false,false,true);
1505 2 : if(TargetGrid_->getDimension()!=getNumberOfArguments()) {
1506 0 : error("mismatch between dimensionality of input grid and number of arguments");
1507 : }
1508 4 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1509 4 : if( getPntrToArgument(i)->isPeriodic()!=TargetGrid_->getIsPeriodic()[i] ) {
1510 0 : error("periodicity mismatch between arguments and input bias");
1511 : }
1512 : }
1513 2 : }
1514 :
1515 148 : if(getRestart()) {
1516 : // if this is a restart the neighbor list should be immediately updated
1517 37 : if(nlist_) {
1518 1 : nlist_update_=true;
1519 : }
1520 : // Calculate the Tiwary-Parrinello reweighting factor if we are restarting from previous hills
1521 37 : if(calc_rct_) {
1522 0 : computeReweightingFactor();
1523 : }
1524 : // Calculate all special bias quantities desired if restarting with nonzero bias.
1525 37 : if(calc_max_bias_) {
1526 0 : max_bias_ = BiasGrid_->getMaxValue();
1527 0 : getPntrToComponent("maxbias")->set(max_bias_);
1528 : }
1529 37 : if(calc_transition_bias_) {
1530 13 : transition_bias_ = getTransitionBarrierBias();
1531 26 : getPntrToComponent("transbias")->set(transition_bias_);
1532 : }
1533 : }
1534 :
1535 : // open grid file for writing
1536 148 : if(wgridstride_>0) {
1537 19 : gridfile_.link(*this);
1538 19 : if(walkers_mpi_) {
1539 0 : int r=0;
1540 0 : if(comm.Get_rank()==0) {
1541 0 : r=multi_sim_comm.Get_rank();
1542 : }
1543 0 : comm.Bcast(r,0);
1544 0 : if(r>0) {
1545 : gridfilename_="/dev/null";
1546 : }
1547 0 : gridfile_.enforceSuffix("");
1548 : }
1549 19 : if(mw_n_>1) {
1550 0 : gridfile_.enforceSuffix("");
1551 : }
1552 19 : gridfile_.open(gridfilename_);
1553 : }
1554 :
1555 : // open hills file for writing
1556 148 : hillsOfile_.link(*this);
1557 148 : if(walkers_mpi_) {
1558 38 : int r=0;
1559 38 : if(comm.Get_rank()==0) {
1560 23 : r=multi_sim_comm.Get_rank();
1561 : }
1562 38 : comm.Bcast(r,0);
1563 38 : if(r>0) {
1564 25 : ifilesnames_[mw_id_]="/dev/null";
1565 : }
1566 76 : hillsOfile_.enforceSuffix("");
1567 : }
1568 148 : if(mw_n_>1) {
1569 12 : hillsOfile_.enforceSuffix("");
1570 : }
1571 148 : hillsOfile_.open(ifilesnames_[mw_id_]);
1572 148 : if(fmt_.length()>0) {
1573 117 : hillsOfile_.fmtField(fmt_);
1574 : }
1575 148 : hillsOfile_.addConstantField("multivariate");
1576 148 : hillsOfile_.addConstantField("kerneltype");
1577 148 : if(doInt_) {
1578 4 : hillsOfile_.addConstantField("lower_int").printField("lower_int",lowI_);
1579 4 : hillsOfile_.addConstantField("upper_int").printField("upper_int",uppI_);
1580 : }
1581 : hillsOfile_.setHeavyFlush();
1582 : // output periodicities of variables
1583 415 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1584 267 : hillsOfile_.setupPrintValue( getPntrToArgument(i) );
1585 : }
1586 :
1587 : bool concurrent=false;
1588 148 : const ActionSet&actionSet(plumed.getActionSet());
1589 1941 : for(const auto & p : actionSet)
1590 1867 : if(dynamic_cast<MetaD*>(p.get())) {
1591 : concurrent=true;
1592 : break;
1593 : }
1594 148 : if(concurrent) {
1595 74 : log<<" You are using concurrent metadynamics\n";
1596 : }
1597 148 : if(rect_biasf_.size()>0) {
1598 18 : if(walkers_mpi_) {
1599 12 : log<<" You are using RECT in its 'altruistic' implementation\n";
1600 : }{
1601 18 : log<<" You are using RECT\n";
1602 : }
1603 : }
1604 :
1605 296 : log<<" Bibliography "<<plumed.cite("Laio and Parrinello, PNAS 99, 12562 (2002)");
1606 148 : if(welltemp_) {
1607 76 : log<<plumed.cite("Barducci, Bussi, and Parrinello, Phys. Rev. Lett. 100, 020603 (2008)");
1608 : }
1609 148 : if(tt_specs_.is_active) {
1610 6 : log << plumed.cite("Dama, Rotskoff, Parrinello, and Voth, J. Chem. Theory Comput. 10, 3626 (2014)");
1611 6 : log << plumed.cite("Dama, Parrinello, and Voth, Phys. Rev. Lett. 112, 240602 (2014)");
1612 : }
1613 148 : if(mw_n_>1||walkers_mpi_) {
1614 88 : log<<plumed.cite("Raiteri, Laio, Gervasio, Micheletti, and Parrinello, J. Phys. Chem. B 110, 3533 (2006)");
1615 : }
1616 148 : if(adaptive_!=FlexibleBin::none) {
1617 42 : log<<plumed.cite("Branduardi, Bussi, and Parrinello, J. Chem. Theory Comput. 8, 2247 (2012)");
1618 : }
1619 148 : if(doInt_) {
1620 4 : log<<plumed.cite("Baftizadeh, Cossio, Pietrucci, and Laio, Curr. Phys. Chem. 2, 79 (2012)");
1621 : }
1622 148 : if(acceleration_) {
1623 8 : log<<plumed.cite("Pratyush and Parrinello, Phys. Rev. Lett. 111, 230602 (2013)");
1624 : }
1625 148 : if(calc_rct_) {
1626 12 : log<<plumed.cite("Pratyush and Parrinello, J. Phys. Chem. B, 119, 736 (2015)");
1627 : }
1628 148 : if(concurrent || rect_biasf_.size()>0) {
1629 160 : log<<plumed.cite("Gil-Ley and Bussi, J. Chem. Theory Comput. 11, 1077 (2015)");
1630 : }
1631 148 : if(rect_biasf_.size()>0 && walkers_mpi_) {
1632 24 : log<<plumed.cite("Hosek, Toulcova, Bortolato, and Spiwok, J. Phys. Chem. B 120, 2209 (2016)");
1633 : }
1634 148 : if(targetfilename_.length()>0) {
1635 4 : log<<plumed.cite("White, Dama, and Voth, J. Chem. Theory Comput. 11, 2451 (2015)");
1636 4 : log<<plumed.cite("Marinelli and Faraldo-Gómez, Biophys. J. 108, 2779 (2015)");
1637 4 : log<<plumed.cite("Gil-Ley, Bottaro, and Bussi, J. Chem. Theory Comput. 12, 2790 (2016)");
1638 : }
1639 148 : if(freq_adaptive_) {
1640 4 : log<<plumed.cite("Wang, Valsson, Tiwary, Parrinello, and Lindorff-Larsen, J. Chem. Phys. 149, 072309 (2018)");
1641 : }
1642 148 : log<<"\n";
1643 396 : }
1644 :
1645 151 : void MetaD::readTemperingSpecs(TemperingSpecs &t_specs) {
1646 : // Set global tempering parameters.
1647 151 : parse(t_specs.name_stem + "BIASFACTOR", t_specs.biasf);
1648 151 : if (t_specs.biasf != -1.0) {
1649 3 : if (kbt_ == 0.0) {
1650 0 : error("Unless the MD engine passes the temperature to plumed, with tempered metad you must specify it using TEMP");
1651 : }
1652 3 : if (t_specs.biasf == 1.0) {
1653 0 : error("A bias factor of 1 corresponds to zero delta T and zero hill size, so it is not allowed.");
1654 : }
1655 3 : t_specs.is_active = true;
1656 3 : parse(t_specs.name_stem + "BIASTHRESHOLD", t_specs.threshold);
1657 3 : if (t_specs.threshold < 0.0) {
1658 0 : error(t_specs.name + " bias threshold is nonsensical");
1659 : }
1660 3 : parse(t_specs.name_stem + "ALPHA", t_specs.alpha);
1661 3 : if (t_specs.alpha <= 0.0 || t_specs.alpha > 1.0) {
1662 0 : error(t_specs.name + " decay shape parameter alpha is nonsensical");
1663 : }
1664 : }
1665 151 : }
1666 :
1667 3 : void MetaD::logTemperingSpecs(const TemperingSpecs &t_specs) {
1668 3 : log.printf(" %s bias factor %f\n", t_specs.name.c_str(), t_specs.biasf);
1669 3 : log.printf(" KbT %f\n", kbt_);
1670 3 : if (t_specs.threshold != 0.0) {
1671 2 : log.printf(" %s bias threshold %f\n", t_specs.name.c_str(), t_specs.threshold);
1672 : }
1673 3 : if (t_specs.alpha != 1.0) {
1674 1 : log.printf(" %s decay shape parameter alpha %f\n", t_specs.name.c_str(), t_specs.alpha);
1675 : }
1676 3 : }
1677 :
1678 6037 : void MetaD::readGaussians(IFile *ifile) {
1679 6037 : unsigned ncv=getNumberOfArguments();
1680 6037 : std::vector<double> center(ncv);
1681 6037 : std::vector<double> sigma(ncv);
1682 : double height;
1683 : int nhills=0;
1684 6037 : bool multivariate=false;
1685 :
1686 : std::vector<Value> tmpvalues;
1687 18124 : for(unsigned j=0; j<getNumberOfArguments(); ++j) {
1688 24174 : tmpvalues.push_back( Value( this, getPntrToArgument(j)->getName(), false ) );
1689 : }
1690 :
1691 8278 : while(scanOneHill(ifile,tmpvalues,center,sigma,height,multivariate)) {
1692 2241 : nhills++;
1693 : // note that for gamma=1 we store directly -F
1694 2241 : if(welltemp_ && biasf_>1.0) {
1695 41 : height*=(biasf_-1.0)/biasf_;
1696 : }
1697 2241 : addGaussian(Gaussian(multivariate,height,center,sigma));
1698 : }
1699 6037 : log.printf(" %d Gaussians read\n",nhills);
1700 12074 : }
1701 :
1702 2922 : void MetaD::writeGaussian(const Gaussian& hill, OFile&file) {
1703 2922 : unsigned ncv=getNumberOfArguments();
1704 2922 : file.printField("time",getTimeStep()*getStep());
1705 8194 : for(unsigned i=0; i<ncv; ++i) {
1706 5272 : file.printField(getPntrToArgument(i),hill.center[i]);
1707 : }
1708 5844 : hillsOfile_.printField("kerneltype","stretched-gaussian");
1709 2922 : if(hill.multivariate) {
1710 892 : hillsOfile_.printField("multivariate","true");
1711 : Matrix<double> mymatrix(ncv,ncv);
1712 : unsigned k=0;
1713 1047 : for(unsigned i=0; i<ncv; i++) {
1714 1357 : for(unsigned j=i; j<ncv; j++) {
1715 : // recompose the full inverse matrix
1716 756 : mymatrix(i,j)=mymatrix(j,i)=hill.sigma[k];
1717 756 : k++;
1718 : }
1719 : }
1720 : // invert it
1721 : Matrix<double> invmatrix(ncv,ncv);
1722 446 : Invert(mymatrix,invmatrix);
1723 : // enforce symmetry
1724 1047 : for(unsigned i=0; i<ncv; i++) {
1725 1357 : for(unsigned j=i; j<ncv; j++) {
1726 756 : invmatrix(i,j)=invmatrix(j,i);
1727 : }
1728 : }
1729 :
1730 : // do cholesky so to have a "sigma like" number
1731 : Matrix<double> lower(ncv,ncv);
1732 446 : cholesky(invmatrix,lower);
1733 : // loop in band form
1734 1047 : for(unsigned i=0; i<ncv; i++) {
1735 1357 : for(unsigned j=0; j<ncv-i; j++) {
1736 1512 : file.printField("sigma_"+getPntrToArgument(j+i)->getName()+"_"+getPntrToArgument(j)->getName(),lower(j+i,j));
1737 : }
1738 : }
1739 : } else {
1740 4952 : hillsOfile_.printField("multivariate","false");
1741 7147 : for(unsigned i=0; i<ncv; ++i) {
1742 9342 : file.printField("sigma_"+getPntrToArgument(i)->getName(),hill.sigma[i]);
1743 : }
1744 : }
1745 2922 : double height=hill.height;
1746 : // note that for gamma=1 we store directly -F
1747 2922 : if(welltemp_ && biasf_>1.0) {
1748 339 : height*=biasf_/(biasf_-1.0);
1749 : }
1750 5844 : file.printField("height",height).printField("biasf",biasf_);
1751 2922 : if(mw_n_>1) {
1752 3018 : file.printField("clock",int(std::time(0)));
1753 : }
1754 2922 : file.printField();
1755 2922 : }
1756 :
1757 5325 : void MetaD::addGaussian(const Gaussian& hill) {
1758 5325 : if(grid_) {
1759 640 : size_t ncv=getNumberOfArguments();
1760 640 : std::vector<unsigned> nneighb=getGaussianSupport(hill);
1761 640 : std::vector<Grid::index_t> neighbors=BiasGrid_->getNeighbors(hill.center,nneighb);
1762 640 : std::vector<double> der(ncv);
1763 640 : std::vector<double> xx(ncv);
1764 640 : if(comm.Get_size()==1) {
1765 : // for performance reasons and thread safety
1766 544 : std::vector<double> dp(ncv);
1767 55324 : for(size_t i=0; i<neighbors.size(); ++i) {
1768 54780 : Grid::index_t ineigh=neighbors[i];
1769 158922 : for(unsigned j=0; j<ncv; ++j) {
1770 104142 : der[j]=0.0;
1771 : }
1772 54780 : BiasGrid_->getPoint(ineigh,xx);
1773 54780 : double bias=evaluateGaussianAndDerivatives(xx,hill,der,dp);
1774 54780 : BiasGrid_->addValueAndDerivatives(ineigh,bias,der);
1775 : }
1776 : } else {
1777 96 : unsigned stride=comm.Get_size();
1778 96 : unsigned rank=comm.Get_rank();
1779 96 : std::vector<double> allder(ncv*neighbors.size(),0.0);
1780 96 : std::vector<double> n_der(ncv,0.0);
1781 96 : std::vector<double> allbias(neighbors.size(),0.0);
1782 : // for performance reasons and thread safety
1783 96 : std::vector<double> dp(ncv);
1784 27148 : for(unsigned i=rank; i<neighbors.size(); i+=stride) {
1785 27052 : Grid::index_t ineigh=neighbors[i];
1786 81156 : for(unsigned j=0; j<ncv; ++j) {
1787 54104 : n_der[j]=0.0;
1788 : }
1789 27052 : BiasGrid_->getPoint(ineigh,xx);
1790 27052 : allbias[i]=evaluateGaussianAndDerivatives(xx,hill,n_der,dp);
1791 81156 : for(unsigned j=0; j<ncv; j++) {
1792 54104 : allder[ncv*i+j]=n_der[j];
1793 : }
1794 : }
1795 96 : comm.Sum(allbias);
1796 96 : comm.Sum(allder);
1797 103200 : for(unsigned i=0; i<neighbors.size(); ++i) {
1798 103104 : Grid::index_t ineigh=neighbors[i];
1799 309312 : for(unsigned j=0; j<ncv; ++j) {
1800 206208 : der[j]=allder[ncv*i+j];
1801 : }
1802 103104 : BiasGrid_->addValueAndDerivatives(ineigh,allbias[i],der);
1803 : }
1804 : }
1805 : } else {
1806 4685 : hills_.push_back(hill);
1807 : }
1808 5325 : }
1809 :
1810 640 : std::vector<unsigned> MetaD::getGaussianSupport(const Gaussian& hill) {
1811 : std::vector<unsigned> nneigh;
1812 : std::vector<double> cutoff;
1813 640 : unsigned ncv=getNumberOfArguments();
1814 :
1815 : // traditional or flexible hill?
1816 640 : if(hill.multivariate) {
1817 : unsigned k=0;
1818 : Matrix<double> mymatrix(ncv,ncv);
1819 0 : for(unsigned i=0; i<ncv; i++) {
1820 0 : for(unsigned j=i; j<ncv; j++) {
1821 : // recompose the full inverse matrix
1822 0 : mymatrix(i,j)=mymatrix(j,i)=hill.sigma[k];
1823 0 : k++;
1824 : }
1825 : }
1826 : // Reinvert so to have the ellipses
1827 : Matrix<double> myinv(ncv,ncv);
1828 0 : Invert(mymatrix,myinv);
1829 : Matrix<double> myautovec(ncv,ncv);
1830 0 : std::vector<double> myautoval(ncv); //should I take this or their square root?
1831 0 : diagMat(myinv,myautoval,myautovec);
1832 : double maxautoval=0.;
1833 : unsigned ind_maxautoval;
1834 : ind_maxautoval=ncv;
1835 0 : for(unsigned i=0; i<ncv; i++) {
1836 0 : if(myautoval[i]>maxautoval) {
1837 : maxautoval=myautoval[i];
1838 : ind_maxautoval=i;
1839 : }
1840 : }
1841 0 : for(unsigned i=0; i<ncv; i++) {
1842 0 : cutoff.push_back(std::sqrt(2.0*dp2cutoff)*std::abs(std::sqrt(maxautoval)*myautovec(i,ind_maxautoval)));
1843 : }
1844 : } else {
1845 1618 : for(unsigned i=0; i<ncv; ++i) {
1846 978 : cutoff.push_back(std::sqrt(2.0*dp2cutoff)*hill.sigma[i]);
1847 : }
1848 : }
1849 :
1850 640 : if(doInt_) {
1851 2 : if(hill.center[0]+cutoff[0] > uppI_ || hill.center[0]-cutoff[0] < lowI_) {
1852 : // in this case, we updated the entire grid to avoid problems
1853 2 : return BiasGrid_->getNbin();
1854 : } else {
1855 0 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff[0]/BiasGrid_->getDx()[0])) );
1856 : return nneigh;
1857 : }
1858 : } else {
1859 1614 : for(unsigned i=0; i<ncv; i++) {
1860 976 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff[i]/BiasGrid_->getDx()[i])) );
1861 : }
1862 : }
1863 :
1864 : return nneigh;
1865 : }
1866 :
1867 285 : double MetaD::getBias(const std::vector<double>& cv) {
1868 285 : double bias=0.0;
1869 285 : if(grid_) {
1870 203 : bias = BiasGrid_->getValue(cv);
1871 : } else {
1872 82 : unsigned nt=OpenMP::getNumThreads();
1873 82 : unsigned stride=comm.Get_size();
1874 82 : unsigned rank=comm.Get_rank();
1875 :
1876 82 : if(!nlist_) {
1877 82 : #pragma omp parallel num_threads(nt)
1878 : {
1879 : #pragma omp for reduction(+:bias) nowait
1880 : for(unsigned i=rank; i<hills_.size(); i+=stride) {
1881 : bias+=evaluateGaussian(cv,hills_[i]);
1882 : }
1883 : }
1884 : } else {
1885 0 : #pragma omp parallel num_threads(nt)
1886 : {
1887 : #pragma omp for reduction(+:bias) nowait
1888 : for(unsigned i=rank; i<nlist_hills_.size(); i+=stride) {
1889 : bias+=evaluateGaussian(cv,nlist_hills_[i]);
1890 : }
1891 : }
1892 : }
1893 82 : comm.Sum(bias);
1894 : }
1895 :
1896 285 : return bias;
1897 : }
1898 :
1899 8395 : double MetaD::getBiasAndDerivatives(const std::vector<double>& cv, std::vector<double>& der) {
1900 8395 : unsigned ncv=getNumberOfArguments();
1901 8395 : double bias=0.0;
1902 8395 : if(grid_) {
1903 1506 : std::vector<double> vder(ncv);
1904 1506 : bias=BiasGrid_->getValueAndDerivatives(cv,vder);
1905 3498 : for(unsigned i=0; i<ncv; i++) {
1906 1992 : der[i]=vder[i];
1907 : }
1908 : } else {
1909 6889 : unsigned nt=OpenMP::getNumThreads();
1910 6889 : unsigned stride=comm.Get_size();
1911 6889 : unsigned rank=comm.Get_rank();
1912 :
1913 6889 : if(!nlist_) {
1914 6884 : if(hills_.size()<2*nt*stride||nt==1) {
1915 : // for performance reasons and thread safety
1916 2586 : std::vector<double> dp(ncv);
1917 6703 : for(unsigned i=rank; i<hills_.size(); i+=stride) {
1918 4117 : bias+=evaluateGaussianAndDerivatives(cv,hills_[i],der,dp);
1919 : }
1920 : } else {
1921 4298 : #pragma omp parallel num_threads(nt)
1922 : {
1923 : std::vector<double> omp_deriv(ncv,0.);
1924 : // for performance reasons and thread safety
1925 : std::vector<double> dp(ncv);
1926 : #pragma omp for reduction(+:bias) nowait
1927 : for(unsigned i=rank; i<hills_.size(); i+=stride) {
1928 : bias+=evaluateGaussianAndDerivatives(cv,hills_[i],omp_deriv,dp);
1929 : }
1930 : #pragma omp critical
1931 : for(unsigned i=0; i<ncv; i++) {
1932 : der[i]+=omp_deriv[i];
1933 : }
1934 : }
1935 : }
1936 : } else {
1937 5 : if(hills_.size()<2*nt*stride||nt==1) {
1938 : // for performance reasons and thread safety
1939 0 : std::vector<double> dp(ncv);
1940 0 : for(unsigned i=rank; i<nlist_hills_.size(); i+=stride) {
1941 0 : bias+=evaluateGaussianAndDerivatives(cv,nlist_hills_[i],der,dp);
1942 : }
1943 : } else {
1944 5 : #pragma omp parallel num_threads(nt)
1945 : {
1946 : std::vector<double> omp_deriv(ncv,0.);
1947 : // for performance reasons and thread safety
1948 : std::vector<double> dp(ncv);
1949 : #pragma omp for reduction(+:bias) nowait
1950 : for(unsigned i=rank; i<nlist_hills_.size(); i+=stride) {
1951 : bias+=evaluateGaussianAndDerivatives(cv,nlist_hills_[i],omp_deriv,dp);
1952 : }
1953 : #pragma omp critical
1954 : for(unsigned i=0; i<ncv; i++) {
1955 : der[i]+=omp_deriv[i];
1956 : }
1957 : }
1958 : }
1959 : }
1960 6889 : comm.Sum(bias);
1961 6889 : comm.Sum(der);
1962 : }
1963 :
1964 8395 : return bias;
1965 : }
1966 :
1967 0 : double MetaD::getGaussianNormalization(const Gaussian& hill) {
1968 : double norm=1;
1969 0 : unsigned ncv=hill.center.size();
1970 :
1971 0 : if(hill.multivariate) {
1972 : // recompose the full sigma from the upper diag cholesky
1973 : unsigned k=0;
1974 : Matrix<double> mymatrix(ncv,ncv);
1975 0 : for(unsigned i=0; i<ncv; i++) {
1976 0 : for(unsigned j=i; j<ncv; j++) {
1977 0 : mymatrix(i,j)=mymatrix(j,i)=hill.sigma[k]; // recompose the full inverse matrix
1978 0 : k++;
1979 : }
1980 : double ldet;
1981 0 : logdet( mymatrix, ldet );
1982 0 : norm = std::exp( ldet ); // Not sure here if mymatrix is sigma or inverse
1983 : }
1984 : } else {
1985 0 : for(unsigned i=0; i<hill.sigma.size(); i++) {
1986 0 : norm*=hill.sigma[i];
1987 : }
1988 : }
1989 :
1990 0 : return norm*std::pow(2*pi,static_cast<double>(ncv)/2.0);
1991 : }
1992 :
1993 192 : double MetaD::evaluateGaussian(const std::vector<double>& cv, const Gaussian& hill) {
1994 192 : unsigned ncv=cv.size();
1995 :
1996 : // I use a pointer here because cv is const (and should be const)
1997 : // but when using doInt it is easier to locally replace cv[0] with
1998 : // the upper/lower limit in case it is out of range
1999 : double tmpcv[1];
2000 : const double *pcv=NULL; // pointer to cv
2001 192 : if(ncv>0) {
2002 : pcv=&cv[0];
2003 : }
2004 192 : if(doInt_) {
2005 0 : plumed_assert(ncv==1);
2006 0 : tmpcv[0]=cv[0];
2007 0 : if(cv[0]<lowI_) {
2008 0 : tmpcv[0]=lowI_;
2009 : }
2010 0 : if(cv[0]>uppI_) {
2011 0 : tmpcv[0]=uppI_;
2012 : }
2013 : pcv=&(tmpcv[0]);
2014 : }
2015 :
2016 : double dp2=0.0;
2017 192 : if(hill.multivariate) {
2018 : unsigned k=0;
2019 : // recompose the full sigma from the upper diag cholesky
2020 : Matrix<double> mymatrix(ncv,ncv);
2021 0 : for(unsigned i=0; i<ncv; i++) {
2022 0 : for(unsigned j=i; j<ncv; j++) {
2023 0 : mymatrix(i,j)=mymatrix(j,i)=hill.sigma[k]; // recompose the full inverse matrix
2024 0 : k++;
2025 : }
2026 : }
2027 0 : for(unsigned i=0; i<ncv; i++) {
2028 0 : double dp_i=difference(i,hill.center[i],pcv[i]);
2029 0 : for(unsigned j=i; j<ncv; j++) {
2030 0 : if(i==j) {
2031 0 : dp2+=dp_i*dp_i*mymatrix(i,j)*0.5;
2032 : } else {
2033 0 : double dp_j=difference(j,hill.center[j],pcv[j]);
2034 0 : dp2+=dp_i*dp_j*mymatrix(i,j);
2035 : }
2036 : }
2037 : }
2038 : } else {
2039 576 : for(unsigned i=0; i<ncv; i++) {
2040 384 : double dp=difference(i,hill.center[i],pcv[i])*hill.invsigma[i];
2041 384 : dp2+=dp*dp;
2042 : }
2043 192 : dp2*=0.5;
2044 : }
2045 :
2046 : double bias=0.0;
2047 192 : if(dp2<dp2cutoff) {
2048 154 : bias=hill.height*(stretchA*std::exp(-dp2)+stretchB);
2049 : }
2050 :
2051 192 : return bias;
2052 : }
2053 :
2054 2408957 : double MetaD::evaluateGaussianAndDerivatives(const std::vector<double>& cv, const Gaussian& hill, std::vector<double>& der, std::vector<double>& dp_) {
2055 2408957 : unsigned ncv=cv.size();
2056 :
2057 : // I use a pointer here because cv is const (and should be const)
2058 : // but when using doInt it is easier to locally replace cv[0] with
2059 : // the upper/lower limit in case it is out of range
2060 : const double *pcv=NULL; // pointer to cv
2061 : double tmpcv[1]; // tmp array with cv (to be used with doInt_)
2062 2408957 : if(ncv>0) {
2063 : pcv=&cv[0];
2064 : }
2065 2408957 : if(doInt_) {
2066 602 : plumed_assert(ncv==1);
2067 602 : tmpcv[0]=cv[0];
2068 602 : if(cv[0]<lowI_) {
2069 118 : tmpcv[0]=lowI_;
2070 : }
2071 602 : if(cv[0]>uppI_) {
2072 360 : tmpcv[0]=uppI_;
2073 : }
2074 : pcv=&(tmpcv[0]);
2075 : }
2076 :
2077 : bool int_der=false;
2078 2408957 : if(doInt_) {
2079 602 : if(cv[0]<lowI_ || cv[0]>uppI_) {
2080 : int_der=true;
2081 : }
2082 : }
2083 :
2084 : double dp2=0.0;
2085 : double bias=0.0;
2086 2408957 : if(hill.multivariate) {
2087 : unsigned k=0;
2088 : // recompose the full sigma from the upper diag cholesky
2089 : Matrix<double> mymatrix(ncv,ncv);
2090 161635 : for(unsigned i=0; i<ncv; i++) {
2091 162513 : for(unsigned j=i; j<ncv; j++) {
2092 81476 : mymatrix(i,j)=mymatrix(j,i)=hill.sigma[k]; // recompose the full inverse matrix
2093 81476 : k++;
2094 : }
2095 : }
2096 161635 : for(unsigned i=0; i<ncv; i++) {
2097 81037 : dp_[i]=difference(i,hill.center[i],pcv[i]);
2098 162513 : for(unsigned j=i; j<ncv; j++) {
2099 81476 : if(i==j) {
2100 81037 : dp2+=dp_[i]*dp_[i]*mymatrix(i,j)*0.5;
2101 : } else {
2102 439 : double dp_j=difference(j,hill.center[j],pcv[j]);
2103 439 : dp2+=dp_[i]*dp_j*mymatrix(i,j);
2104 : }
2105 : }
2106 : }
2107 80598 : if(dp2<dp2cutoff) {
2108 77683 : bias=hill.height*std::exp(-dp2);
2109 77683 : if(!int_der) {
2110 155673 : for(unsigned i=0; i<ncv; i++) {
2111 : double tmp=0.0;
2112 156594 : for(unsigned j=0; j<ncv; j++) {
2113 78604 : tmp += dp_[j]*mymatrix(i,j)*bias;
2114 : }
2115 77990 : der[i]-=tmp*stretchA;
2116 : }
2117 : } else {
2118 0 : for(unsigned i=0; i<ncv; i++) {
2119 0 : der[i]=0.;
2120 : }
2121 : }
2122 77683 : bias=stretchA*bias+hill.height*stretchB;
2123 : }
2124 : } else {
2125 6975435 : for(unsigned i=0; i<ncv; i++) {
2126 4647076 : dp_[i]=difference(i,hill.center[i],pcv[i])*hill.invsigma[i];
2127 4647076 : dp2+=dp_[i]*dp_[i];
2128 : }
2129 2328359 : dp2*=0.5;
2130 2328359 : if(dp2<dp2cutoff) {
2131 1356330 : bias=hill.height*std::exp(-dp2);
2132 1356330 : if(!int_der) {
2133 4060211 : for(unsigned i=0; i<ncv; i++) {
2134 2704120 : der[i]-=bias*dp_[i]*hill.invsigma[i]*stretchA;
2135 : }
2136 : } else {
2137 478 : for(unsigned i=0; i<ncv; i++) {
2138 239 : der[i]=0.;
2139 : }
2140 : }
2141 1356330 : bias=stretchA*bias+hill.height*stretchB;
2142 : }
2143 : }
2144 :
2145 2408957 : return bias;
2146 : }
2147 :
2148 2736 : double MetaD::getHeight(const std::vector<double>& cv) {
2149 2736 : double height=height0_;
2150 2736 : if(welltemp_) {
2151 275 : double vbias = getBias(cv);
2152 275 : if(biasf_>1.0) {
2153 259 : height = height0_*std::exp(-vbias/(kbt_*(biasf_-1.0)));
2154 : } else {
2155 : // notice that if gamma=1 we store directly -F
2156 16 : height = height0_*std::exp(-vbias/kbt_);
2157 : }
2158 : }
2159 2736 : if(dampfactor_>0.0) {
2160 18 : plumed_assert(BiasGrid_);
2161 18 : double m=BiasGrid_->getMaxValue();
2162 18 : height*=std::exp(-m/(kbt_*(dampfactor_)));
2163 : }
2164 2736 : if (tt_specs_.is_active) {
2165 60 : double vbarrier = transition_bias_;
2166 60 : temperHeight(height, tt_specs_, vbarrier);
2167 : }
2168 2736 : if(TargetGrid_) {
2169 18 : double f=TargetGrid_->getValue(cv)-TargetGrid_->getMaxValue();
2170 18 : height*=std::exp(f/kbt_);
2171 : }
2172 2736 : return height;
2173 : }
2174 :
2175 60 : void MetaD::temperHeight(double& height, const TemperingSpecs& t_specs, const double tempering_bias) {
2176 60 : if (t_specs.alpha == 1.0) {
2177 80 : height *= std::exp(-std::max(0.0, tempering_bias - t_specs.threshold) / (kbt_ * (t_specs.biasf - 1.0)));
2178 : } else {
2179 40 : height *= std::pow(1 + (1 - t_specs.alpha) / t_specs.alpha * std::max(0.0, tempering_bias - t_specs.threshold) / (kbt_ * (t_specs.biasf - 1.0)), - t_specs.alpha / (1 - t_specs.alpha));
2180 : }
2181 60 : }
2182 :
2183 8435 : void MetaD::calculate() {
2184 : // this is because presently there is no way to properly pass information
2185 : // on adaptive hills (diff) after exchanges:
2186 8435 : if(adaptive_==FlexibleBin::diffusion && getExchangeStep()) {
2187 0 : error("ADAPTIVE=DIFF is not compatible with replica exchange");
2188 : }
2189 :
2190 8435 : const unsigned ncv=getNumberOfArguments();
2191 8435 : std::vector<double> cv(ncv);
2192 21082 : for(unsigned i=0; i<ncv; ++i) {
2193 12647 : cv[i]=getArgument(i);
2194 : }
2195 :
2196 8435 : if(nlist_) {
2197 5 : nlist_steps_++;
2198 5 : if(getExchangeStep()) {
2199 0 : nlist_update_=true;
2200 : } else {
2201 11 : for(unsigned i=0; i<ncv; ++i) {
2202 8 : double d = difference(i, cv[i], nlist_center_[i]);
2203 8 : double nk_dist2 = d*d/nlist_dev2_[i];
2204 8 : if(nk_dist2>nlist_param_[1]) {
2205 2 : nlist_update_=true;
2206 2 : break;
2207 : }
2208 : }
2209 : }
2210 5 : if(nlist_update_) {
2211 4 : updateNlist();
2212 : }
2213 : }
2214 :
2215 : double ene = 0.;
2216 8435 : std::vector<double> der(ncv,0.);
2217 8435 : if(biasf_!=1.0) {
2218 8395 : ene = getBiasAndDerivatives(cv,der);
2219 : }
2220 8435 : setBias(ene);
2221 21082 : for(unsigned i=0; i<ncv; i++) {
2222 12647 : setOutputForce(i,-der[i]);
2223 : }
2224 :
2225 8435 : if(calc_work_) {
2226 10 : getPntrToComponent("work")->set(work_);
2227 : }
2228 8435 : if(calc_rct_) {
2229 220 : getPntrToComponent("rbias")->set(ene - reweight_factor_);
2230 : }
2231 : // calculate the acceleration factor
2232 8435 : if(acceleration_&&!isFirstStep_) {
2233 329 : acc_ += static_cast<double>(getStride()) * std::exp(ene/(kbt_));
2234 329 : const double mean_acc = acc_/((double) getStep());
2235 329 : getPntrToComponent("acc")->set(mean_acc);
2236 8435 : } else if (acceleration_ && isFirstStep_ && acc_restart_mean_ > 0.0) {
2237 2 : acc_ = acc_restart_mean_ * static_cast<double>(getStep());
2238 2 : if(freq_adaptive_) {
2239 : // has to be done here if restarting, as the acc is not defined before
2240 1 : updateFrequencyAdaptiveStride();
2241 : }
2242 : }
2243 8435 : }
2244 :
2245 6239 : void MetaD::update() {
2246 : // adding hills criteria (could be more complex though)
2247 : bool nowAddAHill;
2248 6239 : if(getStep()%current_stride_==0 && !isFirstStep_) {
2249 : nowAddAHill=true;
2250 : } else {
2251 : nowAddAHill=false;
2252 3503 : isFirstStep_=false;
2253 : }
2254 :
2255 6239 : unsigned ncv=getNumberOfArguments();
2256 6239 : std::vector<double> cv(ncv);
2257 16690 : for(unsigned i=0; i<ncv; ++i) {
2258 10451 : cv[i] = getArgument(i);
2259 : }
2260 :
2261 : double vbias=0.;
2262 6239 : if(calc_work_) {
2263 5 : vbias=getBias(cv);
2264 : }
2265 :
2266 : // if you use adaptive, call the FlexibleBin
2267 : bool multivariate=false;
2268 6239 : if(adaptive_!=FlexibleBin::none) {
2269 778 : flexbin_->update(nowAddAHill);
2270 : multivariate=true;
2271 : }
2272 :
2273 : std::vector<double> thissigma;
2274 6239 : if(nowAddAHill) {
2275 : // add a Gaussian
2276 2736 : double height=getHeight(cv);
2277 : // returns upper diagonal inverse
2278 2736 : if(adaptive_!=FlexibleBin::none) {
2279 748 : thissigma=flexbin_->getInverseMatrix();
2280 : }
2281 : // returns normal sigma
2282 : else {
2283 2362 : thissigma=sigma0_;
2284 : }
2285 :
2286 : // In case we use walkers_mpi, it is now necessary to communicate with other replicas.
2287 2736 : if(walkers_mpi_) {
2288 : // Allocate arrays to store all walkers hills
2289 174 : std::vector<double> all_cv(mpi_nw_*ncv,0.0);
2290 174 : std::vector<double> all_sigma(mpi_nw_*thissigma.size(),0.0);
2291 174 : std::vector<double> all_height(mpi_nw_,0.0);
2292 174 : std::vector<int> all_multivariate(mpi_nw_,0);
2293 174 : if(comm.Get_rank()==0) {
2294 : // Communicate (only root)
2295 99 : multi_sim_comm.Allgather(cv,all_cv);
2296 99 : multi_sim_comm.Allgather(thissigma,all_sigma);
2297 : // notice that if gamma=1 we store directly -F so this scaling is not necessary:
2298 99 : multi_sim_comm.Allgather(height*(biasf_>1.0?biasf_/(biasf_-1.0):1.0),all_height);
2299 99 : multi_sim_comm.Allgather(int(multivariate),all_multivariate);
2300 : }
2301 : // Share info with group members
2302 174 : comm.Bcast(all_cv,0);
2303 174 : comm.Bcast(all_sigma,0);
2304 174 : comm.Bcast(all_height,0);
2305 174 : comm.Bcast(all_multivariate,0);
2306 :
2307 : // Flying Gaussian
2308 174 : if (flying_) {
2309 54 : hills_.clear();
2310 54 : comm.Barrier();
2311 : }
2312 :
2313 696 : for(unsigned i=0; i<mpi_nw_; i++) {
2314 : // actually add hills one by one
2315 522 : std::vector<double> cv_now(ncv);
2316 522 : std::vector<double> sigma_now(thissigma.size());
2317 1566 : for(unsigned j=0; j<ncv; j++) {
2318 1044 : cv_now[j]=all_cv[i*ncv+j];
2319 : }
2320 1674 : for(unsigned j=0; j<thissigma.size(); j++) {
2321 1152 : sigma_now[j]=all_sigma[i*thissigma.size()+j];
2322 : }
2323 : // notice that if gamma=1 we store directly -F so this scaling is not necessary:
2324 522 : double fact=(biasf_>1.0?(biasf_-1.0)/biasf_:1.0);
2325 522 : Gaussian newhill=Gaussian(all_multivariate[i],all_height[i]*fact,cv_now,sigma_now);
2326 522 : addGaussian(newhill);
2327 522 : if(!flying_) {
2328 360 : writeGaussian(newhill,hillsOfile_);
2329 : }
2330 522 : }
2331 : } else {
2332 2562 : Gaussian newhill=Gaussian(multivariate,height,cv,thissigma);
2333 2562 : addGaussian(newhill);
2334 2562 : writeGaussian(newhill,hillsOfile_);
2335 2562 : }
2336 :
2337 : // this is to update the hills neighbor list
2338 2736 : if(nlist_) {
2339 4 : nlist_update_=true;
2340 : }
2341 : }
2342 :
2343 : // this should be outside of the if block in case
2344 : // mw_rstride_ is not a multiple of stride_
2345 6239 : if(mw_n_>1 && getStep()%mw_rstride_==0) {
2346 3012 : hillsOfile_.flush();
2347 : }
2348 :
2349 6239 : if(calc_work_) {
2350 5 : if(nlist_) {
2351 0 : updateNlist();
2352 : }
2353 5 : double vbias1=getBias(cv);
2354 5 : work_+=vbias1-vbias;
2355 : }
2356 :
2357 : // dump grid on file
2358 6239 : if(wgridstride_>0&&(getStep()%wgridstride_==0||getCPT())) {
2359 : // in case old grids are stored, a sequence of grids should appear
2360 : // this call results in a repetition of the header:
2361 91 : if(storeOldGrids_) {
2362 40 : gridfile_.clearFields();
2363 : }
2364 : // in case only latest grid is stored, file should be rewound
2365 : // this will overwrite previously written grids
2366 : else {
2367 51 : int r = 0;
2368 51 : if(walkers_mpi_) {
2369 0 : if(comm.Get_rank()==0) {
2370 0 : r=multi_sim_comm.Get_rank();
2371 : }
2372 0 : comm.Bcast(r,0);
2373 : }
2374 51 : if(r==0) {
2375 51 : gridfile_.rewind();
2376 : }
2377 : }
2378 91 : BiasGrid_->writeToFile(gridfile_);
2379 : // if a single grid is stored, it is necessary to flush it, otherwise
2380 : // the file might stay empty forever (when a single grid is not large enough to
2381 : // trigger flushing from the operating system).
2382 : // on the other hand, if grids are stored one after the other this is
2383 : // no necessary, and we leave the flushing control to the user as usual
2384 : // (with FLUSH keyword)
2385 91 : if(!storeOldGrids_) {
2386 51 : gridfile_.flush();
2387 : }
2388 : }
2389 :
2390 : // if multiple walkers and time to read Gaussians
2391 6239 : if(mw_n_>1 && getStep()%mw_rstride_==0) {
2392 12048 : for(int i=0; i<mw_n_; ++i) {
2393 : // don't read your own Gaussians
2394 9036 : if(i==mw_id_) {
2395 3012 : continue;
2396 : }
2397 : // if the file is not open yet
2398 6024 : if(!(ifiles_[i]->isOpen())) {
2399 : // check if it exists now and open it!
2400 6 : if(ifiles_[i]->FileExist(ifilesnames_[i])) {
2401 6 : ifiles_[i]->open(ifilesnames_[i]);
2402 6 : ifiles_[i]->reset(false);
2403 : }
2404 : // otherwise read the new Gaussians
2405 : } else {
2406 6018 : log.printf(" Reading hills from %s:",ifilesnames_[i].c_str());
2407 6018 : readGaussians(ifiles_[i].get());
2408 6018 : ifiles_[i]->reset(false);
2409 : }
2410 : }
2411 : // this is to update the hills neighbor list
2412 3012 : if(nlist_) {
2413 0 : nlist_update_=true;
2414 : }
2415 : }
2416 :
2417 : // Recalculate special bias quantities whenever the bias has been changed by the update.
2418 6239 : bool bias_has_changed = (nowAddAHill || (mw_n_ > 1 && getStep() % mw_rstride_ == 0));
2419 6239 : if (calc_rct_ && bias_has_changed && getStep()%(stride_*rct_ustride_)==0) {
2420 102 : computeReweightingFactor();
2421 : }
2422 6239 : if (calc_max_bias_ && bias_has_changed) {
2423 0 : max_bias_ = BiasGrid_->getMaxValue();
2424 0 : getPntrToComponent("maxbias")->set(max_bias_);
2425 : }
2426 6239 : if (calc_transition_bias_ && bias_has_changed) {
2427 260 : transition_bias_ = getTransitionBarrierBias();
2428 520 : getPntrToComponent("transbias")->set(transition_bias_);
2429 : }
2430 :
2431 : // Frequency adaptive metadynamics - update hill addition frequency
2432 6239 : if(freq_adaptive_ && getStep()%fa_update_frequency_==0) {
2433 151 : updateFrequencyAdaptiveStride();
2434 : }
2435 6239 : }
2436 :
2437 : /// takes a pointer to the file and a template std::string with values v and gives back the next center, sigma and height
2438 8278 : bool MetaD::scanOneHill(IFile* ifile, std::vector<Value>& tmpvalues, std::vector<double>& center, std::vector<double>& sigma, double& height, bool& multivariate) {
2439 : double dummy;
2440 8278 : multivariate=false;
2441 16556 : if(ifile->scanField("time",dummy)) {
2442 2241 : unsigned ncv=tmpvalues.size();
2443 6677 : for(unsigned i=0; i<ncv; ++i) {
2444 4436 : ifile->scanField( &tmpvalues[i] );
2445 4436 : if( tmpvalues[i].isPeriodic() && ! getPntrToArgument(i)->isPeriodic() ) {
2446 0 : error("in hills file periodicity for variable " + tmpvalues[i].getName() + " does not match periodicity in input");
2447 4436 : } else if( tmpvalues[i].isPeriodic() ) {
2448 : std::string imin, imax;
2449 0 : tmpvalues[i].getDomain( imin, imax );
2450 : std::string rmin, rmax;
2451 0 : getPntrToArgument(i)->getDomain( rmin, rmax );
2452 0 : if( imin!=rmin || imax!=rmax ) {
2453 0 : error("in hills file periodicity for variable " + tmpvalues[i].getName() + " does not match periodicity in input");
2454 : }
2455 : }
2456 4436 : center[i]=tmpvalues[i].get();
2457 : }
2458 : // scan for kerneltype
2459 2241 : std::string ktype="stretched-gaussian";
2460 4482 : if( ifile->FieldExist("kerneltype") ) {
2461 4464 : ifile->scanField("kerneltype",ktype);
2462 : }
2463 2241 : if( ktype=="gaussian" ) {
2464 12 : noStretchWarning();
2465 2229 : } else if( ktype!="stretched-gaussian") {
2466 0 : error("non Gaussian kernels are not supported in MetaD");
2467 : }
2468 : // scan for multivariate label: record the actual file position so to eventually rewind
2469 : std::string sss;
2470 4482 : ifile->scanField("multivariate",sss);
2471 2241 : if(sss=="true") {
2472 0 : multivariate=true;
2473 2241 : } else if(sss=="false") {
2474 2241 : multivariate=false;
2475 : } else {
2476 0 : plumed_merror("cannot parse multivariate = "+ sss);
2477 : }
2478 2241 : if(multivariate) {
2479 0 : sigma.resize(ncv*(ncv+1)/2);
2480 : Matrix<double> upper(ncv,ncv);
2481 : Matrix<double> lower(ncv,ncv);
2482 0 : for(unsigned i=0; i<ncv; i++) {
2483 0 : for(unsigned j=0; j<ncv-i; j++) {
2484 0 : ifile->scanField("sigma_"+getPntrToArgument(j+i)->getName()+"_"+getPntrToArgument(j)->getName(),lower(j+i,j));
2485 0 : upper(j,j+i)=lower(j+i,j);
2486 : }
2487 : }
2488 : Matrix<double> mymult(ncv,ncv);
2489 : Matrix<double> invmatrix(ncv,ncv);
2490 0 : mult(lower,upper,mymult);
2491 : // now invert and get the sigmas
2492 0 : Invert(mymult,invmatrix);
2493 : // put the sigmas in the usual order: upper diagonal (this time in normal form and not in band form)
2494 : unsigned k=0;
2495 0 : for(unsigned i=0; i<ncv; i++) {
2496 0 : for(unsigned j=i; j<ncv; j++) {
2497 0 : sigma[k]=invmatrix(i,j);
2498 0 : k++;
2499 : }
2500 : }
2501 : } else {
2502 6677 : for(unsigned i=0; i<ncv; ++i) {
2503 8872 : ifile->scanField("sigma_"+getPntrToArgument(i)->getName(),sigma[i]);
2504 : }
2505 : }
2506 :
2507 2241 : ifile->scanField("height",height);
2508 2241 : ifile->scanField("biasf",dummy);
2509 4482 : if(ifile->FieldExist("clock")) {
2510 4116 : ifile->scanField("clock",dummy);
2511 : }
2512 4482 : if(ifile->FieldExist("lower_int")) {
2513 0 : ifile->scanField("lower_int",dummy);
2514 : }
2515 4482 : if(ifile->FieldExist("upper_int")) {
2516 0 : ifile->scanField("upper_int",dummy);
2517 : }
2518 2241 : ifile->scanField();
2519 : return true;
2520 : } else {
2521 : return false;
2522 : }
2523 : }
2524 :
2525 102 : void MetaD::computeReweightingFactor() {
2526 102 : if(biasf_==1.0) { // in this case we have no bias, so reweight factor is 0.0
2527 0 : getPntrToComponent("rct")->set(0.0);
2528 0 : return;
2529 : }
2530 :
2531 102 : double Z_0=0; //proportional to the integral of exp(-beta*F)
2532 102 : double Z_V=0; //proportional to the integral of exp(-beta*(F+V))
2533 102 : double minusBetaF=biasf_/(biasf_-1.)/kbt_;
2534 102 : double minusBetaFplusV=1./(biasf_-1.)/kbt_;
2535 102 : if (biasf_==-1.0) { //non well-tempered case
2536 0 : minusBetaF=1./kbt_;
2537 : minusBetaFplusV=0;
2538 : }
2539 102 : max_bias_=BiasGrid_->getMaxValue(); //to avoid exp overflow
2540 :
2541 102 : const unsigned rank=comm.Get_rank();
2542 102 : const unsigned stride=comm.Get_size();
2543 920504 : for (Grid::index_t t=rank; t<BiasGrid_->getSize(); t+=stride) {
2544 920402 : const double val=BiasGrid_->getValue(t);
2545 920402 : Z_0+=std::exp(minusBetaF*(val-max_bias_));
2546 920402 : Z_V+=std::exp(minusBetaFplusV*(val-max_bias_));
2547 : }
2548 102 : comm.Sum(Z_0);
2549 102 : comm.Sum(Z_V);
2550 :
2551 102 : reweight_factor_=kbt_*std::log(Z_0/Z_V)+max_bias_;
2552 204 : getPntrToComponent("rct")->set(reweight_factor_);
2553 : }
2554 :
2555 273 : double MetaD::getTransitionBarrierBias() {
2556 : // If there is only one well of interest, return the bias at that well point.
2557 273 : if (transitionwells_.size() == 1) {
2558 0 : double tb_bias = getBias(transitionwells_[0]);
2559 0 : return tb_bias;
2560 :
2561 : // Otherwise, check for the least barrier bias between all pairs of wells.
2562 : // Note that because the paths can be considered edges between the wells' nodes
2563 : // to make a graph and the path barriers satisfy certain cycle inequalities, it
2564 : // is sufficient to look at paths corresponding to a minimal spanning tree of the
2565 : // overall graph rather than examining every edge in the graph.
2566 : // For simplicity, I chose the star graph with center well 0 as the spanning tree.
2567 : // It is most efficient to start the path searches from the wells that are
2568 : // expected to be sampled last, so transitionwell_[0] should correspond to the
2569 : // starting well. With this choice the searches will terminate in one step until
2570 : // transitionwell_[1] is sampled.
2571 : } else {
2572 : double least_transition_bias;
2573 273 : std::vector<double> sink = transitionwells_[0];
2574 273 : std::vector<double> source = transitionwells_[1];
2575 273 : least_transition_bias = BiasGrid_->findMaximalPathMinimum(source, sink);
2576 273 : for (unsigned i = 2; i < transitionwells_.size(); i++) {
2577 0 : if (least_transition_bias == 0.0) {
2578 : break;
2579 : }
2580 0 : source = transitionwells_[i];
2581 0 : double curr_transition_bias = BiasGrid_->findMaximalPathMinimum(source, sink);
2582 0 : least_transition_bias = fmin(curr_transition_bias, least_transition_bias);
2583 : }
2584 : return least_transition_bias;
2585 : }
2586 : }
2587 :
2588 154 : void MetaD::updateFrequencyAdaptiveStride() {
2589 154 : plumed_massert(freq_adaptive_,"should only be used if frequency adaptive metadynamics is enabled");
2590 154 : plumed_massert(acceleration_,"frequency adaptive metadynamics can only be used if the acceleration factor is calculated");
2591 154 : const double mean_acc = acc_/((double) getStep());
2592 154 : int tmp_stride= stride_*floor((mean_acc/fa_min_acceleration_)+0.5);
2593 154 : if(mean_acc >= fa_min_acceleration_) {
2594 129 : if(tmp_stride > current_stride_) {
2595 6 : current_stride_ = tmp_stride;
2596 : }
2597 : }
2598 154 : if(fa_max_stride_!=0 && current_stride_>fa_max_stride_) {
2599 0 : current_stride_=fa_max_stride_;
2600 : }
2601 154 : getPntrToComponent("pace")->set(current_stride_);
2602 154 : }
2603 :
2604 8435 : bool MetaD::checkNeedsGradients()const {
2605 8435 : if(adaptive_==FlexibleBin::geometry) {
2606 192 : if(getStep()%stride_==0 && !isFirstStep_) {
2607 : return true;
2608 : } else {
2609 109 : return false;
2610 : }
2611 : } else {
2612 : return false;
2613 : }
2614 : }
2615 :
2616 4 : void MetaD::updateNlist() {
2617 : // no need to check for neighbors
2618 4 : if(hills_.size()==0) {
2619 0 : return;
2620 : }
2621 :
2622 : // here we generate the neighbor list
2623 4 : nlist_hills_.clear();
2624 : std::vector<Gaussian> local_flat_nl;
2625 4 : unsigned nt=OpenMP::getNumThreads();
2626 4 : if(hills_.size()<2*nt) {
2627 : nt=1;
2628 : }
2629 4 : #pragma omp parallel num_threads(nt)
2630 : {
2631 : std::vector<Gaussian> private_flat_nl;
2632 : #pragma omp for nowait
2633 : for(unsigned k=0; k<hills_.size(); k++) {
2634 : double dist2=0;
2635 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
2636 : const double d=difference(i,getArgument(i),hills_[k].center[i])/hills_[k].sigma[i];
2637 : dist2+=d*d;
2638 : }
2639 : if(dist2<=nlist_param_[0]*dp2cutoff) {
2640 : private_flat_nl.push_back(hills_[k]);
2641 : }
2642 : }
2643 : #pragma omp critical
2644 : local_flat_nl.insert(local_flat_nl.end(), private_flat_nl.begin(), private_flat_nl.end());
2645 : }
2646 4 : nlist_hills_ = local_flat_nl;
2647 :
2648 : // here we set some properties that are used to decide when to update it again
2649 12 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
2650 8 : nlist_center_[i]=getArgument(i);
2651 : }
2652 : std::vector<double> dev2;
2653 4 : dev2.resize(getNumberOfArguments(),0);
2654 46 : for(unsigned k=0; k<nlist_hills_.size(); k++) {
2655 126 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
2656 84 : const double d=difference(i,getArgument(i),nlist_hills_[k].center[i]);
2657 84 : dev2[i]+=d*d;
2658 : }
2659 : }
2660 12 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
2661 8 : if(dev2[i]>0.) {
2662 8 : nlist_dev2_[i]=dev2[i]/static_cast<double>(nlist_hills_.size());
2663 : } else {
2664 0 : nlist_dev2_[i]=hills_.back().sigma[i]*hills_.back().sigma[i];
2665 : }
2666 : }
2667 :
2668 : // we are done
2669 4 : getPntrToComponent("nlker")->set(nlist_hills_.size());
2670 4 : getPntrToComponent("nlsteps")->set(nlist_steps_);
2671 4 : nlist_steps_=0;
2672 4 : nlist_update_=false;
2673 4 : }
2674 :
2675 : }
2676 : }
|