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