Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2015-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 PBMETAD
41 : /*
42 : Used to performed Parallel Bias metadynamics.
43 :
44 : This action activate Parallel Bias Metadynamics (PBMetaD) \cite pbmetad, a version of metadynamics \cite metad in which
45 : multiple low-dimensional bias potentials are applied in parallel.
46 : In the current implementation, these have the form of mono-dimensional metadynamics bias
47 : potentials:
48 :
49 : \f[
50 : {V(s_1,t), ..., V(s_N,t)}
51 : \f]
52 :
53 : where:
54 :
55 : \f[
56 : V(s_i,t) = \sum_{ k \tau < t} W_i(k \tau)
57 : \exp\left(
58 : - \frac{(s_i-s_i^{(0)}(k \tau))^2}{2\sigma_i^2}
59 : \right).
60 : \f]
61 :
62 : To ensure the convergence of each mono-dimensional bias potential to the corresponding free energy,
63 : at each deposition step the Gaussian heights are multiplied by the so-called conditional term:
64 :
65 : \f[
66 : W_i(k \tau)=W_0 \frac{\exp\left(
67 : - \frac{V(s_i,k \tau)}{k_B T}
68 : \right)}{\sum_{i=1}^N
69 : \exp\left(
70 : - \frac{V(s_i,k \tau)}{k_B T}
71 : \right)}
72 : \f]
73 :
74 : where \f$W_0\f$ is the initial Gaussian height.
75 :
76 : The PBMetaD bias potential is defined by:
77 :
78 : \f[
79 : V_{PB}(\vec{s},t) = -k_B T \log{\sum_{i=1}^N
80 : \exp\left(
81 : - \frac{V(s_i,t)}{k_B T}
82 : \right)}.
83 : \f]
84 :
85 :
86 : Information on the Gaussian functions that build each bias potential are printed to
87 : multiple HILLS files, which
88 : are used both to restart the calculation and to reconstruct the mono-dimensional
89 : free energies as a function of the corresponding CVs.
90 : These can be reconstructed using the \ref sum_hills utility because the final bias is given
91 : by:
92 :
93 : \f[
94 : V(s_i) = -F(s_i)
95 : \f]
96 :
97 : Currently, only a subset of the \ref METAD options are available in PBMetaD.
98 :
99 : The bias potentials can be stored on a grid to increase performances of long PBMetaD simulations.
100 : You should
101 : provide either the number of bins for every collective variable (GRID_BIN) or
102 : the desired grid spacing (GRID_SPACING). In case you provide both PLUMED will use
103 : the most conservative choice (highest number of bins) for each dimension.
104 : In case you do not provide any information about bin size (neither GRID_BIN nor GRID_SPACING)
105 : and if Gaussian width is fixed PLUMED will use 1/5 of the Gaussian width as grid spacing.
106 : This default choice should be reasonable for most applications.
107 :
108 : Another option that is available is well-tempered metadynamics \cite Barducci:2008. In this
109 : variant of PBMetaD the heights of the Gaussian hills are scaled at each step by the
110 : additional well-tempered metadynamics term.
111 : This ensures that each bias converges more smoothly. It should be noted that, in the case of well-tempered metadynamics, in
112 : the output printed the Gaussian height is re-scaled using the bias factor.
113 : Also notice that with well-tempered metadynamics the HILLS files do not contain the bias,
114 : but the negative of the free-energy estimate. This choice has the advantage that
115 : one can restart a simulation using a different value for the \f$\Delta T\f$. The applied bias will be scaled accordingly.
116 :
117 : Note that you can use here also the flexible Gaussian approach \cite Branduardi:2012dl
118 : in which you can adapt the Gaussian to the extent of Cartesian space covered by a variable or
119 : to the space in collective variable covered in a given time. In this case the width of the deposited
120 : Gaussian potential is denoted by one value only that is a Cartesian space (ADAPTIVE=GEOM) or a time
121 : (ADAPTIVE=DIFF). Note that a specific integration technique for the deposited Gaussian kernels
122 : should be used in this case. Check the documentation for utility sum_hills.
123 :
124 : With the keyword INTERVAL one changes the metadynamics algorithm setting the bias force equal to zero
125 : outside boundary \cite baftizadeh2012protein. If, for example, metadynamics is performed on a CV s and one is interested only
126 : to the free energy for s > boundary, the history dependent potential is still updated according to the above
127 : equations but the metadynamics force is set to zero for s < boundary. Notice that Gaussian kernels are added also
128 : if s < boundary, as the tails of these Gaussian kernels influence VG in the relevant region s > boundary. In this way, the
129 : force on the system in the region s > boundary comes from both metadynamics and the force field, in the region
130 : s < boundary only from the latter. This approach allows obtaining a history-dependent bias potential VG that
131 : fluctuates around a stable estimator, equal to the negative of the free energy far enough from the
132 : boundaries. Note that:
133 : - It works only for one-dimensional biases;
134 : - It works both with and without GRID;
135 : - The interval limit boundary in a region where the free energy derivative is not large;
136 : - If in the region outside the limit boundary the system has a free energy minimum, the INTERVAL keyword should
137 : be used together with a \ref UPPER_WALLS or \ref LOWER_WALLS at boundary.
138 :
139 : For systems with multiple CVs that share identical properties, PBMetaD with partitioned families can be used
140 : to group them under one bias potential that each contributes to \cite Prakash2018PF. This is done with a list
141 : of PF keywords, where each PF* argument contains the list of CVs from ARG to be placed in that family. Once
142 : invoked, each CV in ARG must be placed in exactly one PF, even if it results in families containing only one CV.
143 : Additionally, in cases where each of SIGMA or GRID entry would correspond to each ARG entry, they now correspond to
144 : each PF and must be adjusted accordingly.
145 :
146 : Multiple walkers \cite multiplewalkers can also be used. See below the examples.
147 :
148 : \par Examples
149 :
150 : The following input is for PBMetaD calculation using as
151 : collective variables the distance between atoms 3 and 5
152 : and the distance between atoms 2 and 4. The value of the CVs and
153 : the PBMetaD bias potential are written to the COLVAR file every 100 steps.
154 : \plumedfile
155 : DISTANCE ATOMS=3,5 LABEL=d1
156 : DISTANCE ATOMS=2,4 LABEL=d2
157 : PBMETAD ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3 PACE=500 LABEL=pb FILE=HILLS_d1,HILLS_d2
158 : PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
159 : \endplumedfile
160 : (See also \ref DISTANCE and \ref PRINT).
161 :
162 : \par
163 : If you use well-tempered metadynamics, you should specify a single bias factor and initial
164 : Gaussian height.
165 : \plumedfile
166 : DISTANCE ATOMS=3,5 LABEL=d1
167 : DISTANCE ATOMS=2,4 LABEL=d2
168 : PBMETAD ...
169 : ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
170 : PACE=500 BIASFACTOR=8 LABEL=pb
171 : FILE=HILLS_d1,HILLS_d2
172 : ... PBMETAD
173 : PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
174 : \endplumedfile
175 :
176 : \par
177 : Using partitioned families, each CV in ARG must be in exactly one family. Here,
178 : the distance between atoms 1,2 is degenerate with 2,4, but not with the
179 : distance between 3,5. Note that two SIGMA are provided to match the two PF.
180 : \plumedfile
181 : DISTANCE ATOMS=3,5 LABEL=d1
182 : DISTANCE ATOMS=2,4 LABEL=d2
183 : DISTANCE ATOMS=1,2 LABEL=d3
184 : PBMETAD ...
185 : ARG=d1,d2,d3 SIGMA=0.2,0.2 HEIGHT=0.3
186 : PF0=d1 PF1=d2,d3
187 : PACE=500 BIASFACTOR=8 LABEL=pb
188 : FILE=HILLS_d1,HILLS_d2
189 : ... PBMETAD
190 : PRINT ARG=d1,d2,d3,pb.bias STRIDE=100 FILE=COLVAR
191 : \endplumedfile
192 :
193 : \par
194 : The following input enables the MPI version of multiple-walkers.
195 : \plumedfile
196 : DISTANCE ATOMS=3,5 LABEL=d1
197 : DISTANCE ATOMS=2,4 LABEL=d2
198 : PBMETAD ...
199 : ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
200 : PACE=500 BIASFACTOR=8 LABEL=pb
201 : FILE=HILLS_d1,HILLS_d2
202 : WALKERS_MPI
203 : ... PBMETAD
204 : PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
205 : \endplumedfile
206 :
207 : \par
208 : The disk version of multiple-walkers can be
209 : enabled by setting the number of walker used, the id of the
210 : current walker which interprets the input file, the directory where the
211 : hills containing files resides, and the frequency to read the other walkers.
212 : Here is an example
213 : \plumedfile
214 : DISTANCE ATOMS=3,5 LABEL=d1
215 : DISTANCE ATOMS=2,4 LABEL=d2
216 : PBMETAD ...
217 : ARG=d1,d2 SIGMA=0.2,0.2 HEIGHT=0.3
218 : PACE=500 BIASFACTOR=8 LABEL=pb
219 : FILE=HILLS_d1,HILLS_d2
220 : WALKERS_N=10
221 : WALKERS_ID=3
222 : WALKERS_DIR=../
223 : WALKERS_RSTRIDE=100
224 : ... PBMETAD
225 : PRINT ARG=d1,d2,pb.bias STRIDE=100 FILE=COLVAR
226 : \endplumedfile
227 : where WALKERS_N is the total number of walkers, WALKERS_ID is the
228 : id of the present walker (starting from 0 ) and the WALKERS_DIR is the directory
229 : where all the walkers are located. WALKERS_RSTRIDE is the number of step between
230 : one update and the other.
231 :
232 : */
233 : //+ENDPLUMEDOC
234 :
235 : class PBMetaD : public Bias {
236 :
237 : private:
238 : struct Gaussian {
239 : std::vector<double> center;
240 : std::vector<double> sigma;
241 : double height;
242 : bool multivariate; // this is required to discriminate the one dimensional case
243 : std::vector<double> invsigma;
244 1120 : Gaussian(const std::vector<double> & center,const std::vector<double> & sigma, double height, bool multivariate):
245 1120 : center(center),sigma(sigma),height(height),multivariate(multivariate),invsigma(sigma) {
246 : // to avoid troubles from zero element in flexible hills
247 2240 : for(unsigned i=0; i<invsigma.size(); ++i)
248 1120 : if(std::abs(invsigma[i])>1.e-20) {
249 1120 : invsigma[i]=1.0/invsigma[i] ;
250 : } else {
251 0 : invsigma[i]=0.0;
252 : }
253 1120 : }
254 : };
255 : // general setup
256 : double kbt_;
257 : int stride_;
258 : // well-tempered MetaD
259 : bool welltemp_;
260 : double biasf_;
261 : // output files format
262 : std::string fmt_;
263 : // first step?
264 : bool isFirstStep_;
265 : // Gaussian starting parameters
266 : double height0_;
267 : std::vector<double> sigma0_;
268 : std::vector<double> sigma0min_;
269 : std::vector<double> sigma0max_;
270 : // Gaussians
271 : std::vector<std::vector<Gaussian> > hills_;
272 : std::vector<FlexibleBin> flexbin_;
273 : int adaptive_;
274 : std::vector<std::unique_ptr<OFile>> hillsOfiles_;
275 : std::vector<std::unique_ptr<IFile>> ifiles_;
276 : std::vector<std::string> ifilesnames_;
277 : // Grids
278 : bool grid_;
279 : std::vector<std::unique_ptr<GridBase>> BiasGrids_;
280 : std::vector<std::unique_ptr<OFile>> gridfiles_;
281 : int wgridstride_;
282 : // Partitioned Families
283 : unsigned int pf_n_; // initialize number of partitioned families
284 : std::vector<int> pfs_; //std::vector length of arguments that holds which pf# each cv belongs in
285 : std::vector<Value*> pfhold_; // std::vector length of pf_n which stores a pointer to the first argument fed to each family
286 : bool do_pf_; // if partitioned families are enabled
287 : // multiple walkers
288 : int mw_n_;
289 : std::string mw_dir_;
290 : int mw_id_;
291 : int mw_rstride_;
292 : bool walkers_mpi_;
293 : size_t mpi_nw_;
294 : unsigned mpi_id_;
295 : std::vector<std::string> hillsfname_;
296 : // intervals
297 : std::vector<double> uppI_;
298 : std::vector<double> lowI_;
299 : std::vector<bool> doInt_;
300 : // variable for selector
301 : std::string selector_;
302 : bool do_select_;
303 : unsigned select_value_;
304 : unsigned current_value_;
305 :
306 : double stretchA=1.0;
307 : double stretchB=0.0;
308 :
309 : bool noStretchWarningDone=false;
310 :
311 0 : void noStretchWarning() {
312 0 : if(!noStretchWarningDone) {
313 0 : log<<"\nWARNING: you are using a HILLS file with Gaussian kernels, PLUMED 2.8 uses stretched Gaussians by default\n";
314 : }
315 0 : noStretchWarningDone=true;
316 0 : }
317 :
318 : void readGaussians(unsigned iarg, IFile*);
319 : void writeGaussian(unsigned iarg, const Gaussian&, OFile*);
320 : void addGaussian(unsigned iarg, const Gaussian&);
321 : double getBiasAndDerivatives(unsigned iarg, const std::vector<double>&, double* der=NULL);
322 : double evaluateGaussian(unsigned iarg, const std::vector<double>&, const Gaussian&,double* der=NULL);
323 : std::vector<unsigned> getGaussianSupport(unsigned iarg, const Gaussian&);
324 : bool scanOneHill(unsigned iarg, IFile *ifile, std::vector<Value> &v, std::vector<double> ¢er, std::vector<double> &sigma, double &height, bool &multivariate);
325 :
326 : public:
327 : explicit PBMetaD(const ActionOptions&);
328 : void calculate() override;
329 : void update() override;
330 : static void registerKeywords(Keywords& keys);
331 : bool checkNeedsGradients()const override;
332 : };
333 :
334 : PLUMED_REGISTER_ACTION(PBMetaD,"PBMETAD")
335 :
336 44 : void PBMetaD::registerKeywords(Keywords& keys) {
337 44 : Bias::registerKeywords(keys);
338 44 : keys.add("compulsory","SIGMA","the widths of the Gaussian hills");
339 44 : keys.add("compulsory","PACE","the frequency for hill addition, one for all biases");
340 44 : keys.add("optional","FILE","files in which the lists of added hills are stored, default names are assigned using arguments if FILE is not found");
341 44 : keys.add("optional","HEIGHT","the height of the Gaussian hills, one for all biases. Compulsory unless TAU, TEMP and BIASFACTOR are given");
342 44 : keys.add("optional","FMT","specify format for HILLS files (useful for decrease the number of digits in regtests)");
343 44 : keys.add("optional","BIASFACTOR","use well tempered metadynamics with this bias factor, one for all biases. Please note you must also specify temp");
344 44 : keys.add("optional","TEMP","the system temperature - this is only needed if you are doing well-tempered metadynamics");
345 44 : keys.add("optional","TAU","in well tempered metadynamics, sets height to (k_B Delta T*pace*timestep)/tau");
346 44 : keys.add("optional","GRID_MIN","the lower bounds for the grid");
347 44 : keys.add("optional","GRID_MAX","the upper bounds for the grid");
348 44 : keys.add("optional","GRID_BIN","the number of bins for the grid");
349 44 : keys.add("optional","GRID_SPACING","the approximate grid spacing (to be used as an alternative or together with GRID_BIN)");
350 44 : keys.addFlag("GRID_SPARSE",false,"use a sparse grid to store hills");
351 44 : keys.addFlag("GRID_NOSPLINE",false,"don't use spline interpolation with grids");
352 44 : keys.add("optional","GRID_WSTRIDE", "frequency for dumping the grid");
353 44 : keys.add("optional","GRID_WFILES", "dump grid for the bias, default names are used if GRID_WSTRIDE is used without GRID_WFILES.");
354 44 : keys.add("optional","GRID_RFILES", "read grid for the bias");
355 44 : keys.add("optional","ADAPTIVE","use a geometric (=GEOM) or diffusion (=DIFF) based hills width scheme. Sigma is one number that has distance units or timestep dimensions");
356 44 : keys.add("optional","SIGMA_MAX","the upper bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
357 44 : keys.add("optional","SIGMA_MIN","the lower bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
358 88 : keys.addInputKeyword("numbered","PF", "scalar", "specify which CVs belong in a partitioned family. Once a PF is specified, all CVs in ARG must be placed in a PF even if there is one CV per PF”");
359 44 : keys.add("optional","SELECTOR", "add forces and do update based on the value of SELECTOR");
360 44 : keys.add("optional","SELECTOR_ID", "value of SELECTOR");
361 44 : keys.add("optional","WALKERS_ID", "walker id");
362 44 : keys.add("optional","WALKERS_N", "number of walkers");
363 44 : keys.add("optional","WALKERS_DIR", "shared directory with the hills files from all the walkers");
364 44 : keys.add("optional","WALKERS_RSTRIDE","stride for reading hills files");
365 44 : keys.addFlag("WALKERS_MPI",false,"Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR");
366 44 : keys.add("optional","INTERVAL_MIN","one dimensional lower limits, outside the limits the system will not feel the biasing force.");
367 44 : keys.add("optional","INTERVAL_MAX","one dimensional upper limits, outside the limits the system will not feel the biasing force.");
368 44 : keys.use("RESTART");
369 44 : keys.use("UPDATE_FROM");
370 44 : keys.use("UPDATE_UNTIL");
371 44 : }
372 :
373 42 : PBMetaD::PBMetaD(const ActionOptions& ao):
374 : PLUMED_BIAS_INIT(ao),
375 42 : kbt_(0.0),
376 42 : stride_(0),
377 42 : welltemp_(false),
378 42 : biasf_(1.0),
379 42 : isFirstStep_(true),
380 42 : height0_(std::numeric_limits<double>::max()),
381 42 : adaptive_(FlexibleBin::none),
382 42 : grid_(false),
383 42 : wgridstride_(0),
384 42 : pf_n_(0), do_pf_(false),
385 42 : mw_n_(1), mw_dir_(""), mw_id_(0), mw_rstride_(1),
386 42 : walkers_mpi_(false), mpi_nw_(0),
387 42 : do_select_(false) {
388 :
389 : // parse the flexible hills
390 : std::string adaptiveoption;
391 : adaptiveoption="NONE";
392 84 : parse("ADAPTIVE",adaptiveoption);
393 42 : if(adaptiveoption=="GEOM") {
394 0 : log.printf(" Uses Geometry-based hills width: sigma must be in distance units and only one sigma is needed\n");
395 0 : adaptive_=FlexibleBin::geometry;
396 42 : } else if(adaptiveoption=="DIFF") {
397 4 : log.printf(" Uses Diffusion-based hills width: sigma must be in time steps and only one sigma is needed\n");
398 4 : adaptive_=FlexibleBin::diffusion;
399 38 : } else if(adaptiveoption=="NONE") {
400 38 : adaptive_=FlexibleBin::none;
401 : } else {
402 0 : error("I do not know this type of adaptive scheme");
403 : }
404 :
405 42 : parse("FMT",fmt_);
406 :
407 : // Partitioned Families - fill with -1 to mark as invalid
408 42 : pfs_.assign(getNumberOfArguments(), -1);
409 42 : pfhold_.resize(getNumberOfArguments());
410 : std::vector<Value*> familyargs;
411 42 : for(int i = 0;; i++) {
412 100 : parseArgumentList("PF", i, familyargs);
413 50 : if (familyargs.empty()) {
414 : break;
415 : }
416 :
417 8 : do_pf_ = true;
418 8 : log << " Identified Partitioned Family " << i << ":";
419 20 : for (unsigned j = 0; j < familyargs.size(); j++) {
420 12 : log << " " << familyargs[j]->getName();
421 : // loop through the argument list to make sure it exists and assign it
422 : bool foundArg = false;
423 48 : for (unsigned argnum = 0; argnum < getNumberOfArguments(); argnum++) {
424 36 : if (familyargs[j]->getName() == getPntrToArgument(argnum)->getName()) {
425 : foundArg = true;
426 12 : if (pfs_[argnum] != -1) {
427 0 : error(familyargs[j]->getName() + " already present in PF" + std::to_string(pfs_[argnum]));
428 : }
429 12 : pfs_[argnum] = i; // store the pf# for each cv
430 12 : if (pfhold_[i] == nullptr) {
431 : // if this is the first argument in the family, store a pointer for it (this is for HILLS & GRID files)
432 8 : pfhold_[i] = getPntrToArgument(argnum);
433 : }
434 : }
435 : }
436 12 : if (!foundArg) {
437 0 : error(familyargs[j]->getName() + " in PF" + std::to_string(i) + " not found in ARG");
438 : }
439 : }
440 8 : log << "\n";
441 8 : pf_n_++;
442 8 : }
443 :
444 : // if PF were specified, every argument gets treated as its own PF
445 42 : if (!do_pf_) {
446 38 : pf_n_ = getNumberOfArguments();
447 114 : for(unsigned i=0; i < pf_n_; i++) {
448 76 : pfhold_[i] = getPntrToArgument(i);
449 76 : pfs_[i] = i;
450 : }
451 : } else {
452 : // If we are doing PF, make sure each argument got assigned to a family.
453 16 : for (unsigned i = 0; i < getNumberOfArguments(); i++) {
454 12 : if (pfs_[i] == -1) {
455 0 : error(getPntrToArgument(i)->getName() + " was not assigned a PF");
456 : }
457 : }
458 : }
459 :
460 : // parse the sigma
461 42 : parseVector("SIGMA",sigma0_);
462 42 : if(adaptive_==FlexibleBin::none) {
463 : // if you use normal sigma you need one sigma per argument
464 38 : if( sigma0_.size()!=pf_n_ ) {
465 0 : std::string fields = do_pf_ ? "PFs" : "arguments";
466 0 : error("number of " + fields + " does not match number of SIGMA parameters");
467 : }
468 : } else {
469 : // if you use flexible hills you need one sigma
470 4 : if(sigma0_.size()!=1) {
471 0 : error("If you choose ADAPTIVE you need only one sigma according to your choice of type (GEOM/DIFF)");
472 : }
473 : // if adaptive then the number must be an integer
474 4 : if(adaptive_==FlexibleBin::diffusion) {
475 4 : if(int(sigma0_[0])-sigma0_[0]>1.e-9 || int(sigma0_[0])-sigma0_[0] <-1.e-9 || int(sigma0_[0])<1 ) {
476 0 : error("In case of adaptive hills with diffusion, the sigma must be an integer which is the number of time steps\n");
477 : }
478 : }
479 : // here evtl parse the sigma min and max values
480 8 : parseVector("SIGMA_MIN",sigma0min_);
481 4 : if(sigma0min_.size()>0 && sigma0min_.size()!=pf_n_) {
482 0 : error("the number of SIGMA_MIN values be the same of the number of the arguments/PF");
483 4 : } else if(sigma0min_.size()==0) {
484 0 : sigma0min_.resize(pf_n_);
485 0 : for(unsigned i=0; i<pf_n_; i++) {
486 0 : sigma0min_[i]=-1.;
487 : }
488 : }
489 :
490 8 : parseVector("SIGMA_MAX",sigma0max_);
491 4 : if(sigma0max_.size()>0 && sigma0max_.size()!=pf_n_) {
492 0 : error("the number of SIGMA_MAX values be the same of the number of the arguments/PF");
493 4 : } else if(sigma0max_.size()==0) {
494 4 : sigma0max_.resize(pf_n_);
495 12 : for(unsigned i=0; i<pf_n_; i++) {
496 8 : sigma0max_[i]=-1.;
497 : }
498 : }
499 :
500 12 : for(unsigned i=0; i<pf_n_; i++) {
501 : std::vector<double> tmp_smin, tmp_smax;
502 8 : tmp_smin.resize(1,sigma0min_[i]);
503 8 : tmp_smax.resize(1,sigma0max_[i]);
504 16 : flexbin_.push_back(FlexibleBin(adaptive_,this,i,sigma0_[0],tmp_smin,tmp_smax));
505 : }
506 : }
507 :
508 : // note: HEIGHT is not compulsory, since one could use the TAU keyword, see below
509 42 : parse("HEIGHT",height0_);
510 42 : parse("PACE",stride_);
511 42 : if(stride_<=0) {
512 0 : error("frequency for hill addition is nonsensical");
513 : }
514 :
515 :
516 84 : parseVector("FILE",hillsfname_);
517 42 : if(hillsfname_.size()==0) {
518 30 : for(unsigned i=0; i< pf_n_; i++) {
519 20 : std::string name = do_pf_ ? "HILLS.PF"+std::to_string(i) : "HILLS."+getPntrToArgument(i)->getName();
520 20 : hillsfname_.push_back(name);
521 : }
522 : }
523 42 : if( hillsfname_.size()!=pf_n_ ) {
524 0 : error("number of FILE arguments does not match number of HILLS files");
525 : }
526 :
527 42 : parse("BIASFACTOR",biasf_);
528 42 : if( biasf_<1.0 ) {
529 0 : error("well tempered bias factor is nonsensical");
530 : }
531 42 : kbt_=getkBT();
532 42 : if(biasf_>1.0) {
533 41 : if(kbt_==0.0) {
534 0 : error("Unless the MD engine passes the temperature to plumed, with well-tempered metad you must specify it using TEMP");
535 : }
536 41 : welltemp_=true;
537 : }
538 42 : double tau=0.0;
539 42 : parse("TAU",tau);
540 42 : if(tau==0.0) {
541 42 : if(height0_==std::numeric_limits<double>::max()) {
542 0 : error("At least one between HEIGHT and TAU should be specified");
543 : }
544 : // if tau is not set, we compute it here from the other input parameters
545 42 : if(welltemp_) {
546 41 : tau=(kbt_*(biasf_-1.0))/height0_*getTimeStep()*stride_;
547 : }
548 : } else {
549 0 : if(!welltemp_) {
550 0 : error("TAU only makes sense in well-tempered metadynamics");
551 : }
552 0 : if(height0_!=std::numeric_limits<double>::max()) {
553 0 : error("At most one between HEIGHT and TAU should be specified");
554 : }
555 0 : height0_=(kbt_*(biasf_-1.0))/tau*getTimeStep()*stride_;
556 : }
557 :
558 :
559 : // Multiple walkers
560 42 : parse("WALKERS_N",mw_n_);
561 42 : parse("WALKERS_ID",mw_id_);
562 42 : if(mw_n_<=mw_id_) {
563 0 : error("walker ID should be a numerical value less than the total number of walkers");
564 : }
565 42 : parse("WALKERS_DIR",mw_dir_);
566 42 : parse("WALKERS_RSTRIDE",mw_rstride_);
567 :
568 : // MPI version
569 42 : parseFlag("WALKERS_MPI",walkers_mpi_);
570 :
571 : // Grid file
572 84 : parse("GRID_WSTRIDE",wgridstride_);
573 : std::vector<std::string> gridfilenames_;
574 42 : parseVector("GRID_WFILES",gridfilenames_);
575 42 : if (wgridstride_ == 0 && gridfilenames_.size() > 0) {
576 0 : error("frequency with which to output grid not specified use GRID_WSTRIDE");
577 : }
578 42 : if(gridfilenames_.size()==0 && wgridstride_ > 0) {
579 12 : for(unsigned i=0; i<pf_n_; i++) {
580 8 : std::string name = do_pf_ ? "GRID.PF"+std::to_string(i) : "GRID."+getPntrToArgument(i)->getName();
581 8 : gridfilenames_.push_back(name);
582 : }
583 : }
584 42 : if(gridfilenames_.size() > 0 && hillsfname_.size() > 0 && gridfilenames_.size() != hillsfname_.size()) {
585 0 : error("number of GRID_WFILES arguments does not match number of HILLS files");
586 : }
587 :
588 : // Read grid
589 : std::vector<std::string> gridreadfilenames_;
590 42 : parseVector("GRID_RFILES",gridreadfilenames_);
591 :
592 : // Grid Stuff
593 42 : std::vector<std::string> gmin(pf_n_);
594 84 : parseVector("GRID_MIN",gmin);
595 42 : if(gmin.size()!=pf_n_ && gmin.size()!=0) {
596 0 : error("not enough values for GRID_MIN");
597 : }
598 42 : std::vector<std::string> gmax(pf_n_);
599 84 : parseVector("GRID_MAX",gmax);
600 42 : if(gmax.size()!=pf_n_ && gmax.size()!=0) {
601 0 : error("not enough values for GRID_MAX");
602 : }
603 42 : std::vector<unsigned> gbin(pf_n_);
604 : std::vector<double> gspacing;
605 84 : parseVector("GRID_BIN",gbin);
606 42 : if(gbin.size()!=pf_n_ && gbin.size()!=0) {
607 0 : error("not enough values for GRID_BIN");
608 : }
609 84 : parseVector("GRID_SPACING",gspacing);
610 42 : if(gspacing.size()!=pf_n_ && gspacing.size()!=0) {
611 0 : error("not enough values for GRID_SPACING");
612 : }
613 42 : if(gmin.size()!=gmax.size()) {
614 0 : error("GRID_MAX and GRID_MIN should be either present or absent");
615 : }
616 42 : if(gspacing.size()!=0 && gmin.size()==0) {
617 0 : error("If GRID_SPACING is present also GRID_MIN and GRID_MAX should be present");
618 : }
619 42 : if(gbin.size()!=0 && gmin.size()==0) {
620 0 : error("If GRID_BIN is present also GRID_MIN and GRID_MAX should be present");
621 : }
622 42 : if(gmin.size()!=0) {
623 10 : if(gbin.size()==0 && gspacing.size()==0) {
624 10 : if(adaptive_==FlexibleBin::none) {
625 6 : log<<" Binsize not specified, 1/5 of sigma will be be used\n";
626 6 : plumed_assert(sigma0_.size()==pf_n_);
627 6 : gspacing.resize(pf_n_);
628 18 : for(unsigned i=0; i<gspacing.size(); i++) {
629 12 : gspacing[i]=0.2*sigma0_[i];
630 : }
631 : } else {
632 : // with adaptive hills and grid a sigma min must be specified
633 12 : for(unsigned i=0; i<sigma0min_.size(); i++)
634 8 : if(sigma0min_[i]<=0) {
635 0 : error("When using ADAPTIVE Gaussians on a grid SIGMA_MIN must be specified");
636 : }
637 4 : log<<" Binsize not specified, 1/5 of sigma_min will be be used\n";
638 4 : gspacing.resize(pf_n_);
639 12 : for(unsigned i=0; i<gspacing.size(); i++) {
640 8 : gspacing[i]=0.2*sigma0min_[i];
641 : }
642 : }
643 0 : } else if(gspacing.size()!=0 && gbin.size()==0) {
644 0 : log<<" The number of bins will be estimated from GRID_SPACING\n";
645 0 : } else if(gspacing.size()!=0 && gbin.size()!=0) {
646 0 : log<<" You specified both GRID_BIN and GRID_SPACING\n";
647 0 : log<<" The more conservative (highest) number of bins will be used for each variable\n";
648 : }
649 10 : if(gbin.size()==0) {
650 10 : gbin.assign(pf_n_,1);
651 : }
652 10 : if(gspacing.size()!=0)
653 30 : for(unsigned i=0; i<pf_n_; i++) {
654 : double a,b;
655 20 : Tools::convert(gmin[i],a);
656 20 : Tools::convert(gmax[i],b);
657 20 : unsigned n=std::ceil(((b-a)/gspacing[i]));
658 20 : if(gbin[i]<n) {
659 20 : gbin[i]=n;
660 : }
661 : }
662 : }
663 42 : if(gbin.size()>0) {
664 10 : grid_=true;
665 : }
666 :
667 42 : bool sparsegrid=false;
668 42 : parseFlag("GRID_SPARSE",sparsegrid);
669 42 : bool nospline=false;
670 42 : parseFlag("GRID_NOSPLINE",nospline);
671 42 : bool spline=!nospline;
672 42 : if(!grid_&&gridfilenames_.size() > 0) {
673 0 : error("To write a grid you need first to define it!");
674 : }
675 42 : if(!grid_&&gridreadfilenames_.size() > 0) {
676 0 : error("To read a grid you need first to define it!");
677 : }
678 :
679 42 : doInt_.resize(pf_n_,false);
680 : // Interval keyword
681 42 : parseVector("INTERVAL_MIN",lowI_);
682 84 : parseVector("INTERVAL_MAX",uppI_);
683 : // various checks
684 42 : if(lowI_.size()!=uppI_.size()) {
685 0 : error("both a lower and an upper limits must be provided with INTERVAL");
686 : }
687 42 : if(lowI_.size()!=0 && lowI_.size()!=pf_n_) {
688 0 : error("check number of argument of INTERVAL");
689 : }
690 50 : for(unsigned i=0; i<lowI_.size(); ++i) {
691 8 : if(uppI_[i]<lowI_[i]) {
692 0 : error("The Upper limit must be greater than the Lower limit!");
693 : }
694 8 : if(pfhold_[i]->isPeriodic()) {
695 0 : warning("INTERVAL is not used for periodic variables");
696 : } else {
697 : doInt_[i]=true;
698 : }
699 : }
700 :
701 : // parse selector stuff
702 84 : parse("SELECTOR", selector_);
703 42 : if(selector_.length()>0) {
704 1 : do_select_ = true;
705 1 : select_value_ = 0; // set defalt value or it might be not initialized if the user does not pass SELECTOR_ID
706 2 : parse("SELECTOR_ID", select_value_);
707 : }
708 :
709 42 : checkRead();
710 :
711 42 : log.printf(" Gaussian width ");
712 42 : if (adaptive_==FlexibleBin::diffusion) {
713 4 : log.printf(" (Note: The units of sigma are in timesteps) ");
714 : }
715 42 : if (adaptive_==FlexibleBin::geometry) {
716 0 : log.printf(" (Note: The units of sigma are in dist units) ");
717 : }
718 122 : for(unsigned i=0; i<sigma0_.size(); ++i) {
719 80 : log.printf(" %f",sigma0_[i]);
720 : }
721 42 : log.printf(" Gaussian height %f\n",height0_);
722 42 : log.printf(" Gaussian deposition pace %d\n",stride_);
723 42 : log.printf(" Gaussian files ");
724 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
725 84 : log.printf("%s ",hillsfname_[i].c_str());
726 : }
727 42 : log.printf("\n");
728 42 : if(welltemp_) {
729 41 : log.printf(" Well-Tempered Bias Factor %f\n",biasf_);
730 41 : log.printf(" Hills relaxation time (tau) %f\n",tau);
731 41 : log.printf(" KbT %f\n",kbt_);
732 : }
733 :
734 42 : if(do_select_) {
735 1 : log.printf(" Add forces and update bias based on the value of SELECTOR %s\n",selector_.c_str());
736 1 : log.printf(" Id of the SELECTOR for this action %u\n", select_value_);
737 : }
738 :
739 42 : if(mw_n_>1) {
740 0 : if(walkers_mpi_) {
741 0 : error("MPI version of multiple walkers is not compatible with filesystem version of multiple walkers");
742 : }
743 0 : log.printf(" %d multiple walkers active\n",mw_n_);
744 0 : log.printf(" walker id %d\n",mw_id_);
745 0 : log.printf(" reading stride %d\n",mw_rstride_);
746 0 : if(mw_dir_!="") {
747 0 : log.printf(" directory with hills files %s\n",mw_dir_.c_str());
748 : }
749 : } else {
750 42 : if(walkers_mpi_) {
751 34 : log.printf(" Multiple walkers active using MPI communnication\n");
752 34 : if(mw_dir_!="") {
753 0 : log.printf(" directory with hills files %s\n",mw_dir_.c_str());
754 : }
755 34 : if(comm.Get_rank()==0) {
756 : // Only root of group can communicate with other walkers
757 18 : mpi_nw_ = multi_sim_comm.Get_size();
758 18 : mpi_id_ = multi_sim_comm.Get_rank();
759 : }
760 : // Communicate to the other members of the same group
761 : // info abount number of walkers and walker index
762 34 : comm.Bcast(mpi_nw_,0);
763 34 : comm.Bcast(mpi_id_,0);
764 : }
765 : }
766 :
767 126 : for(unsigned i=0; i<doInt_.size(); i++) {
768 84 : if(doInt_[i]) {
769 8 : log.printf(" Upper and Lower limits boundaries for the bias of CV %u are activated\n", i);
770 : }
771 : }
772 42 : if(grid_) {
773 10 : log.printf(" Grid min");
774 30 : for(unsigned i=0; i<gmin.size(); ++i) {
775 20 : log.printf(" %s",gmin[i].c_str() );
776 : }
777 10 : log.printf("\n");
778 10 : log.printf(" Grid max");
779 30 : for(unsigned i=0; i<gmax.size(); ++i) {
780 20 : log.printf(" %s",gmax[i].c_str() );
781 : }
782 10 : log.printf("\n");
783 10 : log.printf(" Grid bin");
784 30 : for(unsigned i=0; i<gbin.size(); ++i) {
785 20 : log.printf(" %u",gbin[i]);
786 : }
787 10 : log.printf("\n");
788 10 : if(spline) {
789 10 : log.printf(" Grid uses spline interpolation\n");
790 : }
791 10 : if(sparsegrid) {
792 0 : log.printf(" Grid uses sparse grid\n");
793 : }
794 10 : if(wgridstride_>0) {
795 18 : for(unsigned i=0; i<gridfilenames_.size(); ++i) {
796 12 : log.printf(" Grid is written on file %s with stride %d\n",gridfilenames_[i].c_str(),wgridstride_);
797 : }
798 : }
799 10 : if(gridreadfilenames_.size()>0) {
800 3 : for(unsigned i=0; i<gridreadfilenames_.size(); ++i) {
801 2 : log.printf(" Reading bias from grid in file %s \n",gridreadfilenames_[i].c_str());
802 : }
803 : }
804 : }
805 :
806 : // initializing vector of hills
807 42 : hills_.resize(pf_n_);
808 :
809 : // restart from external grid
810 : bool restartedFromGrid=false;
811 :
812 : // initializing and checking grid
813 42 : if(grid_) {
814 : // check for mesh and sigma size
815 30 : for(unsigned i=0; i<pf_n_; i++) {
816 : double a,b;
817 20 : int family = pfs_[i]; // point to families instead of arguments
818 20 : Tools::convert(gmin[family],a);
819 20 : Tools::convert(gmax[family],b);
820 20 : double mesh=(b-a)/((double)gbin[family]);
821 20 : if(adaptive_==FlexibleBin::none) {
822 12 : if(mesh>0.5*sigma0_[i]) {
823 0 : log<<" WARNING: Using a PBMETAD with a Grid Spacing larger than half of the Gaussians width can produce artifacts\n";
824 : }
825 : } else {
826 8 : if(mesh>0.5*sigma0min_[i]||sigma0min_[i]<0.) {
827 0 : log<<" WARNING: to use a PBMETAD with a GRID and ADAPTIVE you need to set a Grid Spacing larger than half of the Gaussians \n";
828 : }
829 : }
830 : }
831 10 : std::string funcl=getLabel() + ".bias";
832 30 : for(unsigned i=0; i<pf_n_; ++i) {
833 20 : std::vector<Value*> args(1);
834 20 : args[0] = pfhold_[i]; //Use first argument in family for interactions.
835 20 : std::vector<std::string> gmin_t(1);
836 20 : std::vector<std::string> gmax_t(1);
837 20 : std::vector<unsigned> gbin_t(1);
838 : gmin_t[0] = gmin[i];
839 : gmax_t[0] = gmax[i];
840 20 : gbin_t[0] = gbin[i];
841 20 : std::unique_ptr<GridBase> BiasGrid_;
842 : // Read grid from file
843 20 : if(gridreadfilenames_.size()>0) {
844 2 : IFile gridfile;
845 2 : gridfile.link(*this);
846 2 : if(gridfile.FileExist(gridreadfilenames_[i])) {
847 2 : gridfile.open(gridreadfilenames_[i]);
848 : } else {
849 0 : error("The GRID file you want to read: " + gridreadfilenames_[i] + ", cannot be found!");
850 : }
851 2 : std::string funcl = getLabel() + ".bias";
852 4 : BiasGrid_=GridBase::create(funcl, args, gridfile, gmin_t, gmax_t, gbin_t, sparsegrid, spline, true);
853 2 : if(BiasGrid_->getDimension() != args.size()) {
854 0 : error("mismatch between dimensionality of input grid and number of arguments");
855 : }
856 4 : if(pfhold_[i]->isPeriodic() != BiasGrid_->getIsPeriodic()[0]) {
857 0 : error("periodicity mismatch between arguments and input bias");
858 : }
859 2 : log.printf(" Restarting from %s:\n",gridreadfilenames_[i].c_str());
860 2 : if(getRestart()) {
861 : restartedFromGrid=true;
862 : }
863 2 : } else {
864 18 : if(!sparsegrid) {
865 18 : BiasGrid_=Tools::make_unique<Grid>(funcl,args,gmin_t,gmax_t,gbin_t,spline,true);
866 : } else {
867 0 : BiasGrid_=Tools::make_unique<SparseGrid>(funcl,args,gmin_t,gmax_t,gbin_t,spline,true);
868 : }
869 18 : std::vector<std::string> actualmin=BiasGrid_->getMin();
870 18 : std::vector<std::string> actualmax=BiasGrid_->getMax();
871 : std::string is;
872 18 : Tools::convert(i,is);
873 18 : if(gmin_t[0]!=actualmin[0]) {
874 0 : error("GRID_MIN["+is+"] must be adjusted to "+actualmin[0]+" to fit periodicity");
875 : }
876 18 : if(gmax_t[0]!=actualmax[0]) {
877 0 : error("GRID_MAX["+is+"] must be adjusted to "+actualmax[0]+" to fit periodicity");
878 : }
879 18 : }
880 20 : BiasGrids_.emplace_back(std::move(BiasGrid_));
881 40 : }
882 : }
883 :
884 :
885 :
886 : // creating vector of ifile* for hills reading
887 : // open all files at the beginning and read Gaussians if restarting
888 :
889 84 : for(int j=0; j<mw_n_; ++j) {
890 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
891 84 : unsigned k=j*hillsfname_.size()+i;
892 : std::string fname;
893 84 : if(mw_dir_!="") {
894 0 : if(mw_n_>1) {
895 0 : std::stringstream out;
896 0 : out << j;
897 0 : fname = mw_dir_+"/"+hillsfname_[i]+"."+out.str();
898 0 : } else if(walkers_mpi_) {
899 0 : fname = mw_dir_+"/"+hillsfname_[i];
900 : } else {
901 : fname = hillsfname_[i];
902 : }
903 : } else {
904 84 : if(mw_n_>1) {
905 0 : std::stringstream out;
906 0 : out << j;
907 0 : fname = hillsfname_[i]+"."+out.str();
908 0 : } else {
909 : fname = hillsfname_[i];
910 : }
911 : }
912 84 : ifiles_.emplace_back(Tools::make_unique<IFile>());
913 : // this is just a shortcut pointer to the last element:
914 : IFile *ifile = ifiles_.back().get();
915 84 : ifile->link(*this);
916 84 : ifilesnames_.push_back(fname);
917 84 : if(ifile->FileExist(fname)) {
918 18 : ifile->open(fname);
919 18 : if(getRestart()&&!restartedFromGrid) {
920 2 : log.printf(" Restarting from %s:",ifilesnames_[k].c_str());
921 2 : readGaussians(i,ifiles_[k].get());
922 : }
923 18 : ifiles_[k]->reset(false);
924 : // close only the walker own hills file for later writing
925 18 : if(j==mw_id_) {
926 18 : ifiles_[k]->close();
927 : }
928 : } else {
929 : // in case a file does not exist and we are restarting, complain that the file was not found
930 66 : if(getRestart()) {
931 2 : log<<" WARNING: restart file "<<fname<<" not found\n";
932 : }
933 : }
934 : }
935 : }
936 :
937 42 : comm.Barrier();
938 42 : if(comm.Get_rank()==0 && walkers_mpi_) {
939 18 : multi_sim_comm.Barrier();
940 : }
941 :
942 : // open hills files for writing
943 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
944 : auto ofile=Tools::make_unique<OFile>();
945 84 : ofile->link(*this);
946 : // if MPI multiple walkers, only rank 0 will write to file
947 84 : if(walkers_mpi_) {
948 68 : int r=0;
949 68 : if(comm.Get_rank()==0) {
950 36 : r=multi_sim_comm.Get_rank();
951 : }
952 68 : comm.Bcast(r,0);
953 68 : if(r>0) {
954 34 : ifilesnames_[mw_id_*hillsfname_.size()+i]="/dev/null";
955 : }
956 136 : ofile->enforceSuffix("");
957 : }
958 84 : if(mw_n_>1) {
959 0 : ofile->enforceSuffix("");
960 : }
961 84 : ofile->open(ifilesnames_[mw_id_*hillsfname_.size()+i]);
962 84 : if(fmt_.length()>0) {
963 24 : ofile->fmtField(fmt_);
964 : }
965 168 : ofile->addConstantField("multivariate");
966 168 : ofile->addConstantField("kerneltype");
967 84 : if(doInt_[i]) {
968 16 : ofile->addConstantField("lower_int").printField("lower_int",lowI_[i]);
969 16 : ofile->addConstantField("upper_int").printField("upper_int",uppI_[i]);
970 : }
971 84 : ofile->setHeavyFlush();
972 : // output periodicities of variables
973 84 : ofile->setupPrintValue( pfhold_[i] ); //assuming cvs in the same family have the same periodicity and boundaries.
974 : // push back
975 84 : hillsOfiles_.emplace_back(std::move(ofile));
976 84 : }
977 :
978 : // Dump grid to files
979 42 : if(wgridstride_ > 0) {
980 18 : for(unsigned i = 0; i < gridfilenames_.size(); ++i) {
981 : auto ofile=Tools::make_unique<OFile>();
982 12 : ofile->link(*this);
983 12 : std::string gridfname_tmp = gridfilenames_[i];
984 12 : if(walkers_mpi_) {
985 8 : int r = 0;
986 8 : if(comm.Get_rank() == 0) {
987 4 : r = multi_sim_comm.Get_rank();
988 : }
989 8 : comm.Bcast(r, 0);
990 8 : if(r>0) {
991 : gridfname_tmp = "/dev/null";
992 : }
993 16 : ofile->enforceSuffix("");
994 : }
995 12 : if(mw_n_>1) {
996 0 : ofile->enforceSuffix("");
997 : }
998 12 : ofile->open(gridfname_tmp);
999 12 : ofile->setHeavyFlush();
1000 12 : gridfiles_.emplace_back(std::move(ofile));
1001 12 : }
1002 : }
1003 :
1004 84 : log<<" Bibliography "<<plumed.cite("Pfaendtner and Bonomi. J. Chem. Theory Comput. 11, 5062 (2015)");
1005 42 : if(doInt_[0])
1006 8 : log<<plumed.cite(
1007 8 : "Baftizadeh, Cossio, Pietrucci, and Laio, Curr. Phys. Chem. 2, 79 (2012)");
1008 42 : if(mw_n_>1||walkers_mpi_)
1009 68 : log<<plumed.cite(
1010 68 : "Raiteri, Laio, Gervasio, Micheletti, and Parrinello, J. Phys. Chem. B 110, 3533 (2006)");
1011 42 : if(adaptive_!=FlexibleBin::none)
1012 8 : log<<plumed.cite(
1013 8 : "Branduardi, Bussi, and Parrinello, J. Chem. Theory Comput. 8, 2247 (2012)");
1014 42 : if (do_pf_) {
1015 8 : log<<plumed.cite("Prakash, Fu, Bonomi, and Pfaendtner, J. Chem. Theory Comput. 14, 4985 (2018)");
1016 : }
1017 42 : log<<"\n";
1018 :
1019 :
1020 84 : }
1021 :
1022 2 : void PBMetaD::readGaussians(unsigned iarg, IFile *ifile) {
1023 2 : std::vector<double> center(1);
1024 2 : std::vector<double> sigma(1);
1025 : double height;
1026 : int nhills=0;
1027 2 : bool multivariate=false;
1028 2 : int family=pfs_[iarg];
1029 :
1030 : std::vector<Value> tmpvalues;
1031 4 : tmpvalues.push_back( Value( this, pfhold_[family]->getName(), false ) );
1032 :
1033 10 : while(scanOneHill(iarg,ifile,tmpvalues,center,sigma,height,multivariate)) {
1034 : ;
1035 8 : nhills++;
1036 8 : if(welltemp_) {
1037 8 : height*=(biasf_-1.0)/biasf_;
1038 : }
1039 8 : addGaussian(family, Gaussian(center,sigma,height,multivariate));
1040 : }
1041 2 : log.printf(" %d Gaussians read\n",nhills);
1042 4 : }
1043 :
1044 1112 : void PBMetaD::writeGaussian(unsigned iarg, const Gaussian& hill, OFile *ofile) {
1045 1112 : int family=pfs_[iarg];
1046 2224 : ofile->printField("time",getTimeStep()*getStep());
1047 1112 : ofile->printField(pfhold_[family],hill.center[0]);
1048 :
1049 2224 : ofile->printField("kerneltype","stretched-gaussian");
1050 1112 : if(hill.multivariate) {
1051 288 : ofile->printField("multivariate","true");
1052 144 : double lower = std::sqrt(1./hill.sigma[0]);
1053 288 : ofile->printField("sigma_"+pfhold_[family]->getName()+"_"+
1054 144 : pfhold_[family]->getName(),lower);
1055 : } else {
1056 1936 : ofile->printField("multivariate","false");
1057 1936 : ofile->printField("sigma_"+pfhold_[family]->getName(),hill.sigma[0]);
1058 : }
1059 1112 : double height=hill.height;
1060 1112 : if(welltemp_) {
1061 1104 : height *= biasf_/(biasf_-1.0);
1062 : }
1063 1112 : ofile->printField("height",height);
1064 1112 : ofile->printField("biasf",biasf_);
1065 1112 : if(mw_n_>1) {
1066 0 : ofile->printField("clock",int(std::time(0)));
1067 : }
1068 1112 : ofile->printField();
1069 1112 : }
1070 :
1071 1120 : void PBMetaD::addGaussian(unsigned iarg, const Gaussian& hill) {
1072 1120 : if(!grid_) {
1073 912 : hills_[iarg].push_back(hill);
1074 : } else {
1075 208 : std::vector<unsigned> nneighb=getGaussianSupport(iarg, hill);
1076 208 : std::vector<Grid::index_t> neighbors=BiasGrids_[iarg]->getNeighbors(hill.center,nneighb);
1077 208 : std::vector<double> der(1);
1078 208 : std::vector<double> xx(1);
1079 208 : if(comm.Get_size()==1) {
1080 1520 : for(unsigned i=0; i<neighbors.size(); ++i) {
1081 1480 : Grid::index_t ineigh=neighbors[i];
1082 1480 : der[0]=0.0;
1083 1480 : BiasGrids_[iarg]->getPoint(ineigh,xx);
1084 1480 : double bias=evaluateGaussian(iarg,xx,hill,&der[0]);
1085 1480 : BiasGrids_[iarg]->addValueAndDerivatives(ineigh,bias,der);
1086 : }
1087 : } else {
1088 168 : unsigned stride=comm.Get_size();
1089 168 : unsigned rank=comm.Get_rank();
1090 168 : std::vector<double> allder(neighbors.size(),0.0);
1091 168 : std::vector<double> allbias(neighbors.size(),0.0);
1092 3276 : for(unsigned i=rank; i<neighbors.size(); i+=stride) {
1093 3108 : Grid::index_t ineigh=neighbors[i];
1094 3108 : BiasGrids_[iarg]->getPoint(ineigh,xx);
1095 3108 : allbias[i]=evaluateGaussian(iarg,xx,hill,&allder[i]);
1096 : }
1097 168 : comm.Sum(allbias);
1098 168 : comm.Sum(allder);
1099 6384 : for(unsigned i=0; i<neighbors.size(); ++i) {
1100 6216 : Grid::index_t ineigh=neighbors[i];
1101 6216 : der[0]=allder[i];
1102 6216 : BiasGrids_[iarg]->addValueAndDerivatives(ineigh,allbias[i],der);
1103 : }
1104 : }
1105 : }
1106 1120 : }
1107 :
1108 208 : std::vector<unsigned> PBMetaD::getGaussianSupport(unsigned iarg, const Gaussian& hill) {
1109 : std::vector<unsigned> nneigh;
1110 : double cutoff;
1111 208 : if(hill.multivariate) {
1112 144 : double maxautoval=1./hill.sigma[0];
1113 144 : cutoff=std::sqrt(2.0*dp2cutoff*maxautoval);
1114 : } else {
1115 64 : cutoff=std::sqrt(2.0*dp2cutoff)*hill.sigma[0];
1116 : }
1117 :
1118 208 : if(doInt_[iarg]) {
1119 144 : if(hill.center[0]+cutoff > uppI_[iarg] || hill.center[0]-cutoff < lowI_[iarg]) {
1120 : // in this case, we updated the entire grid to avoid problems
1121 0 : return BiasGrids_[iarg]->getNbin();
1122 : } else {
1123 144 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff/BiasGrids_[iarg]->getDx()[0])));
1124 : return nneigh;
1125 : }
1126 : }
1127 :
1128 64 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff/BiasGrids_[iarg]->getDx()[0])) );
1129 :
1130 : return nneigh;
1131 : }
1132 :
1133 1296 : double PBMetaD::getBiasAndDerivatives(unsigned iarg, const std::vector<double>& cv, double* der) {
1134 1296 : double bias=0.0;
1135 1296 : int family = pfs_[iarg];
1136 1296 : if(!grid_) {
1137 1000 : unsigned stride=comm.Get_size();
1138 1000 : unsigned rank=comm.Get_rank();
1139 5520 : for(unsigned i=rank; i<hills_[family].size(); i+=stride) {
1140 4520 : bias += evaluateGaussian(iarg,cv,hills_[family][i],der);
1141 : }
1142 1000 : comm.Sum(bias);
1143 1000 : if(der) {
1144 540 : comm.Sum(der,1);
1145 : }
1146 : } else {
1147 296 : if(der) {
1148 160 : std::vector<double> vder(1);
1149 160 : bias = BiasGrids_[family]->getValueAndDerivatives(cv,vder);
1150 160 : der[0] = vder[0];
1151 : } else {
1152 136 : bias = BiasGrids_[family]->getValue(cv);
1153 : }
1154 : }
1155 :
1156 1296 : return bias;
1157 : }
1158 :
1159 9108 : double PBMetaD::evaluateGaussian(unsigned iarg, const std::vector<double>& cv, const Gaussian& hill, double* der) {
1160 : double bias=0.0;
1161 : // I use a pointer here because cv is const (and should be const)
1162 : // but when using doInt it is easier to locally replace cv[0] with
1163 : // the upper/lower limit in case it is out of range
1164 : const double *pcv=NULL;
1165 : double tmpcv[1]; // tmp array with cv (to be used with doInt_)
1166 9108 : tmpcv[0]=cv[0];
1167 : bool isOutOfInt = false;
1168 9108 : if(doInt_[iarg]) {
1169 2664 : if(cv[0]<lowI_[iarg]) {
1170 : tmpcv[0]=lowI_[iarg];
1171 : isOutOfInt = true;
1172 2664 : } else if(cv[0]>uppI_[iarg]) {
1173 : tmpcv[0]=uppI_[iarg];
1174 : isOutOfInt = true;
1175 : }
1176 : }
1177 : pcv=&(tmpcv[0]);
1178 :
1179 9108 : if(hill.multivariate) {
1180 2664 : double dp = difference(iarg, hill.center[0], pcv[0]);
1181 2664 : double dp2 = 0.5 * dp * dp * hill.sigma[0];
1182 2664 : if(dp2<dp2cutoff) {
1183 2534 : bias = hill.height*std::exp(-dp2);
1184 2534 : if(der && !isOutOfInt) {
1185 2534 : der[0] += -bias * dp * hill.sigma[0] * stretchA;
1186 : }
1187 2534 : bias=stretchA*bias+hill.height*stretchB;
1188 : }
1189 : } else {
1190 6444 : double dp = difference(iarg, hill.center[0], pcv[0]) * hill.invsigma[0];
1191 6444 : double dp2 = 0.5 * dp * dp;
1192 6444 : if(dp2<dp2cutoff) {
1193 6363 : bias = hill.height*std::exp(-dp2);
1194 6363 : if(der && !isOutOfInt) {
1195 4107 : der[0] += -bias * dp * hill.invsigma[0] * stretchA;
1196 : }
1197 6363 : bias=stretchA*bias+hill.height*stretchB;
1198 : }
1199 : }
1200 :
1201 9108 : return bias;
1202 : }
1203 :
1204 340 : void PBMetaD::calculate() {
1205 : // this is because presently there is no way to properly pass information
1206 : // on adaptive hills (diff) after exchanges:
1207 340 : if(adaptive_==FlexibleBin::diffusion && getExchangeStep()) {
1208 0 : error("ADAPTIVE=DIFF is not compatible with replica exchange");
1209 : }
1210 :
1211 340 : std::vector<double> cv(1);
1212 : double der[1];
1213 340 : std::vector<double> bias(getNumberOfArguments());
1214 340 : std::vector<double> deriv(getNumberOfArguments());
1215 :
1216 340 : double ncv = (double) getNumberOfArguments();
1217 : double bmin = 1.0e+19;
1218 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1219 700 : cv[0] = getArgument(i);
1220 700 : der[0] = 0.0;
1221 700 : bias[i] = getBiasAndDerivatives(i, cv, der);
1222 700 : deriv[i] = der[0];
1223 700 : if(bias[i] < bmin) {
1224 : bmin = bias[i];
1225 : }
1226 : }
1227 : double ene = 0.;
1228 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1229 700 : ene += std::exp((-bias[i]+bmin)/kbt_);
1230 : }
1231 :
1232 : // set Forces - set them to zero if SELECTOR is active
1233 340 : if(do_select_) {
1234 5 : current_value_ = static_cast<unsigned>(plumed.passMap[selector_]);
1235 : }
1236 :
1237 340 : if(!do_select_ || select_value_==current_value_) {
1238 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1239 700 : const double f = - std::exp((-bias[i]+bmin)/kbt_) / (ene) * deriv[i];
1240 700 : setOutputForce(i, f);
1241 : }
1242 : }
1243 :
1244 340 : if(do_select_ && select_value_!=current_value_) {
1245 0 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1246 0 : setOutputForce(i, 0.0);
1247 : }
1248 : }
1249 :
1250 : // set bias
1251 340 : ene = -kbt_ * (std::log(ene) - std::log(ncv)) + bmin;
1252 340 : setBias(ene);
1253 340 : }
1254 :
1255 340 : void PBMetaD::update() {
1256 : bool multivariate;
1257 : // adding hills criteria
1258 : bool nowAddAHill;
1259 340 : if(getStep()%stride_==0 && !isFirstStep_) {
1260 : nowAddAHill=true;
1261 : } else {
1262 : nowAddAHill=false;
1263 50 : isFirstStep_=false;
1264 : }
1265 :
1266 : // if you use adaptive, call the FlexibleBin
1267 340 : if(adaptive_!=FlexibleBin::none) {
1268 120 : for(unsigned i=0; i<getNumberOfArguments(); i++) {
1269 80 : flexbin_[i].update(nowAddAHill,i);
1270 : }
1271 : multivariate=true;
1272 : } else {
1273 : multivariate=false;
1274 : }
1275 :
1276 340 : if(nowAddAHill && (!do_select_ || select_value_==current_value_)) {
1277 : // get all biases and heights
1278 290 : std::vector<double> cv(getNumberOfArguments());
1279 290 : std::vector<double> bias(getNumberOfArguments());
1280 290 : std::vector<double> thissigma(getNumberOfArguments());
1281 290 : std::vector<double> height(getNumberOfArguments());
1282 290 : std::vector<double> cv_tmp(1);
1283 290 : std::vector<double> sigma_tmp(1);
1284 : double norm = 0.0;
1285 : double bmin = 1.0e+19;
1286 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1287 596 : int family=pfs_[i];
1288 : // get flex/sigmas for each family and assign them to this args sigma
1289 596 : if(adaptive_!=FlexibleBin::none) {
1290 72 : thissigma[i]=flexbin_[family].getInverseMatrix(i)[0];
1291 : } else {
1292 524 : thissigma[i]=sigma0_[family];
1293 : }
1294 596 : cv[i] = getArgument(i);
1295 596 : cv_tmp[0] = getArgument(i);
1296 596 : bias[i] = getBiasAndDerivatives(i, cv_tmp);
1297 596 : if(bias[i] < bmin) {
1298 : bmin = bias[i];
1299 : }
1300 : }
1301 : // calculate heights and norm
1302 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1303 596 : double h = std::exp((-bias[i]+bmin)/kbt_);
1304 596 : norm += h;
1305 596 : height[i] = h;
1306 : }
1307 : // normalize and apply welltemp correction
1308 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1309 596 : height[i] *= height0_ / norm;
1310 596 : if(welltemp_) {
1311 588 : height[i] *= std::exp(-bias[i]/(kbt_*(biasf_-1.0)));
1312 : }
1313 : }
1314 :
1315 : // MPI Multiple walkers: share hills and add them all
1316 290 : if(walkers_mpi_) {
1317 : // Allocate arrays to store all walkers hills
1318 258 : std::vector<double> all_cv(mpi_nw_*cv.size(), 0.0);
1319 258 : std::vector<double> all_sigma(mpi_nw_*getNumberOfArguments(), 0.0);
1320 258 : std::vector<double> all_height(mpi_nw_*height.size(), 0.0);
1321 258 : if(comm.Get_rank()==0) {
1322 : // fill in value
1323 390 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1324 260 : unsigned j = mpi_id_ * getNumberOfArguments() + i;
1325 260 : all_cv[j] = cv[i];
1326 260 : all_sigma[j] = thissigma[i];
1327 260 : all_height[j] = height[i];
1328 : }
1329 : // Communicate (only root)
1330 130 : multi_sim_comm.Sum(&all_cv[0], all_cv.size());
1331 130 : multi_sim_comm.Sum(&all_sigma[0], all_sigma.size());
1332 130 : multi_sim_comm.Sum(&all_height[0], all_height.size());
1333 : }
1334 : // Share info with group members
1335 258 : comm.Sum(&all_cv[0], all_cv.size());
1336 258 : comm.Sum(&all_sigma[0], all_sigma.size());
1337 258 : comm.Sum(&all_height[0], all_height.size());
1338 : // now add hills one by one
1339 774 : for(unsigned j=0; j<mpi_nw_; ++j) {
1340 1548 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1341 : // Add CVs of same family together and write to same file
1342 1032 : int family = pfs_[i];
1343 1032 : cv_tmp[0] = all_cv[j*cv.size()+i];
1344 1032 : double height_tmp = all_height[j*cv.size()+i];
1345 1032 : sigma_tmp[0] = all_sigma[j*cv.size()+i];
1346 1032 : Gaussian newhill = Gaussian(cv_tmp, sigma_tmp, height_tmp, multivariate);
1347 1032 : addGaussian(family, newhill);
1348 1032 : writeGaussian(i, newhill, hillsOfiles_[family].get());
1349 1032 : }
1350 : }
1351 : // just add your own hills
1352 : } else {
1353 112 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1354 : // Add CVs of same family together and write to same file
1355 80 : int family = pfs_[i];
1356 80 : cv_tmp[0] = cv[i];
1357 80 : if(adaptive_!=FlexibleBin::none) {
1358 0 : sigma_tmp[0]=thissigma[i];
1359 : } else {
1360 80 : sigma_tmp[0] = sigma0_[family];
1361 : }
1362 80 : Gaussian newhill = Gaussian(cv_tmp, sigma_tmp, height[i], multivariate);
1363 80 : addGaussian(family, newhill);
1364 80 : writeGaussian(i, newhill, hillsOfiles_[family].get());
1365 80 : }
1366 : }
1367 : }
1368 :
1369 : // write grid files
1370 340 : if(wgridstride_>0 && (getStep()%wgridstride_==0 || getCPT())) {
1371 14 : int r = 0;
1372 14 : if(walkers_mpi_) {
1373 4 : if(comm.Get_rank()==0) {
1374 2 : r=multi_sim_comm.Get_rank();
1375 : }
1376 4 : comm.Bcast(r,0);
1377 : }
1378 14 : if(r==0) {
1379 36 : for(unsigned i=0; i<gridfiles_.size(); ++i) {
1380 24 : gridfiles_[i]->rewind();
1381 24 : BiasGrids_[i]->writeToFile(*gridfiles_[i]);
1382 24 : gridfiles_[i]->flush();
1383 : }
1384 : }
1385 : }
1386 :
1387 : // if multiple walkers and time to read Gaussians
1388 340 : if(mw_n_>1 && getStep()%mw_rstride_==0) {
1389 0 : for(int j=0; j<mw_n_; ++j) {
1390 0 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
1391 0 : unsigned k=j*hillsfname_.size()+i;
1392 : // don't read your own Gaussians
1393 0 : if(j==mw_id_) {
1394 0 : continue;
1395 : }
1396 : // if the file is not open yet
1397 0 : if(!(ifiles_[k]->isOpen())) {
1398 : // check if it exists now and open it!
1399 0 : if(ifiles_[k]->FileExist(ifilesnames_[k])) {
1400 0 : ifiles_[k]->open(ifilesnames_[k]);
1401 0 : ifiles_[k]->reset(false);
1402 : }
1403 : // otherwise read the new Gaussians
1404 : } else {
1405 0 : log.printf(" Reading hills from %s:",ifilesnames_[k].c_str());
1406 0 : readGaussians(i,ifiles_[k].get());
1407 0 : ifiles_[k]->reset(false);
1408 : }
1409 : }
1410 : }
1411 : }
1412 :
1413 340 : }
1414 :
1415 : /// takes a pointer to the file and a template string with values v and gives back the next center, sigma and height
1416 10 : bool PBMetaD::scanOneHill(unsigned iarg, IFile *ifile, std::vector<Value> &tmpvalues, std::vector<double> ¢er, std::vector<double> &sigma, double &height, bool &multivariate) {
1417 : double dummy;
1418 10 : multivariate=false;
1419 10 : Value* argPtr = pfhold_[pfs_[iarg]];
1420 20 : if(ifile->scanField("time",dummy)) {
1421 8 : ifile->scanField( &tmpvalues[0] );
1422 8 : if( tmpvalues[0].isPeriodic() && ! argPtr->isPeriodic() ) {
1423 0 : error("in hills file periodicity for variable " + tmpvalues[0].getName() + " does not match periodicity in input");
1424 8 : } else if( tmpvalues[0].isPeriodic() ) {
1425 : std::string imin, imax;
1426 0 : tmpvalues[0].getDomain( imin, imax );
1427 : std::string rmin, rmax;
1428 0 : argPtr->getDomain( rmin, rmax );
1429 0 : if( imin!=rmin || imax!=rmax ) {
1430 0 : error("in hills file periodicity for variable " + tmpvalues[0].getName() + " does not match periodicity in input");
1431 : }
1432 : }
1433 8 : center[0]=tmpvalues[0].get();
1434 8 : std::string ktype="stretched-gaussian";
1435 16 : if( ifile->FieldExist("kerneltype") ) {
1436 16 : ifile->scanField("kerneltype",ktype);
1437 : }
1438 :
1439 8 : if( ktype=="gaussian" ) {
1440 0 : noStretchWarning();
1441 8 : } else if( ktype!="stretched-gaussian") {
1442 0 : error("non Gaussian kernels are not supported in MetaD");
1443 : }
1444 :
1445 : std::string sss;
1446 16 : ifile->scanField("multivariate",sss);
1447 8 : if(sss=="true") {
1448 0 : multivariate=true;
1449 8 : } else if(sss=="false") {
1450 8 : multivariate=false;
1451 : } else {
1452 0 : plumed_merror("cannot parse multivariate = "+ sss);
1453 : }
1454 8 : if(multivariate) {
1455 0 : ifile->scanField("sigma_"+argPtr->getName()+"_"+
1456 : argPtr->getName(),sigma[0]);
1457 0 : sigma[0] = 1./(sigma[0]*sigma[0]);
1458 : } else {
1459 16 : ifile->scanField("sigma_"+argPtr->getName(),sigma[0]);
1460 : }
1461 8 : ifile->scanField("height",height);
1462 8 : ifile->scanField("biasf",dummy);
1463 16 : if(ifile->FieldExist("clock")) {
1464 0 : ifile->scanField("clock",dummy);
1465 : }
1466 16 : if(ifile->FieldExist("lower_int")) {
1467 0 : ifile->scanField("lower_int",dummy);
1468 : }
1469 16 : if(ifile->FieldExist("upper_int")) {
1470 0 : ifile->scanField("upper_int",dummy);
1471 : }
1472 8 : ifile->scanField();
1473 : return true;
1474 : } else {
1475 : return false;
1476 : }
1477 :
1478 : }
1479 :
1480 340 : bool PBMetaD::checkNeedsGradients()const {
1481 340 : if(adaptive_==FlexibleBin::geometry) {
1482 0 : if(getStep()%stride_==0 && !isFirstStep_) {
1483 : return true;
1484 : } else {
1485 0 : return false;
1486 : }
1487 : } else {
1488 : return false;
1489 : }
1490 : }
1491 :
1492 : }
1493 : }
1494 :
|