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) if(std::abs(invsigma[i])>1.e-20) invsigma[i]=1.0/invsigma[i] ; else invsigma[i]=0.0;
248 1120 : }
249 : };
250 : // general setup
251 : double kbt_;
252 : int stride_;
253 : // well-tempered MetaD
254 : bool welltemp_;
255 : double biasf_;
256 : // output files format
257 : std::string fmt_;
258 : // first step?
259 : bool isFirstStep_;
260 : // Gaussian starting parameters
261 : double height0_;
262 : std::vector<double> sigma0_;
263 : std::vector<double> sigma0min_;
264 : std::vector<double> sigma0max_;
265 : // Gaussians
266 : std::vector<std::vector<Gaussian> > hills_;
267 : std::vector<FlexibleBin> flexbin_;
268 : int adaptive_;
269 : std::vector<std::unique_ptr<OFile>> hillsOfiles_;
270 : std::vector<std::unique_ptr<IFile>> ifiles_;
271 : std::vector<std::string> ifilesnames_;
272 : // Grids
273 : bool grid_;
274 : std::vector<std::unique_ptr<GridBase>> BiasGrids_;
275 : std::vector<std::unique_ptr<OFile>> gridfiles_;
276 : int wgridstride_;
277 : // Partitioned Families
278 : unsigned int pf_n_; // initialize number of partitioned families
279 : std::vector<int> pfs_; //std::vector length of arguments that holds which pf# each cv belongs in
280 : std::vector<Value*> pfhold_; // std::vector length of pf_n which stores a pointer to the first argument fed to each family
281 : bool do_pf_; // if partitioned families are enabled
282 : // multiple walkers
283 : int mw_n_;
284 : std::string mw_dir_;
285 : int mw_id_;
286 : int mw_rstride_;
287 : bool walkers_mpi_;
288 : size_t mpi_nw_;
289 : unsigned mpi_id_;
290 : std::vector<std::string> hillsfname_;
291 : // intervals
292 : std::vector<double> uppI_;
293 : std::vector<double> lowI_;
294 : std::vector<bool> doInt_;
295 : // variable for selector
296 : std::string selector_;
297 : bool do_select_;
298 : unsigned select_value_;
299 : unsigned current_value_;
300 :
301 : double stretchA=1.0;
302 : double stretchB=0.0;
303 :
304 : bool noStretchWarningDone=false;
305 :
306 0 : void noStretchWarning() {
307 0 : if(!noStretchWarningDone) {
308 0 : log<<"\nWARNING: you are using a HILLS file with Gaussian kernels, PLUMED 2.8 uses stretched Gaussians by default\n";
309 : }
310 0 : noStretchWarningDone=true;
311 0 : }
312 :
313 : void readGaussians(unsigned iarg, IFile*);
314 : void writeGaussian(unsigned iarg, const Gaussian&, OFile*);
315 : void addGaussian(unsigned iarg, const Gaussian&);
316 : double getBiasAndDerivatives(unsigned iarg, const std::vector<double>&, double* der=NULL);
317 : double evaluateGaussian(unsigned iarg, const std::vector<double>&, const Gaussian&,double* der=NULL);
318 : std::vector<unsigned> getGaussianSupport(unsigned iarg, const Gaussian&);
319 : bool scanOneHill(unsigned iarg, IFile *ifile, std::vector<Value> &v, std::vector<double> ¢er, std::vector<double> &sigma, double &height, bool &multivariate);
320 :
321 : public:
322 : explicit PBMetaD(const ActionOptions&);
323 : void calculate() override;
324 : void update() override;
325 : static void registerKeywords(Keywords& keys);
326 : bool checkNeedsGradients()const override;
327 : };
328 :
329 : PLUMED_REGISTER_ACTION(PBMetaD,"PBMETAD")
330 :
331 44 : void PBMetaD::registerKeywords(Keywords& keys) {
332 44 : Bias::registerKeywords(keys);
333 44 : keys.use("ARG");
334 88 : keys.add("compulsory","SIGMA","the widths of the Gaussian hills");
335 88 : keys.add("compulsory","PACE","the frequency for hill addition, one for all biases");
336 88 : 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");
337 88 : keys.add("optional","HEIGHT","the height of the Gaussian hills, one for all biases. Compulsory unless TAU, TEMP and BIASFACTOR are given");
338 88 : keys.add("optional","FMT","specify format for HILLS files (useful for decrease the number of digits in regtests)");
339 88 : keys.add("optional","BIASFACTOR","use well tempered metadynamics with this bias factor, one for all biases. Please note you must also specify temp");
340 88 : keys.add("optional","TEMP","the system temperature - this is only needed if you are doing well-tempered metadynamics");
341 88 : keys.add("optional","TAU","in well tempered metadynamics, sets height to (k_B Delta T*pace*timestep)/tau");
342 88 : keys.add("optional","GRID_MIN","the lower bounds for the grid");
343 88 : keys.add("optional","GRID_MAX","the upper bounds for the grid");
344 88 : keys.add("optional","GRID_BIN","the number of bins for the grid");
345 88 : keys.add("optional","GRID_SPACING","the approximate grid spacing (to be used as an alternative or together with GRID_BIN)");
346 88 : keys.addFlag("GRID_SPARSE",false,"use a sparse grid to store hills");
347 88 : keys.addFlag("GRID_NOSPLINE",false,"don't use spline interpolation with grids");
348 88 : keys.add("optional","GRID_WSTRIDE", "frequency for dumping the grid");
349 88 : keys.add("optional","GRID_WFILES", "dump grid for the bias, default names are used if GRID_WSTRIDE is used without GRID_WFILES.");
350 88 : keys.add("optional","GRID_RFILES", "read grid for the bias");
351 88 : 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");
352 88 : keys.add("optional","SIGMA_MAX","the upper bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
353 88 : keys.add("optional","SIGMA_MIN","the lower bounds for the sigmas (in CV units) when using adaptive hills. Negative number means no bounds ");
354 88 : keys.add("numbered","PF", "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”");
355 88 : keys.add("optional","SELECTOR", "add forces and do update based on the value of SELECTOR");
356 88 : keys.add("optional","SELECTOR_ID", "value of SELECTOR");
357 88 : keys.add("optional","WALKERS_ID", "walker id");
358 88 : keys.add("optional","WALKERS_N", "number of walkers");
359 88 : keys.add("optional","WALKERS_DIR", "shared directory with the hills files from all the walkers");
360 88 : keys.add("optional","WALKERS_RSTRIDE","stride for reading hills files");
361 88 : keys.addFlag("WALKERS_MPI",false,"Switch on MPI version of multiple walkers - not compatible with WALKERS_* options other than WALKERS_DIR");
362 88 : keys.add("optional","INTERVAL_MIN","one dimensional lower limits, outside the limits the system will not feel the biasing force.");
363 88 : keys.add("optional","INTERVAL_MAX","one dimensional upper limits, outside the limits the system will not feel the biasing force.");
364 44 : keys.use("RESTART");
365 44 : keys.use("UPDATE_FROM");
366 44 : keys.use("UPDATE_UNTIL");
367 44 : }
368 :
369 42 : PBMetaD::PBMetaD(const ActionOptions& ao):
370 : PLUMED_BIAS_INIT(ao),
371 42 : kbt_(0.0),
372 42 : stride_(0),
373 42 : welltemp_(false),
374 42 : biasf_(1.0),
375 42 : isFirstStep_(true),
376 42 : height0_(std::numeric_limits<double>::max()),
377 42 : adaptive_(FlexibleBin::none),
378 42 : grid_(false),
379 42 : wgridstride_(0),
380 42 : pf_n_(0), do_pf_(false),
381 42 : mw_n_(1), mw_dir_(""), mw_id_(0), mw_rstride_(1),
382 42 : walkers_mpi_(false), mpi_nw_(0),
383 42 : do_select_(false)
384 : {
385 :
386 : // parse the flexible hills
387 : std::string adaptiveoption;
388 : adaptiveoption="NONE";
389 84 : parse("ADAPTIVE",adaptiveoption);
390 42 : if(adaptiveoption=="GEOM") {
391 0 : log.printf(" Uses Geometry-based hills width: sigma must be in distance units and only one sigma is needed\n");
392 0 : adaptive_=FlexibleBin::geometry;
393 42 : } else if(adaptiveoption=="DIFF") {
394 4 : log.printf(" Uses Diffusion-based hills width: sigma must be in time steps and only one sigma is needed\n");
395 4 : adaptive_=FlexibleBin::diffusion;
396 38 : } else if(adaptiveoption=="NONE") {
397 38 : adaptive_=FlexibleBin::none;
398 : } else {
399 0 : error("I do not know this type of adaptive scheme");
400 : }
401 :
402 42 : parse("FMT",fmt_);
403 :
404 : // Partitioned Families - fill with -1 to mark as invalid
405 42 : pfs_.assign(getNumberOfArguments(), -1);
406 42 : pfhold_.resize(getNumberOfArguments());
407 : std::vector<Value*> familyargs;
408 42 : for(int i = 0;; i++) {
409 100 : parseArgumentList("PF", i, familyargs);
410 50 : if (familyargs.empty()) break;
411 :
412 8 : do_pf_ = true;
413 8 : log << " Identified Partitioned Family " << i << ":";
414 20 : for (unsigned j = 0; j < familyargs.size(); j++) {
415 12 : log << " " << familyargs[j]->getName();
416 : // loop through the argument list to make sure it exists and assign it
417 : bool foundArg = false;
418 48 : for (unsigned argnum = 0; argnum < getNumberOfArguments(); argnum++) {
419 36 : if (familyargs[j]->getName() == getPntrToArgument(argnum)->getName()) {
420 : foundArg = true;
421 12 : if (pfs_[argnum] != -1) {
422 0 : error(familyargs[j]->getName() + " already present in PF" + std::to_string(pfs_[argnum]));
423 : }
424 12 : pfs_[argnum] = i; // store the pf# for each cv
425 12 : if (pfhold_[i] == nullptr) {
426 : // if this is the first argument in the family, store a pointer for it (this is for HILLS & GRID files)
427 8 : pfhold_[i] = getPntrToArgument(argnum);
428 : }
429 : }
430 : }
431 12 : if (!foundArg) {
432 0 : error(familyargs[j]->getName() + " in PF" + std::to_string(i) + " not found in ARG");
433 : }
434 : }
435 8 : log << "\n";
436 8 : pf_n_++;
437 8 : }
438 :
439 : // if PF were specified, every argument gets treated as its own PF
440 42 : if (!do_pf_) {
441 38 : pf_n_ = getNumberOfArguments();
442 114 : for(unsigned i=0; i < pf_n_; i++) {
443 76 : pfhold_[i] = getPntrToArgument(i);
444 76 : pfs_[i] = i;
445 : }
446 : } else {
447 : // If we are doing PF, make sure each argument got assigned to a family.
448 16 : for (unsigned i = 0; i < getNumberOfArguments(); i++) {
449 12 : if (pfs_[i] == -1) error(getPntrToArgument(i)->getName() + " was not assigned a PF");
450 : }
451 : }
452 :
453 : // parse the sigma
454 42 : parseVector("SIGMA",sigma0_);
455 42 : if(adaptive_==FlexibleBin::none) {
456 : // if you use normal sigma you need one sigma per argument
457 38 : if( sigma0_.size()!=pf_n_ ) {
458 0 : std::string fields = do_pf_ ? "PFs" : "arguments";
459 0 : error("number of " + fields + " does not match number of SIGMA parameters");
460 : }
461 : } else {
462 : // if you use flexible hills you need one sigma
463 4 : if(sigma0_.size()!=1) {
464 0 : error("If you choose ADAPTIVE you need only one sigma according to your choice of type (GEOM/DIFF)");
465 : }
466 : // if adaptive then the number must be an integer
467 4 : if(adaptive_==FlexibleBin::diffusion) {
468 4 : if(int(sigma0_[0])-sigma0_[0]>1.e-9 || int(sigma0_[0])-sigma0_[0] <-1.e-9 || int(sigma0_[0])<1 ) {
469 0 : error("In case of adaptive hills with diffusion, the sigma must be an integer which is the number of time steps\n");
470 : }
471 : }
472 : // here evtl parse the sigma min and max values
473 8 : parseVector("SIGMA_MIN",sigma0min_);
474 4 : if(sigma0min_.size()>0 && sigma0min_.size()!=pf_n_) {
475 0 : error("the number of SIGMA_MIN values be the same of the number of the arguments/PF");
476 4 : } else if(sigma0min_.size()==0) {
477 0 : sigma0min_.resize(pf_n_);
478 0 : for(unsigned i=0; i<pf_n_; i++) {sigma0min_[i]=-1.;}
479 : }
480 :
481 8 : parseVector("SIGMA_MAX",sigma0max_);
482 4 : if(sigma0max_.size()>0 && sigma0max_.size()!=pf_n_) {
483 0 : error("the number of SIGMA_MAX values be the same of the number of the arguments/PF");
484 4 : } else if(sigma0max_.size()==0) {
485 4 : sigma0max_.resize(pf_n_);
486 12 : for(unsigned i=0; i<pf_n_; i++) {sigma0max_[i]=-1.;}
487 : }
488 :
489 12 : for(unsigned i=0; i<pf_n_; i++) {
490 : std::vector<double> tmp_smin, tmp_smax;
491 8 : tmp_smin.resize(1,sigma0min_[i]);
492 8 : tmp_smax.resize(1,sigma0max_[i]);
493 16 : flexbin_.push_back(FlexibleBin(adaptive_,this,i,sigma0_[0],tmp_smin,tmp_smax));
494 : }
495 : }
496 :
497 : // note: HEIGHT is not compulsory, since one could use the TAU keyword, see below
498 42 : parse("HEIGHT",height0_);
499 42 : parse("PACE",stride_);
500 42 : if(stride_<=0) error("frequency for hill addition is nonsensical");
501 :
502 :
503 84 : parseVector("FILE",hillsfname_);
504 42 : if(hillsfname_.size()==0) {
505 30 : for(unsigned i=0; i< pf_n_; i++) {
506 20 : std::string name = do_pf_ ? "HILLS.PF"+std::to_string(i) : "HILLS."+getPntrToArgument(i)->getName();
507 20 : hillsfname_.push_back(name);
508 : }
509 : }
510 42 : if( hillsfname_.size()!=pf_n_ ) {
511 0 : error("number of FILE arguments does not match number of HILLS files");
512 : }
513 :
514 42 : parse("BIASFACTOR",biasf_);
515 42 : if( biasf_<1.0 ) error("well tempered bias factor is nonsensical");
516 42 : kbt_=getkBT();
517 42 : if(biasf_>1.0) {
518 41 : if(kbt_==0.0) error("Unless the MD engine passes the temperature to plumed, with well-tempered metad you must specify it using TEMP");
519 41 : welltemp_=true;
520 : }
521 42 : double tau=0.0;
522 42 : parse("TAU",tau);
523 42 : if(tau==0.0) {
524 42 : if(height0_==std::numeric_limits<double>::max()) error("At least one between HEIGHT and TAU should be specified");
525 : // if tau is not set, we compute it here from the other input parameters
526 42 : if(welltemp_) tau=(kbt_*(biasf_-1.0))/height0_*getTimeStep()*stride_;
527 : } else {
528 0 : if(!welltemp_)error("TAU only makes sense in well-tempered metadynamics");
529 0 : if(height0_!=std::numeric_limits<double>::max()) error("At most one between HEIGHT and TAU should be specified");
530 0 : height0_=(kbt_*(biasf_-1.0))/tau*getTimeStep()*stride_;
531 : }
532 :
533 :
534 : // Multiple walkers
535 42 : parse("WALKERS_N",mw_n_);
536 42 : parse("WALKERS_ID",mw_id_);
537 42 : if(mw_n_<=mw_id_) error("walker ID should be a numerical value less than the total number of walkers");
538 42 : parse("WALKERS_DIR",mw_dir_);
539 42 : parse("WALKERS_RSTRIDE",mw_rstride_);
540 :
541 : // MPI version
542 42 : parseFlag("WALKERS_MPI",walkers_mpi_);
543 :
544 : // Grid file
545 84 : parse("GRID_WSTRIDE",wgridstride_);
546 : std::vector<std::string> gridfilenames_;
547 42 : parseVector("GRID_WFILES",gridfilenames_);
548 42 : if (wgridstride_ == 0 && gridfilenames_.size() > 0) {
549 0 : error("frequency with which to output grid not specified use GRID_WSTRIDE");
550 : }
551 42 : if(gridfilenames_.size()==0 && wgridstride_ > 0) {
552 12 : for(unsigned i=0; i<pf_n_; i++) {
553 8 : std::string name = do_pf_ ? "GRID.PF"+std::to_string(i) : "GRID."+getPntrToArgument(i)->getName();
554 8 : gridfilenames_.push_back(name);
555 : }
556 : }
557 42 : if(gridfilenames_.size() > 0 && hillsfname_.size() > 0 && gridfilenames_.size() != hillsfname_.size())
558 0 : error("number of GRID_WFILES arguments does not match number of HILLS files");
559 :
560 : // Read grid
561 : std::vector<std::string> gridreadfilenames_;
562 42 : parseVector("GRID_RFILES",gridreadfilenames_);
563 :
564 : // Grid Stuff
565 42 : std::vector<std::string> gmin(pf_n_);
566 84 : parseVector("GRID_MIN",gmin);
567 42 : if(gmin.size()!=pf_n_ && gmin.size()!=0) error("not enough values for GRID_MIN");
568 42 : std::vector<std::string> gmax(pf_n_);
569 84 : parseVector("GRID_MAX",gmax);
570 42 : if(gmax.size()!=pf_n_ && gmax.size()!=0) error("not enough values for GRID_MAX");
571 42 : std::vector<unsigned> gbin(pf_n_);
572 : std::vector<double> gspacing;
573 84 : parseVector("GRID_BIN",gbin);
574 42 : if(gbin.size()!=pf_n_ && gbin.size()!=0) error("not enough values for GRID_BIN");
575 84 : parseVector("GRID_SPACING",gspacing);
576 42 : if(gspacing.size()!=pf_n_ && gspacing.size()!=0) error("not enough values for GRID_SPACING");
577 42 : if(gmin.size()!=gmax.size()) error("GRID_MAX and GRID_MIN should be either present or absent");
578 42 : if(gspacing.size()!=0 && gmin.size()==0) error("If GRID_SPACING is present also GRID_MIN and GRID_MAX should be present");
579 42 : if(gbin.size()!=0 && gmin.size()==0) error("If GRID_BIN is present also GRID_MIN and GRID_MAX should be present");
580 42 : if(gmin.size()!=0) {
581 10 : if(gbin.size()==0 && gspacing.size()==0) {
582 10 : if(adaptive_==FlexibleBin::none) {
583 6 : log<<" Binsize not specified, 1/5 of sigma will be be used\n";
584 6 : plumed_assert(sigma0_.size()==pf_n_);
585 6 : gspacing.resize(pf_n_);
586 18 : for(unsigned i=0; i<gspacing.size(); i++) gspacing[i]=0.2*sigma0_[i];
587 : } else {
588 : // with adaptive hills and grid a sigma min must be specified
589 12 : 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");
590 4 : log<<" Binsize not specified, 1/5 of sigma_min will be be used\n";
591 4 : gspacing.resize(pf_n_);
592 12 : for(unsigned i=0; i<gspacing.size(); i++) gspacing[i]=0.2*sigma0min_[i];
593 : }
594 0 : } else if(gspacing.size()!=0 && gbin.size()==0) {
595 0 : log<<" The number of bins will be estimated from GRID_SPACING\n";
596 0 : } else if(gspacing.size()!=0 && gbin.size()!=0) {
597 0 : log<<" You specified both GRID_BIN and GRID_SPACING\n";
598 0 : log<<" The more conservative (highest) number of bins will be used for each variable\n";
599 : }
600 10 : if(gbin.size()==0) gbin.assign(pf_n_,1);
601 30 : if(gspacing.size()!=0) for(unsigned i=0; i<pf_n_; i++) {
602 : double a,b;
603 20 : Tools::convert(gmin[i],a);
604 20 : Tools::convert(gmax[i],b);
605 20 : unsigned n=std::ceil(((b-a)/gspacing[i]));
606 20 : if(gbin[i]<n) gbin[i]=n;
607 : }
608 : }
609 42 : if(gbin.size()>0) grid_=true;
610 :
611 42 : bool sparsegrid=false;
612 42 : parseFlag("GRID_SPARSE",sparsegrid);
613 42 : bool nospline=false;
614 42 : parseFlag("GRID_NOSPLINE",nospline);
615 42 : bool spline=!nospline;
616 42 : if(!grid_&&gridfilenames_.size() > 0) error("To write a grid you need first to define it!");
617 42 : if(!grid_&&gridreadfilenames_.size() > 0) error("To read a grid you need first to define it!");
618 :
619 42 : doInt_.resize(pf_n_,false);
620 : // Interval keyword
621 42 : parseVector("INTERVAL_MIN",lowI_);
622 84 : parseVector("INTERVAL_MAX",uppI_);
623 : // various checks
624 42 : if(lowI_.size()!=uppI_.size()) error("both a lower and an upper limits must be provided with INTERVAL");
625 42 : if(lowI_.size()!=0 && lowI_.size()!=pf_n_) error("check number of argument of INTERVAL");
626 50 : for(unsigned i=0; i<lowI_.size(); ++i) {
627 8 : if(uppI_[i]<lowI_[i]) error("The Upper limit must be greater than the Lower limit!");
628 8 : if(pfhold_[i]->isPeriodic()) warning("INTERVAL is not used for periodic variables");
629 : else doInt_[i]=true;
630 : }
631 :
632 : // parse selector stuff
633 84 : parse("SELECTOR", selector_);
634 42 : if(selector_.length()>0) {
635 1 : do_select_ = true;
636 1 : select_value_ = 0; // set defalt value or it might be not initialized if the user does not pass SELECTOR_ID
637 2 : parse("SELECTOR_ID", select_value_);
638 : }
639 :
640 42 : checkRead();
641 :
642 42 : log.printf(" Gaussian width ");
643 42 : if (adaptive_==FlexibleBin::diffusion)log.printf(" (Note: The units of sigma are in timesteps) ");
644 42 : if (adaptive_==FlexibleBin::geometry)log.printf(" (Note: The units of sigma are in dist units) ");
645 122 : for(unsigned i=0; i<sigma0_.size(); ++i) log.printf(" %f",sigma0_[i]);
646 42 : log.printf(" Gaussian height %f\n",height0_);
647 42 : log.printf(" Gaussian deposition pace %d\n",stride_);
648 42 : log.printf(" Gaussian files ");
649 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) log.printf("%s ",hillsfname_[i].c_str());
650 42 : log.printf("\n");
651 42 : if(welltemp_) {
652 41 : log.printf(" Well-Tempered Bias Factor %f\n",biasf_);
653 41 : log.printf(" Hills relaxation time (tau) %f\n",tau);
654 41 : log.printf(" KbT %f\n",kbt_);
655 : }
656 :
657 42 : if(do_select_) {
658 1 : log.printf(" Add forces and update bias based on the value of SELECTOR %s\n",selector_.c_str());
659 1 : log.printf(" Id of the SELECTOR for this action %u\n", select_value_);
660 : }
661 :
662 42 : if(mw_n_>1) {
663 0 : if(walkers_mpi_) error("MPI version of multiple walkers is not compatible with filesystem version of multiple walkers");
664 0 : log.printf(" %d multiple walkers active\n",mw_n_);
665 0 : log.printf(" walker id %d\n",mw_id_);
666 0 : log.printf(" reading stride %d\n",mw_rstride_);
667 0 : if(mw_dir_!="")log.printf(" directory with hills files %s\n",mw_dir_.c_str());
668 : } else {
669 42 : if(walkers_mpi_) {
670 34 : log.printf(" Multiple walkers active using MPI communnication\n");
671 34 : if(mw_dir_!="")log.printf(" directory with hills files %s\n",mw_dir_.c_str());
672 34 : if(comm.Get_rank()==0) {
673 : // Only root of group can communicate with other walkers
674 18 : mpi_nw_ = multi_sim_comm.Get_size();
675 18 : mpi_id_ = multi_sim_comm.Get_rank();
676 : }
677 : // Communicate to the other members of the same group
678 : // info abount number of walkers and walker index
679 34 : comm.Bcast(mpi_nw_,0);
680 34 : comm.Bcast(mpi_id_,0);
681 : }
682 : }
683 :
684 126 : for(unsigned i=0; i<doInt_.size(); i++) {
685 84 : if(doInt_[i]) log.printf(" Upper and Lower limits boundaries for the bias of CV %u are activated\n", i);
686 : }
687 42 : if(grid_) {
688 10 : log.printf(" Grid min");
689 30 : for(unsigned i=0; i<gmin.size(); ++i) log.printf(" %s",gmin[i].c_str() );
690 10 : log.printf("\n");
691 10 : log.printf(" Grid max");
692 30 : for(unsigned i=0; i<gmax.size(); ++i) log.printf(" %s",gmax[i].c_str() );
693 10 : log.printf("\n");
694 10 : log.printf(" Grid bin");
695 30 : for(unsigned i=0; i<gbin.size(); ++i) log.printf(" %u",gbin[i]);
696 10 : log.printf("\n");
697 10 : if(spline) {log.printf(" Grid uses spline interpolation\n");}
698 10 : if(sparsegrid) {log.printf(" Grid uses sparse grid\n");}
699 10 : if(wgridstride_>0) {
700 18 : for(unsigned i=0; i<gridfilenames_.size(); ++i) {
701 12 : log.printf(" Grid is written on file %s with stride %d\n",gridfilenames_[i].c_str(),wgridstride_);
702 : }
703 : }
704 10 : if(gridreadfilenames_.size()>0) {
705 3 : for(unsigned i=0; i<gridreadfilenames_.size(); ++i) {
706 2 : log.printf(" Reading bias from grid in file %s \n",gridreadfilenames_[i].c_str());
707 : }
708 : }
709 : }
710 :
711 : // initializing vector of hills
712 42 : hills_.resize(pf_n_);
713 :
714 : // restart from external grid
715 : bool restartedFromGrid=false;
716 :
717 : // initializing and checking grid
718 42 : if(grid_) {
719 : // check for mesh and sigma size
720 30 : for(unsigned i=0; i<pf_n_; i++) {
721 : double a,b;
722 20 : int family = pfs_[i]; // point to families instead of arguments
723 20 : Tools::convert(gmin[family],a);
724 20 : Tools::convert(gmax[family],b);
725 20 : double mesh=(b-a)/((double)gbin[family]);
726 20 : if(adaptive_==FlexibleBin::none) {
727 12 : if(mesh>0.5*sigma0_[i]) log<<" WARNING: Using a PBMETAD with a Grid Spacing larger than half of the Gaussians width can produce artifacts\n";
728 : } else {
729 8 : if(mesh>0.5*sigma0min_[i]||sigma0min_[i]<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";
730 : }
731 : }
732 10 : std::string funcl=getLabel() + ".bias";
733 30 : for(unsigned i=0; i<pf_n_; ++i) {
734 20 : std::vector<Value*> args(1);
735 20 : args[0] = pfhold_[i]; //Use first argument in family for interactions.
736 20 : std::vector<std::string> gmin_t(1);
737 20 : std::vector<std::string> gmax_t(1);
738 20 : std::vector<unsigned> gbin_t(1);
739 : gmin_t[0] = gmin[i];
740 : gmax_t[0] = gmax[i];
741 20 : gbin_t[0] = gbin[i];
742 20 : std::unique_ptr<GridBase> BiasGrid_;
743 : // Read grid from file
744 20 : if(gridreadfilenames_.size()>0) {
745 2 : IFile gridfile;
746 2 : gridfile.link(*this);
747 2 : if(gridfile.FileExist(gridreadfilenames_[i])) {
748 2 : gridfile.open(gridreadfilenames_[i]);
749 : } else {
750 0 : error("The GRID file you want to read: " + gridreadfilenames_[i] + ", cannot be found!");
751 : }
752 2 : std::string funcl = getLabel() + ".bias";
753 4 : BiasGrid_=GridBase::create(funcl, args, gridfile, gmin_t, gmax_t, gbin_t, sparsegrid, spline, true);
754 2 : if(BiasGrid_->getDimension() != args.size()) {
755 0 : error("mismatch between dimensionality of input grid and number of arguments");
756 : }
757 4 : if(pfhold_[i]->isPeriodic() != BiasGrid_->getIsPeriodic()[0]) {
758 0 : error("periodicity mismatch between arguments and input bias");
759 : }
760 2 : log.printf(" Restarting from %s:\n",gridreadfilenames_[i].c_str());
761 2 : if(getRestart()) restartedFromGrid=true;
762 2 : } else {
763 18 : if(!sparsegrid) {BiasGrid_=Tools::make_unique<Grid>(funcl,args,gmin_t,gmax_t,gbin_t,spline,true);}
764 0 : else {BiasGrid_=Tools::make_unique<SparseGrid>(funcl,args,gmin_t,gmax_t,gbin_t,spline,true);}
765 18 : std::vector<std::string> actualmin=BiasGrid_->getMin();
766 18 : std::vector<std::string> actualmax=BiasGrid_->getMax();
767 : std::string is;
768 18 : Tools::convert(i,is);
769 18 : if(gmin_t[0]!=actualmin[0]) error("GRID_MIN["+is+"] must be adjusted to "+actualmin[0]+" to fit periodicity");
770 18 : if(gmax_t[0]!=actualmax[0]) error("GRID_MAX["+is+"] must be adjusted to "+actualmax[0]+" to fit periodicity");
771 18 : }
772 20 : BiasGrids_.emplace_back(std::move(BiasGrid_));
773 40 : }
774 : }
775 :
776 :
777 :
778 : // creating vector of ifile* for hills reading
779 : // open all files at the beginning and read Gaussians if restarting
780 :
781 84 : for(int j=0; j<mw_n_; ++j) {
782 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
783 84 : unsigned k=j*hillsfname_.size()+i;
784 : std::string fname;
785 84 : if(mw_dir_!="") {
786 0 : if(mw_n_>1) {
787 0 : std::stringstream out; out << j;
788 0 : fname = mw_dir_+"/"+hillsfname_[i]+"."+out.str();
789 0 : } else if(walkers_mpi_) {
790 0 : fname = mw_dir_+"/"+hillsfname_[i];
791 : } else {
792 : fname = hillsfname_[i];
793 : }
794 : } else {
795 84 : if(mw_n_>1) {
796 0 : std::stringstream out; out << j;
797 0 : fname = hillsfname_[i]+"."+out.str();
798 0 : } else {
799 : fname = hillsfname_[i];
800 : }
801 : }
802 84 : ifiles_.emplace_back(Tools::make_unique<IFile>());
803 : // this is just a shortcut pointer to the last element:
804 : IFile *ifile = ifiles_.back().get();
805 84 : ifile->link(*this);
806 84 : ifilesnames_.push_back(fname);
807 84 : if(ifile->FileExist(fname)) {
808 18 : ifile->open(fname);
809 18 : if(getRestart()&&!restartedFromGrid) {
810 2 : log.printf(" Restarting from %s:",ifilesnames_[k].c_str());
811 2 : readGaussians(i,ifiles_[k].get());
812 : }
813 18 : ifiles_[k]->reset(false);
814 : // close only the walker own hills file for later writing
815 18 : if(j==mw_id_) ifiles_[k]->close();
816 : } else {
817 : // in case a file does not exist and we are restarting, complain that the file was not found
818 66 : if(getRestart()) log<<" WARNING: restart file "<<fname<<" not found\n";
819 : }
820 : }
821 : }
822 :
823 42 : comm.Barrier();
824 42 : if(comm.Get_rank()==0 && walkers_mpi_) multi_sim_comm.Barrier();
825 :
826 : // open hills files for writing
827 126 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
828 : auto ofile=Tools::make_unique<OFile>();
829 84 : ofile->link(*this);
830 : // if MPI multiple walkers, only rank 0 will write to file
831 84 : if(walkers_mpi_) {
832 68 : int r=0;
833 68 : if(comm.Get_rank()==0) r=multi_sim_comm.Get_rank();
834 68 : comm.Bcast(r,0);
835 68 : if(r>0) ifilesnames_[mw_id_*hillsfname_.size()+i]="/dev/null";
836 136 : ofile->enforceSuffix("");
837 : }
838 84 : if(mw_n_>1) ofile->enforceSuffix("");
839 84 : ofile->open(ifilesnames_[mw_id_*hillsfname_.size()+i]);
840 84 : if(fmt_.length()>0) ofile->fmtField(fmt_);
841 168 : ofile->addConstantField("multivariate");
842 168 : ofile->addConstantField("kerneltype");
843 84 : if(doInt_[i]) {
844 16 : ofile->addConstantField("lower_int").printField("lower_int",lowI_[i]);
845 16 : ofile->addConstantField("upper_int").printField("upper_int",uppI_[i]);
846 : }
847 84 : ofile->setHeavyFlush();
848 : // output periodicities of variables
849 84 : ofile->setupPrintValue( pfhold_[i] ); //assuming cvs in the same family have the same periodicity and boundaries.
850 : // push back
851 84 : hillsOfiles_.emplace_back(std::move(ofile));
852 84 : }
853 :
854 : // Dump grid to files
855 42 : if(wgridstride_ > 0) {
856 18 : for(unsigned i = 0; i < gridfilenames_.size(); ++i) {
857 : auto ofile=Tools::make_unique<OFile>();
858 12 : ofile->link(*this);
859 12 : std::string gridfname_tmp = gridfilenames_[i];
860 12 : if(walkers_mpi_) {
861 8 : int r = 0;
862 8 : if(comm.Get_rank() == 0) {
863 4 : r = multi_sim_comm.Get_rank();
864 : }
865 8 : comm.Bcast(r, 0);
866 8 : if(r>0) {
867 : gridfname_tmp = "/dev/null";
868 : }
869 16 : ofile->enforceSuffix("");
870 : }
871 12 : if(mw_n_>1) ofile->enforceSuffix("");
872 12 : ofile->open(gridfname_tmp);
873 12 : ofile->setHeavyFlush();
874 12 : gridfiles_.emplace_back(std::move(ofile));
875 12 : }
876 : }
877 :
878 84 : log<<" Bibliography "<<plumed.cite("Pfaendtner and Bonomi. J. Chem. Theory Comput. 11, 5062 (2015)");
879 50 : if(doInt_[0]) log<<plumed.cite(
880 8 : "Baftizadeh, Cossio, Pietrucci, and Laio, Curr. Phys. Chem. 2, 79 (2012)");
881 110 : if(mw_n_>1||walkers_mpi_) log<<plumed.cite(
882 68 : "Raiteri, Laio, Gervasio, Micheletti, and Parrinello, J. Phys. Chem. B 110, 3533 (2006)");
883 50 : if(adaptive_!=FlexibleBin::none) log<<plumed.cite(
884 8 : "Branduardi, Bussi, and Parrinello, J. Chem. Theory Comput. 8, 2247 (2012)");
885 46 : if (do_pf_) log<<plumed.cite("Prakash, Fu, Bonomi, and Pfaendtner, J. Chem. Theory Comput. 14, 4985 (2018)");
886 42 : log<<"\n";
887 :
888 :
889 84 : }
890 :
891 2 : void PBMetaD::readGaussians(unsigned iarg, IFile *ifile)
892 : {
893 2 : std::vector<double> center(1);
894 2 : std::vector<double> sigma(1);
895 : double height;
896 : int nhills=0;
897 2 : bool multivariate=false;
898 2 : int family=pfs_[iarg];
899 :
900 : std::vector<Value> tmpvalues;
901 4 : tmpvalues.push_back( Value( this, pfhold_[family]->getName(), false ) );
902 :
903 10 : while(scanOneHill(iarg,ifile,tmpvalues,center,sigma,height,multivariate)) {
904 : ;
905 8 : nhills++;
906 8 : if(welltemp_) {height*=(biasf_-1.0)/biasf_;}
907 8 : addGaussian(family, Gaussian(center,sigma,height,multivariate));
908 : }
909 2 : log.printf(" %d Gaussians read\n",nhills);
910 4 : }
911 :
912 1112 : void PBMetaD::writeGaussian(unsigned iarg, const Gaussian& hill, OFile *ofile)
913 : {
914 1112 : int family=pfs_[iarg];
915 2224 : ofile->printField("time",getTimeStep()*getStep());
916 1112 : ofile->printField(pfhold_[family],hill.center[0]);
917 :
918 2224 : ofile->printField("kerneltype","stretched-gaussian");
919 1112 : if(hill.multivariate) {
920 288 : ofile->printField("multivariate","true");
921 144 : double lower = std::sqrt(1./hill.sigma[0]);
922 288 : ofile->printField("sigma_"+pfhold_[family]->getName()+"_"+
923 144 : pfhold_[family]->getName(),lower);
924 : } else {
925 1936 : ofile->printField("multivariate","false");
926 1936 : ofile->printField("sigma_"+pfhold_[family]->getName(),hill.sigma[0]);
927 : }
928 1112 : double height=hill.height;
929 1112 : if(welltemp_) height *= biasf_/(biasf_-1.0);
930 1112 : ofile->printField("height",height);
931 1112 : ofile->printField("biasf",biasf_);
932 1112 : if(mw_n_>1) ofile->printField("clock",int(std::time(0)));
933 1112 : ofile->printField();
934 1112 : }
935 :
936 1120 : void PBMetaD::addGaussian(unsigned iarg, const Gaussian& hill)
937 : {
938 1120 : if(!grid_) {hills_[iarg].push_back(hill);}
939 : else {
940 208 : std::vector<unsigned> nneighb=getGaussianSupport(iarg, hill);
941 208 : std::vector<Grid::index_t> neighbors=BiasGrids_[iarg]->getNeighbors(hill.center,nneighb);
942 208 : std::vector<double> der(1);
943 208 : std::vector<double> xx(1);
944 208 : if(comm.Get_size()==1) {
945 1520 : for(unsigned i=0; i<neighbors.size(); ++i) {
946 1480 : Grid::index_t ineigh=neighbors[i];
947 1480 : der[0]=0.0;
948 1480 : BiasGrids_[iarg]->getPoint(ineigh,xx);
949 1480 : double bias=evaluateGaussian(iarg,xx,hill,&der[0]);
950 1480 : BiasGrids_[iarg]->addValueAndDerivatives(ineigh,bias,der);
951 : }
952 : } else {
953 168 : unsigned stride=comm.Get_size();
954 168 : unsigned rank=comm.Get_rank();
955 168 : std::vector<double> allder(neighbors.size(),0.0);
956 168 : std::vector<double> allbias(neighbors.size(),0.0);
957 3276 : for(unsigned i=rank; i<neighbors.size(); i+=stride) {
958 3108 : Grid::index_t ineigh=neighbors[i];
959 3108 : BiasGrids_[iarg]->getPoint(ineigh,xx);
960 3108 : allbias[i]=evaluateGaussian(iarg,xx,hill,&allder[i]);
961 : }
962 168 : comm.Sum(allbias);
963 168 : comm.Sum(allder);
964 6384 : for(unsigned i=0; i<neighbors.size(); ++i) {
965 6216 : Grid::index_t ineigh=neighbors[i];
966 6216 : der[0]=allder[i];
967 6216 : BiasGrids_[iarg]->addValueAndDerivatives(ineigh,allbias[i],der);
968 : }
969 : }
970 : }
971 1120 : }
972 :
973 208 : std::vector<unsigned> PBMetaD::getGaussianSupport(unsigned iarg, const Gaussian& hill)
974 : {
975 : std::vector<unsigned> nneigh;
976 : double cutoff;
977 208 : if(hill.multivariate) {
978 144 : double maxautoval=1./hill.sigma[0];
979 144 : cutoff=std::sqrt(2.0*dp2cutoff*maxautoval);
980 : } else {
981 64 : cutoff=std::sqrt(2.0*dp2cutoff)*hill.sigma[0];
982 : }
983 :
984 208 : if(doInt_[iarg]) {
985 144 : if(hill.center[0]+cutoff > uppI_[iarg] || hill.center[0]-cutoff < lowI_[iarg]) {
986 : // in this case, we updated the entire grid to avoid problems
987 0 : return BiasGrids_[iarg]->getNbin();
988 : } else {
989 144 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff/BiasGrids_[iarg]->getDx()[0])));
990 : return nneigh;
991 : }
992 : }
993 :
994 64 : nneigh.push_back( static_cast<unsigned>(ceil(cutoff/BiasGrids_[iarg]->getDx()[0])) );
995 :
996 : return nneigh;
997 : }
998 :
999 1296 : double PBMetaD::getBiasAndDerivatives(unsigned iarg, const std::vector<double>& cv, double* der)
1000 : {
1001 1296 : double bias=0.0;
1002 1296 : int family = pfs_[iarg];
1003 1296 : if(!grid_) {
1004 1000 : unsigned stride=comm.Get_size();
1005 1000 : unsigned rank=comm.Get_rank();
1006 5520 : for(unsigned i=rank; i<hills_[family].size(); i+=stride) {
1007 4520 : bias += evaluateGaussian(iarg,cv,hills_[family][i],der);
1008 : }
1009 1000 : comm.Sum(bias);
1010 1000 : if(der) comm.Sum(der,1);
1011 : } else {
1012 296 : if(der) {
1013 160 : std::vector<double> vder(1);
1014 160 : bias = BiasGrids_[family]->getValueAndDerivatives(cv,vder);
1015 160 : der[0] = vder[0];
1016 : } else {
1017 136 : bias = BiasGrids_[family]->getValue(cv);
1018 : }
1019 : }
1020 :
1021 1296 : return bias;
1022 : }
1023 :
1024 9108 : double PBMetaD::evaluateGaussian(unsigned iarg, const std::vector<double>& cv, const Gaussian& hill, double* der)
1025 : {
1026 : double bias=0.0;
1027 : // I use a pointer here because cv is const (and should be const)
1028 : // but when using doInt it is easier to locally replace cv[0] with
1029 : // the upper/lower limit in case it is out of range
1030 : const double *pcv=NULL;
1031 : double tmpcv[1]; // tmp array with cv (to be used with doInt_)
1032 9108 : tmpcv[0]=cv[0];
1033 : bool isOutOfInt = false;
1034 9108 : if(doInt_[iarg]) {
1035 2664 : if(cv[0]<lowI_[iarg]) { tmpcv[0]=lowI_[iarg]; isOutOfInt = true; }
1036 2664 : else if(cv[0]>uppI_[iarg]) { tmpcv[0]=uppI_[iarg]; isOutOfInt = true; }
1037 : }
1038 : pcv=&(tmpcv[0]);
1039 :
1040 9108 : if(hill.multivariate) {
1041 2664 : double dp = difference(iarg, hill.center[0], pcv[0]);
1042 2664 : double dp2 = 0.5 * dp * dp * hill.sigma[0];
1043 2664 : if(dp2<dp2cutoff) {
1044 2534 : bias = hill.height*std::exp(-dp2);
1045 2534 : if(der && !isOutOfInt) {
1046 2534 : der[0] += -bias * dp * hill.sigma[0] * stretchA;
1047 : }
1048 2534 : bias=stretchA*bias+hill.height*stretchB;
1049 : }
1050 : } else {
1051 6444 : double dp = difference(iarg, hill.center[0], pcv[0]) * hill.invsigma[0];
1052 6444 : double dp2 = 0.5 * dp * dp;
1053 6444 : if(dp2<dp2cutoff) {
1054 6363 : bias = hill.height*std::exp(-dp2);
1055 6363 : if(der && !isOutOfInt) {
1056 4107 : der[0] += -bias * dp * hill.invsigma[0] * stretchA;
1057 : }
1058 6363 : bias=stretchA*bias+hill.height*stretchB;
1059 : }
1060 : }
1061 :
1062 9108 : return bias;
1063 : }
1064 :
1065 340 : void PBMetaD::calculate()
1066 : {
1067 : // this is because presently there is no way to properly pass information
1068 : // on adaptive hills (diff) after exchanges:
1069 340 : if(adaptive_==FlexibleBin::diffusion && getExchangeStep()) error("ADAPTIVE=DIFF is not compatible with replica exchange");
1070 :
1071 340 : std::vector<double> cv(1);
1072 : double der[1];
1073 340 : std::vector<double> bias(getNumberOfArguments());
1074 340 : std::vector<double> deriv(getNumberOfArguments());
1075 :
1076 340 : double ncv = (double) getNumberOfArguments();
1077 : double bmin = 1.0e+19;
1078 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1079 700 : cv[0] = getArgument(i);
1080 700 : der[0] = 0.0;
1081 700 : bias[i] = getBiasAndDerivatives(i, cv, der);
1082 700 : deriv[i] = der[0];
1083 700 : if(bias[i] < bmin) bmin = bias[i];
1084 : }
1085 : double ene = 0.;
1086 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1087 700 : ene += std::exp((-bias[i]+bmin)/kbt_);
1088 : }
1089 :
1090 : // set Forces - set them to zero if SELECTOR is active
1091 340 : if(do_select_) current_value_ = static_cast<unsigned>(plumed.passMap[selector_]);
1092 :
1093 340 : if(!do_select_ || select_value_==current_value_) {
1094 1040 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1095 700 : const double f = - std::exp((-bias[i]+bmin)/kbt_) / (ene) * deriv[i];
1096 700 : setOutputForce(i, f);
1097 : }
1098 : }
1099 :
1100 340 : if(do_select_ && select_value_!=current_value_) {
1101 0 : for(unsigned i=0; i<getNumberOfArguments(); ++i) setOutputForce(i, 0.0);
1102 : }
1103 :
1104 : // set bias
1105 340 : ene = -kbt_ * (std::log(ene) - std::log(ncv)) + bmin;
1106 340 : setBias(ene);
1107 340 : }
1108 :
1109 340 : void PBMetaD::update()
1110 : {
1111 : bool multivariate;
1112 : // adding hills criteria
1113 : bool nowAddAHill;
1114 340 : if(getStep()%stride_==0 && !isFirstStep_) nowAddAHill=true;
1115 : else {
1116 : nowAddAHill=false;
1117 50 : isFirstStep_=false;
1118 : }
1119 :
1120 : // if you use adaptive, call the FlexibleBin
1121 340 : if(adaptive_!=FlexibleBin::none) {
1122 120 : for(unsigned i=0; i<getNumberOfArguments(); i++) flexbin_[i].update(nowAddAHill,i);
1123 : multivariate=true;
1124 : } else {
1125 : multivariate=false;
1126 : }
1127 :
1128 340 : if(nowAddAHill && (!do_select_ || select_value_==current_value_)) {
1129 : // get all biases and heights
1130 290 : std::vector<double> cv(getNumberOfArguments());
1131 290 : std::vector<double> bias(getNumberOfArguments());
1132 290 : std::vector<double> thissigma(getNumberOfArguments());
1133 290 : std::vector<double> height(getNumberOfArguments());
1134 290 : std::vector<double> cv_tmp(1);
1135 290 : std::vector<double> sigma_tmp(1);
1136 : double norm = 0.0;
1137 : double bmin = 1.0e+19;
1138 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1139 596 : int family=pfs_[i];
1140 : // get flex/sigmas for each family and assign them to this args sigma
1141 596 : if(adaptive_!=FlexibleBin::none) thissigma[i]=flexbin_[family].getInverseMatrix(i)[0];
1142 524 : else thissigma[i]=sigma0_[family];
1143 596 : cv[i] = getArgument(i);
1144 596 : cv_tmp[0] = getArgument(i);
1145 596 : bias[i] = getBiasAndDerivatives(i, cv_tmp);
1146 596 : if(bias[i] < bmin) bmin = bias[i];
1147 : }
1148 : // calculate heights and norm
1149 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1150 596 : double h = std::exp((-bias[i]+bmin)/kbt_);
1151 596 : norm += h;
1152 596 : height[i] = h;
1153 : }
1154 : // normalize and apply welltemp correction
1155 886 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1156 596 : height[i] *= height0_ / norm;
1157 596 : if(welltemp_) height[i] *= std::exp(-bias[i]/(kbt_*(biasf_-1.0)));
1158 : }
1159 :
1160 : // MPI Multiple walkers: share hills and add them all
1161 290 : if(walkers_mpi_) {
1162 : // Allocate arrays to store all walkers hills
1163 258 : std::vector<double> all_cv(mpi_nw_*cv.size(), 0.0);
1164 258 : std::vector<double> all_sigma(mpi_nw_*getNumberOfArguments(), 0.0);
1165 258 : std::vector<double> all_height(mpi_nw_*height.size(), 0.0);
1166 258 : if(comm.Get_rank()==0) {
1167 : // fill in value
1168 390 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1169 260 : unsigned j = mpi_id_ * getNumberOfArguments() + i;
1170 260 : all_cv[j] = cv[i];
1171 260 : all_sigma[j] = thissigma[i];
1172 260 : all_height[j] = height[i];
1173 : }
1174 : // Communicate (only root)
1175 130 : multi_sim_comm.Sum(&all_cv[0], all_cv.size());
1176 130 : multi_sim_comm.Sum(&all_sigma[0], all_sigma.size());
1177 130 : multi_sim_comm.Sum(&all_height[0], all_height.size());
1178 : }
1179 : // Share info with group members
1180 258 : comm.Sum(&all_cv[0], all_cv.size());
1181 258 : comm.Sum(&all_sigma[0], all_sigma.size());
1182 258 : comm.Sum(&all_height[0], all_height.size());
1183 : // now add hills one by one
1184 774 : for(unsigned j=0; j<mpi_nw_; ++j) {
1185 1548 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1186 : // Add CVs of same family together and write to same file
1187 1032 : int family = pfs_[i];
1188 1032 : cv_tmp[0] = all_cv[j*cv.size()+i];
1189 1032 : double height_tmp = all_height[j*cv.size()+i];
1190 1032 : sigma_tmp[0] = all_sigma[j*cv.size()+i];
1191 1032 : Gaussian newhill = Gaussian(cv_tmp, sigma_tmp, height_tmp, multivariate);
1192 1032 : addGaussian(family, newhill);
1193 1032 : writeGaussian(i, newhill, hillsOfiles_[family].get());
1194 1032 : }
1195 : }
1196 : // just add your own hills
1197 : } else {
1198 112 : for(unsigned i=0; i<getNumberOfArguments(); ++i) {
1199 : // Add CVs of same family together and write to same file
1200 80 : int family = pfs_[i];
1201 80 : cv_tmp[0] = cv[i];
1202 80 : if(adaptive_!=FlexibleBin::none) sigma_tmp[0]=thissigma[i];
1203 80 : else sigma_tmp[0] = sigma0_[family];
1204 80 : Gaussian newhill = Gaussian(cv_tmp, sigma_tmp, height[i], multivariate);
1205 80 : addGaussian(family, newhill);
1206 80 : writeGaussian(i, newhill, hillsOfiles_[family].get());
1207 80 : }
1208 : }
1209 : }
1210 :
1211 : // write grid files
1212 340 : if(wgridstride_>0 && (getStep()%wgridstride_==0 || getCPT())) {
1213 14 : int r = 0;
1214 14 : if(walkers_mpi_) {
1215 4 : if(comm.Get_rank()==0) r=multi_sim_comm.Get_rank();
1216 4 : comm.Bcast(r,0);
1217 : }
1218 14 : if(r==0) {
1219 36 : for(unsigned i=0; i<gridfiles_.size(); ++i) {
1220 24 : gridfiles_[i]->rewind();
1221 24 : BiasGrids_[i]->writeToFile(*gridfiles_[i]);
1222 24 : gridfiles_[i]->flush();
1223 : }
1224 : }
1225 : }
1226 :
1227 : // if multiple walkers and time to read Gaussians
1228 340 : if(mw_n_>1 && getStep()%mw_rstride_==0) {
1229 0 : for(int j=0; j<mw_n_; ++j) {
1230 0 : for(unsigned i=0; i<hillsfname_.size(); ++i) {
1231 0 : unsigned k=j*hillsfname_.size()+i;
1232 : // don't read your own Gaussians
1233 0 : if(j==mw_id_) continue;
1234 : // if the file is not open yet
1235 0 : if(!(ifiles_[k]->isOpen())) {
1236 : // check if it exists now and open it!
1237 0 : if(ifiles_[k]->FileExist(ifilesnames_[k])) {
1238 0 : ifiles_[k]->open(ifilesnames_[k]);
1239 0 : ifiles_[k]->reset(false);
1240 : }
1241 : // otherwise read the new Gaussians
1242 : } else {
1243 0 : log.printf(" Reading hills from %s:",ifilesnames_[k].c_str());
1244 0 : readGaussians(i,ifiles_[k].get());
1245 0 : ifiles_[k]->reset(false);
1246 : }
1247 : }
1248 : }
1249 : }
1250 :
1251 340 : }
1252 :
1253 : /// takes a pointer to the file and a template string with values v and gives back the next center, sigma and height
1254 10 : bool PBMetaD::scanOneHill(unsigned iarg, IFile *ifile, std::vector<Value> &tmpvalues, std::vector<double> ¢er, std::vector<double> &sigma, double &height, bool &multivariate)
1255 : {
1256 : double dummy;
1257 10 : multivariate=false;
1258 10 : Value* argPtr = pfhold_[pfs_[iarg]];
1259 20 : if(ifile->scanField("time",dummy)) {
1260 8 : ifile->scanField( &tmpvalues[0] );
1261 8 : if( tmpvalues[0].isPeriodic() && ! argPtr->isPeriodic() ) {
1262 0 : error("in hills file periodicity for variable " + tmpvalues[0].getName() + " does not match periodicity in input");
1263 8 : } else if( tmpvalues[0].isPeriodic() ) {
1264 0 : std::string imin, imax; tmpvalues[0].getDomain( imin, imax );
1265 0 : std::string rmin, rmax; argPtr->getDomain( rmin, rmax );
1266 0 : if( imin!=rmin || imax!=rmax ) {
1267 0 : error("in hills file periodicity for variable " + tmpvalues[0].getName() + " does not match periodicity in input");
1268 : }
1269 : }
1270 8 : center[0]=tmpvalues[0].get();
1271 8 : std::string ktype="stretched-gaussian";
1272 24 : if( ifile->FieldExist("kerneltype") ) ifile->scanField("kerneltype",ktype);
1273 :
1274 8 : if( ktype=="gaussian" ) {
1275 0 : noStretchWarning();
1276 8 : } else if( ktype!="stretched-gaussian") {
1277 0 : error("non Gaussian kernels are not supported in MetaD");
1278 : }
1279 :
1280 : std::string sss;
1281 16 : ifile->scanField("multivariate",sss);
1282 8 : if(sss=="true") multivariate=true;
1283 8 : else if(sss=="false") multivariate=false;
1284 0 : else plumed_merror("cannot parse multivariate = "+ sss);
1285 8 : if(multivariate) {
1286 0 : ifile->scanField("sigma_"+argPtr->getName()+"_"+
1287 : argPtr->getName(),sigma[0]);
1288 0 : sigma[0] = 1./(sigma[0]*sigma[0]);
1289 : } else {
1290 16 : ifile->scanField("sigma_"+argPtr->getName(),sigma[0]);
1291 : }
1292 8 : ifile->scanField("height",height);
1293 8 : ifile->scanField("biasf",dummy);
1294 16 : if(ifile->FieldExist("clock")) ifile->scanField("clock",dummy);
1295 16 : if(ifile->FieldExist("lower_int")) ifile->scanField("lower_int",dummy);
1296 16 : if(ifile->FieldExist("upper_int")) ifile->scanField("upper_int",dummy);
1297 8 : ifile->scanField();
1298 : return true;
1299 : } else {
1300 : return false;
1301 : }
1302 :
1303 : }
1304 :
1305 340 : bool PBMetaD::checkNeedsGradients()const
1306 : {
1307 340 : if(adaptive_==FlexibleBin::geometry) {
1308 0 : if(getStep()%stride_==0 && !isFirstStep_) return true;
1309 0 : else return false;
1310 : } else return false;
1311 : }
1312 :
1313 : }
1314 : }
1315 :
|