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