Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2020-2021 of Michele Invernizzi.
3 :
4 : This file is part of the OPES plumed module.
5 :
6 : The OPES plumed module is free software: you can redistribute it and/or modify
7 : it under the terms of the GNU Lesser General Public License as published by
8 : the Free Software Foundation, either version 3 of the License, or
9 : (at your option) any later version.
10 :
11 : The OPES plumed module is distributed in the hope that it will be useful,
12 : but WITHOUT ANY WARRANTY; without even the implied warranty of
13 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 : GNU Lesser General Public License for more details.
15 :
16 : You should have received a copy of the GNU Lesser General Public License
17 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
18 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
19 : #include "bias/Bias.h"
20 : #include "core/PlumedMain.h"
21 : #include "core/ActionRegister.h"
22 : #include "core/ActionSet.h"
23 : #include "tools/Communicator.h"
24 : #include "tools/File.h"
25 : #include "tools/OpenMP.h"
26 :
27 : #include "ExpansionCVs.h"
28 :
29 : namespace PLMD {
30 : namespace opes {
31 :
32 : //+PLUMEDOC OPES_BIAS OPES_EXPANDED
33 : /*
34 : On-the-fly probability enhanced sampling (\ref OPES "OPES") with expanded ensembles target distribution (replica-exchange-like) \cite Invernizzi2020unified.
35 :
36 : An expanded ensemble is obtained by summing a set of ensembles at slightly different termodynamic conditions, or with slightly different Hamiltonians.
37 : Such ensembles can be sampled via methods like replica exchange, or this \ref OPES_EXPANDED bias action.
38 : A typical example is a multicanonical simulation, in which a whole range of temperatures is sampled instead of a single one.
39 :
40 : In oreder to define an expanded target ensemble we use \ref EXPANSION_CV "expansion collective variables" (ECVs), \f$\Delta u_\lambda(\mathbf{x})\f$.
41 : The bias at step \f$n\f$ is
42 : \f[
43 : V_n(\mathbf{x})=-\frac{1}{\beta}\log \left(\frac{1}{N_{\{\lambda\}}}\sum_\lambda e^{-\Delta u_\lambda(\mathbf{x})+\beta\Delta F_n(\lambda)}\right)\, .
44 : \f]
45 : See Ref.\cite Invernizzi2020unified for more details on the method.
46 :
47 : Notice that the estimates in the DELTAFS file are expressed in energy units, and should be multiplied by \f$\beta\f$ to be dimensionless as in Ref.\cite Invernizzi2020unified.
48 : The DELTAFS file also contains an estimate of \f$c(t)=\frac{1}{\beta} \log \langle e^{\beta V}\rangle\f$.
49 : Similarly to \ref OPES_METAD, it is printed only for reference, since \f$c(t)\f$ reaches a fixed value as the bias converges, and should NOT be used for reweighting.
50 : Its value is also needed for restarting a simulation.
51 :
52 : You can store the instantaneous \f$\Delta F_n(\lambda)\f$ estimates also in a more readable format using STATE_WFILE and STATE_WSTRIDE.
53 : Restart can be done either from a DELTAFS file or from a STATE_RFILE, it is equivalent.
54 :
55 : Contrary to \ref OPES_METAD, \ref OPES_EXPANDED does not use kernel density estimation.
56 :
57 : \par Examples
58 :
59 : \plumedfile
60 : # simulate multiple temperatures, as in parallel tempering
61 : ene: ENERGY
62 : ecv: ECV_MULTITHERMAL ARG=ene TEMP_MAX=1000
63 : opes: OPES_EXPANDED ARG=ecv.* PACE=500
64 : PRINT FILE=COLVAR STRIDE=500 ARG=ene,opes.bias
65 : \endplumedfile
66 :
67 : You can easily combine multiple ECVs.
68 : The \ref OPES_EXPANDED bias will create a multidimensional target grid to sample all the combinations.
69 :
70 : \plumedfile
71 : # simulate multiple temperatures while biasing a CV
72 : ene: ENERGY
73 : dst: DISTANCE ATOMS=1,2
74 :
75 : ecv1: ECV_MULTITHERMAL ARG=ene TEMP_SET_ALL=200,300,500,1000
76 : ecv2: ECV_UMBRELLAS_LINE ARG=dst CV_MIN=1.2 CV_MAX=4.3 SIGMA=0.5
77 : opes: OPES_EXPANDED ARG=ecv1.*,ecv2.* PACE=500 OBSERVATION_STEPS=1
78 :
79 : PRINT FILE=COLVAR STRIDE=500 ARG=ene,dst,opes.bias
80 : \endplumedfile
81 :
82 : If an ECV is based on more than one CV you must provide all the output component, in the proper order.
83 : You can use \ref Regex for that, like in the following example.
84 :
85 : \plumedfile
86 : # simulate multiple temperatures and pressures while biasing a two-CVs linear path
87 : ene: ENERGY
88 : vol: VOLUME
89 : ecv_mtp: ECV_MULTITHERMAL_MULTIBARIC ...
90 : ARG=ene,vol
91 : TEMP=300
92 : TEMP_MIN=200
93 : TEMP_MAX=800
94 : PRESSURE=0.06022140857*1000 #1 kbar
95 : PRESSURE_MIN=0
96 : PRESSURE_MAX=0.06022140857*2000 #2 kbar
97 : ...
98 :
99 : cv1: DISTANCE ATOMS=1,2
100 : cv2: DISTANCE ATOMS=3,4
101 : ecv_umb: ECV_UMBRELLAS_LINE ARG=cv1,cv2 TEMP=300 CV_MIN=0.1,0.1 CV_MAX=1.5,1.5 SIGMA=0.2 BARRIER=70
102 :
103 : opes: OPES_EXPANDED ARG=(ecv_.*) PACE=500 WALKERS_MPI PRINT_STRIDE=1000
104 :
105 : PRINT FILE=COLVAR STRIDE=500 ARG=ene,vol,cv1,cv2,opes.bias
106 : \endplumedfile
107 :
108 :
109 : */
110 : //+ENDPLUMEDOC
111 :
112 : class OPESexpanded : public bias::Bias {
113 :
114 : private:
115 : bool isFirstStep_;
116 : unsigned NumOMP_;
117 : unsigned NumParallel_;
118 : unsigned rank_;
119 : unsigned NumWalkers_;
120 : unsigned walker_rank_;
121 : unsigned long counter_;
122 : std::size_t ncv_;
123 :
124 : std::vector<const double *> ECVs_;
125 : std::vector<const double *> derECVs_;
126 : std::vector<opes::ExpansionCVs*> pntrToECVsClass_;
127 : std::vector< std::vector<unsigned> > index_k_;
128 : // A note on indexes usage:
129 : // j -> underlying CVs
130 : // i -> DeltaFs
131 : // k -> single ECVs, which might not be trivially numbered
132 : // l -> groups of ECVs, pntrToECVsClass
133 : // h -> subgroups of ECVs, arguments in ECVsClass
134 : // w -> walkers
135 :
136 : double kbt_;
137 : unsigned stride_;
138 : unsigned deltaF_size_; //different from deltaF_.size() if NumParallel_>1
139 : std::vector<double> deltaF_;
140 : std::vector<double> diff_;
141 : double rct_;
142 :
143 : std::vector<double> all_deltaF_;
144 : std::vector<int> all_size_;
145 : std::vector<int> disp_;
146 :
147 : unsigned obs_steps_;
148 : std::vector<double> obs_cvs_;
149 :
150 : bool calc_work_;
151 : double work_;
152 :
153 : unsigned print_stride_;
154 : OFile deltaFsOfile_;
155 : std::vector<std::string> deltaF_name_;
156 :
157 : OFile stateOfile_;
158 : int wStateStride_;
159 : bool storeOldStates_;
160 :
161 : void init_pntrToECVsClass();
162 : void init_linkECVs();
163 : void init_fromObs();
164 :
165 : void printDeltaF();
166 : void dumpStateToFile();
167 : void updateDeltaF(double);
168 : double getExpansion(const unsigned) const;
169 :
170 : public:
171 : explicit OPESexpanded(const ActionOptions&);
172 : void calculate() override;
173 : void update() override;
174 : static void registerKeywords(Keywords& keys);
175 : };
176 :
177 10479 : PLUMED_REGISTER_ACTION(OPESexpanded,"OPES_EXPANDED")
178 :
179 31 : void OPESexpanded::registerKeywords(Keywords& keys)
180 : {
181 31 : Bias::registerKeywords(keys);
182 31 : keys.remove("ARG");
183 62 : keys.add("compulsory","ARG","the label of the ECVs that define the expansion. You can use an * to make sure all the output components of the ECVs are used, as in the examples above");
184 62 : keys.add("compulsory","PACE","how often the bias is updated");
185 93 : keys.add("compulsory","OBSERVATION_STEPS","100","number of unbiased initial PACE steps to collect statistics for initialization");
186 : //DeltaFs and state files
187 62 : keys.add("compulsory","FILE","DELTAFS","a file with the estimate of the relative \\f$\\Delta F\\f$ for each component of the target and of the global \\f$c(t)\\f$");
188 62 : keys.add("compulsory","PRINT_STRIDE","100","stride for printing to DELTAFS file, in units of PACE");
189 62 : keys.add("optional","FMT","specify format for DELTAFS file");
190 62 : keys.add("optional","STATE_RFILE","read from this file the \\f$\\Delta F\\f$ estimates and all the info needed to RESTART the simulation");
191 62 : keys.add("optional","STATE_WFILE","write to this file the \\f$\\Delta F\\f$ estimates and all the info needed to RESTART the simulation");
192 62 : keys.add("optional","STATE_WSTRIDE","number of MD steps between writing the STATE_WFILE. Default is only on CPT events (but not all MD codes set them)");
193 62 : keys.addFlag("STORE_STATES",false,"append to STATE_WFILE instead of ovewriting it each time");
194 : //miscellaneous
195 62 : keys.addFlag("CALC_WORK",false,"calculate the total accumulated work done by the bias since last restart");
196 62 : keys.addFlag("WALKERS_MPI",false,"switch on MPI version of multiple walkers");
197 62 : keys.addFlag("SERIAL",false,"perform calculations in serial");
198 31 : keys.use("RESTART");
199 31 : keys.use("UPDATE_FROM");
200 31 : keys.use("UPDATE_UNTIL");
201 :
202 : //output components
203 31 : componentsAreNotOptional(keys);
204 62 : keys.addOutputComponent("work","CALC_WORK","total accumulated work done by the bias");
205 31 : }
206 :
207 30 : OPESexpanded::OPESexpanded(const ActionOptions&ao)
208 : : PLUMED_BIAS_INIT(ao)
209 30 : , isFirstStep_(true)
210 30 : , counter_(0)
211 30 : , ncv_(getNumberOfArguments())
212 30 : , deltaF_size_(0)
213 30 : , rct_(0)
214 30 : , work_(0)
215 : {
216 : //set pace
217 30 : parse("PACE",stride_);
218 30 : parse("OBSERVATION_STEPS",obs_steps_);
219 30 : plumed_massert(obs_steps_!=0,"minimum is OBSERVATION_STEPS=1");
220 30 : obs_cvs_.resize(obs_steps_*ncv_);
221 :
222 : //deltaFs file
223 : std::string deltaFsFileName;
224 30 : parse("FILE",deltaFsFileName);
225 60 : parse("PRINT_STRIDE",print_stride_);
226 : std::string fmt;
227 60 : parse("FMT",fmt);
228 : //output checkpoint of current state
229 : std::string restartFileName;
230 60 : parse("STATE_RFILE",restartFileName);
231 : std::string stateFileName;
232 30 : parse("STATE_WFILE",stateFileName);
233 30 : wStateStride_=0;
234 30 : parse("STATE_WSTRIDE",wStateStride_);
235 30 : storeOldStates_=false;
236 30 : parseFlag("STORE_STATES",storeOldStates_);
237 30 : if(wStateStride_!=0 || storeOldStates_)
238 5 : plumed_massert(stateFileName.length()>0,"filename for storing simulation status not specified, use STATE_WFILE");
239 30 : if(wStateStride_>0)
240 5 : plumed_massert(wStateStride_>=(int)stride_,"STATE_WSTRIDE is in units of MD steps, thus should be a multiple of PACE");
241 30 : if(stateFileName.length()>0 && wStateStride_==0)
242 1 : wStateStride_=-1;//will print only on CPT events (checkpoints set by some MD engines, like gromacs)
243 :
244 : //work flag
245 30 : parseFlag("CALC_WORK",calc_work_);
246 :
247 : //multiple walkers //external MW for cp2k not supported, but anyway cp2k cannot put bias on energy!
248 30 : bool walkers_mpi=false;
249 30 : parseFlag("WALKERS_MPI",walkers_mpi);
250 30 : if(walkers_mpi)
251 : {
252 4 : if(comm.Get_rank()==0) //multi_sim_comm works on first rank only
253 : {
254 4 : NumWalkers_=multi_sim_comm.Get_size();
255 4 : walker_rank_=multi_sim_comm.Get_rank();
256 : }
257 4 : comm.Bcast(NumWalkers_,0); //if each walker has more than one processor update them all
258 4 : comm.Bcast(walker_rank_,0);
259 : }
260 : else
261 : {
262 26 : NumWalkers_=1;
263 26 : walker_rank_=0;
264 : }
265 :
266 : //parallelization stuff
267 30 : NumOMP_=OpenMP::getNumThreads();
268 30 : NumParallel_=comm.Get_size();
269 30 : rank_=comm.Get_rank();
270 30 : bool serial=false;
271 30 : parseFlag("SERIAL",serial);
272 30 : if(serial)
273 : {
274 5 : NumOMP_=1;
275 5 : NumParallel_=1;
276 5 : rank_=0;
277 : }
278 :
279 30 : checkRead();
280 :
281 : //check ECVs and link them
282 30 : init_pntrToECVsClass();
283 : //set kbt_
284 30 : kbt_=pntrToECVsClass_[0]->getKbT();
285 67 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
286 37 : plumed_massert(std::abs(kbt_-pntrToECVsClass_[l]->getKbT())<1e-4,"must set same TEMP for each ECV");
287 :
288 : //restart if needed
289 30 : if(getRestart())
290 : {
291 : bool stateRestart=true;
292 10 : if(restartFileName.length()==0)
293 : {
294 : stateRestart=false;
295 : restartFileName=deltaFsFileName;
296 : }
297 10 : IFile ifile;
298 10 : ifile.link(*this);
299 10 : if(ifile.FileExist(restartFileName))
300 : {
301 10 : log.printf(" RESTART - make sure all ECVs used are the same as before\n");
302 10 : log.printf(" restarting from: %s\n",restartFileName.c_str());
303 10 : ifile.open(restartFileName);
304 10 : if(stateRestart) //get all info
305 : {
306 2 : log.printf(" it should be a STATE file (not a DELTAFS file)\n");
307 : double time; //not used
308 2 : ifile.scanField("time",time);
309 2 : ifile.scanField("counter",counter_);
310 4 : ifile.scanField("rct",rct_);
311 : std::string tmp_lambda;
312 66 : while(ifile.scanField(getPntrToArgument(0)->getName(),tmp_lambda))
313 : {
314 64 : std::string subs="DeltaF_"+tmp_lambda;
315 128 : for(unsigned jj=1; jj<ncv_; jj++)
316 : {
317 : tmp_lambda.clear();
318 64 : ifile.scanField(getPntrToArgument(jj)->getName(),tmp_lambda);
319 128 : subs+="_"+tmp_lambda;
320 : }
321 64 : deltaF_name_.push_back(subs);
322 : double tmp_deltaF;
323 64 : ifile.scanField("DeltaF",tmp_deltaF);
324 64 : deltaF_.push_back(tmp_deltaF);
325 64 : ifile.scanField();
326 : tmp_lambda.clear();
327 : }
328 2 : log.printf(" successfully read %lu DeltaF values\n",deltaF_name_.size());
329 2 : if(NumParallel_>1)
330 2 : all_deltaF_=deltaF_;
331 : }
332 : else //get just deltaFs names
333 : {
334 8 : ifile.scanFieldList(deltaF_name_);
335 8 : plumed_massert(deltaF_name_.size()>=4,"RESTART - fewer than expected FIELDS found in '"+deltaFsFileName+"' file");
336 8 : plumed_massert(deltaF_name_[deltaF_name_.size()-1]=="print_stride","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
337 8 : plumed_massert(deltaF_name_[0]=="time","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
338 8 : plumed_massert(deltaF_name_[1]=="rct","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
339 8 : deltaF_name_.pop_back();
340 : deltaF_name_.erase(deltaF_name_.begin(),deltaF_name_.begin()+2);
341 : std::size_t pos=5; //each name starts with "DeltaF"
342 22 : for(unsigned j=0; j<ncv_; j++)
343 14 : pos=deltaF_name_[0].find("_",pos+1); //checking only first one, hopefully is enough
344 8 : plumed_massert(pos<deltaF_name_[0].length(),"RESTART - fewer '_' than expected in DeltaF fields: did you remove any CV?");
345 8 : pos=deltaF_name_[0].find("_",pos+1);
346 8 : plumed_massert(pos>deltaF_name_[0].length(),"RESTART - more '_' than expected in DeltaF fields: did you add new CV?");
347 : }
348 : //get lambdas, init ECVs and Link them
349 10 : deltaF_size_=deltaF_name_.size();
350 525 : auto getLambdaName=[](const std::string& name,const unsigned start,const unsigned dim)
351 : {
352 : std::size_t pos_start=5; //each name starts with "DeltaF"
353 1068 : for(unsigned j=0; j<=start; j++)
354 543 : pos_start=name.find("_",pos_start+1);
355 : std::size_t pos_end=pos_start;
356 1527 : for(unsigned j=0; j<dim; j++)
357 1002 : pos_end=name.find("_",pos_end+1);
358 525 : pos_start++; //do not include heading "_"
359 525 : return name.substr(pos_start,pos_end-pos_start);
360 : };
361 10 : unsigned index_j=ncv_;
362 : unsigned sizeSkip=1;
363 22 : for(int l=pntrToECVsClass_.size()-1; l>=0; l--)
364 : {
365 12 : const unsigned dim_l=pntrToECVsClass_[l]->getNumberOfArguments();
366 12 : index_j-=dim_l;
367 12 : std::vector<std::string> lambdas_l(1);
368 12 : lambdas_l[0]=getLambdaName(deltaF_name_[0],index_j,dim_l);
369 523 : for(unsigned i=sizeSkip; i<deltaF_size_; i+=sizeSkip)
370 : {
371 513 : std::string tmp_lambda=getLambdaName(deltaF_name_[i],index_j,dim_l);
372 513 : if(tmp_lambda==lambdas_l[0])
373 : break;
374 511 : lambdas_l.push_back(tmp_lambda);
375 : }
376 12 : pntrToECVsClass_[l]->initECVs_restart(lambdas_l);
377 12 : sizeSkip*=lambdas_l.size();
378 12 : }
379 10 : plumed_massert(sizeSkip==deltaF_size_,"RESTART - this should not happen");
380 10 : init_linkECVs(); //link ECVs and initializes index_k_
381 10 : log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
382 10 : obs_steps_=0; //avoid initializing again
383 10 : if(stateRestart)
384 : {
385 2 : if(NumParallel_>1)
386 : {
387 2 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
388 : unsigned iter=0;
389 34 : for(unsigned i=start; i<start+deltaF_.size(); i++)
390 32 : deltaF_[iter++]=all_deltaF_[i];
391 : }
392 : }
393 : else //read each step
394 : {
395 8 : counter_=1;
396 : unsigned count_lines=0;
397 8 : ifile.allowIgnoredFields(); //this allows for multiple restart, but without checking for consistency between them!
398 : double time;
399 48 : while(ifile.scanField("time",time)) //only number of lines and last line is important
400 : {
401 : unsigned restart_stride;
402 16 : ifile.scanField("print_stride",restart_stride);
403 16 : ifile.scanField("rct",rct_);
404 16 : if(NumParallel_==1)
405 : {
406 1014 : for(unsigned i=0; i<deltaF_size_; i++)
407 998 : ifile.scanField(deltaF_name_[i],deltaF_[i]);
408 : }
409 : else
410 : {
411 0 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
412 : unsigned iter=0;
413 0 : for(unsigned i=start; i<start+deltaF_.size(); i++)
414 0 : ifile.scanField(deltaF_name_[i],deltaF_[iter++]);
415 : }
416 16 : ifile.scanField();
417 16 : if(count_lines>0)
418 8 : counter_+=restart_stride;
419 16 : count_lines++;
420 : }
421 8 : counter_*=NumWalkers_;
422 8 : log.printf(" successfully read %u lines, up to t=%g\n",count_lines,time);
423 : }
424 10 : ifile.reset(false);
425 10 : ifile.close();
426 : }
427 : else //same behaviour as METAD
428 0 : plumed_merror("RESTART requested, but file '"+restartFileName+"' was not found!\n Set RESTART=NO or provide a restart file");
429 10 : if(NumWalkers_>1) //make sure that all walkers are doing the same thing
430 : {
431 2 : std::vector<unsigned long> all_counter(NumWalkers_);
432 2 : if(comm.Get_rank()==0)
433 2 : multi_sim_comm.Allgather(counter_,all_counter);
434 2 : comm.Bcast(all_counter,0);
435 : bool same_number_of_steps=true;
436 4 : for(unsigned w=1; w<NumWalkers_; w++)
437 2 : if(all_counter[0]!=all_counter[w])
438 : same_number_of_steps=false;
439 2 : plumed_massert(same_number_of_steps,"RESTART - not all walkers are reading the same file!");
440 : }
441 10 : }
442 20 : else if(restartFileName.length()>0)
443 0 : log.printf(" +++ WARNING +++ the provided STATE_RFILE will be ignored, since RESTART was not requested\n");
444 :
445 : //sync all walkers to avoid opening files before reding is over (see also METAD)
446 30 : comm.Barrier();
447 30 : if(comm.Get_rank()==0 && walkers_mpi)
448 4 : multi_sim_comm.Barrier();
449 :
450 : //setup DeltaFs file
451 30 : deltaFsOfile_.link(*this);
452 30 : if(NumWalkers_>1)
453 : {
454 4 : if(walker_rank_>0)
455 : deltaFsFileName="/dev/null"; //only first walker writes on file
456 8 : deltaFsOfile_.enforceSuffix("");
457 : }
458 30 : deltaFsOfile_.open(deltaFsFileName);
459 30 : if(fmt.length()>0)
460 60 : deltaFsOfile_.fmtField(" "+fmt);
461 : deltaFsOfile_.setHeavyFlush(); //do I need it?
462 30 : deltaFsOfile_.addConstantField("print_stride");
463 30 : deltaFsOfile_.printField("print_stride",print_stride_);
464 :
465 : //open file for storing state
466 30 : if(wStateStride_!=0)
467 : {
468 6 : stateOfile_.link(*this);
469 6 : if(NumWalkers_>1)
470 : {
471 0 : if(walker_rank_>0)
472 : stateFileName="/dev/null"; //only first walker writes on file
473 0 : stateOfile_.enforceSuffix("");
474 : }
475 6 : stateOfile_.open(stateFileName);
476 6 : if(fmt.length()>0)
477 12 : stateOfile_.fmtField(" "+fmt);
478 : }
479 :
480 : //add output components
481 30 : if(calc_work_)
482 : {
483 6 : addComponent("work");
484 12 : componentIsNotPeriodic("work");
485 : }
486 :
487 : //printing some info
488 30 : log.printf(" updating the bias with PACE = %u\n",stride_);
489 30 : log.printf(" initial unbiased OBSERVATION_STEPS = %u (in units of PACE)\n",obs_steps_);
490 30 : if(wStateStride_>0)
491 5 : log.printf(" state checkpoints are written on file %s every %d MD steps\n",stateFileName.c_str(),wStateStride_);
492 30 : if(wStateStride_==-1)
493 1 : log.printf(" state checkpoints are written on file '%s' only on CPT events (or never if MD code does define them!)\n",stateFileName.c_str());
494 30 : if(walkers_mpi)
495 4 : log.printf(" -- WALKERS_MPI: if multiple replicas are present, they will share the same bias via MPI\n");
496 30 : if(NumWalkers_>1)
497 : {
498 4 : log.printf(" using multiple walkers\n");
499 4 : log.printf(" number of walkers: %u\n",NumWalkers_);
500 4 : log.printf(" walker rank: %u\n",walker_rank_);
501 : }
502 30 : int mw_warning=0;
503 30 : if(!walkers_mpi && comm.Get_rank()==0 && multi_sim_comm.Get_size()>(int)NumWalkers_)
504 0 : mw_warning=1;
505 30 : comm.Bcast(mw_warning,0);
506 30 : if(mw_warning) //log.printf messes up with comm, so never use it without Bcast!
507 0 : log.printf(" +++ WARNING +++ multiple replicas will NOT communicate unless the flag WALKERS_MPI is used\n");
508 30 : if(NumParallel_>1)
509 2 : log.printf(" using multiple MPI processes per simulation: %u\n",NumParallel_);
510 30 : if(NumOMP_>1)
511 25 : log.printf(" using multiple OpenMP threads per simulation: %u\n",NumOMP_);
512 30 : if(serial)
513 5 : log.printf(" -- SERIAL: no loop parallelization, despite %d MPI processes and %u OpenMP threads available\n",comm.Get_size(),OpenMP::getNumThreads());
514 30 : log.printf(" Bibliography: ");
515 60 : log<<plumed.cite("M. Invernizzi, P.M. Piaggi, and M. Parrinello, Phys. Rev. X 10, 041034 (2020)");
516 30 : log.printf("\n");
517 30 : }
518 :
519 1490 : void OPESexpanded::calculate()
520 : {
521 1490 : if(deltaF_size_==0) //no bias before initialization
522 325 : return;
523 :
524 : //get diffMax, to avoid over/underflow
525 1165 : double diffMax=-std::numeric_limits<double>::max();
526 1165 : #pragma omp parallel num_threads(NumOMP_)
527 : {
528 : #pragma omp for reduction(max:diffMax)
529 : for(unsigned i=0; i<deltaF_.size(); i++)
530 : {
531 : diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
532 : if(diff_[i]>diffMax)
533 : diffMax=diff_[i];
534 : }
535 : }
536 1165 : if(NumParallel_>1)
537 102 : comm.Max(diffMax);
538 :
539 : //calculate the bias and the forces
540 1165 : double sum=0;
541 1165 : std::vector<double> der_sum_cv(ncv_,0);
542 1165 : if(NumOMP_==1)
543 : {
544 2730 : for(unsigned i=0; i<deltaF_.size(); i++)
545 : {
546 2520 : double add_i=std::exp(diff_[i]-diffMax);
547 2520 : sum+=add_i;
548 : //set derivatives
549 6960 : for(unsigned j=0; j<ncv_; j++)
550 4440 : der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
551 : }
552 : }
553 : else
554 : {
555 955 : #pragma omp parallel num_threads(NumOMP_)
556 : {
557 : std::vector<double> omp_der_sum_cv(ncv_,0);
558 : #pragma omp for reduction(+:sum) nowait
559 : for(unsigned i=0; i<deltaF_.size(); i++)
560 : {
561 : double add_i=std::exp(diff_[i]-diffMax);
562 : sum+=add_i;
563 : //set derivatives
564 : for(unsigned j=0; j<ncv_; j++)
565 : omp_der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
566 : }
567 : #pragma omp critical
568 : for(unsigned j=0; j<ncv_; j++)
569 : der_sum_cv[j]+=omp_der_sum_cv[j];
570 : }
571 : }
572 1165 : if(NumParallel_>1)
573 : { //each MPI process has part of the full deltaF_ vector, so must Sum
574 102 : comm.Sum(sum);
575 102 : comm.Sum(der_sum_cv);
576 : }
577 :
578 : //set bias and forces
579 1165 : const double bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
580 : setBias(bias);
581 3163 : for(unsigned j=0; j<ncv_; j++)
582 1998 : setOutputForce(j,kbt_*der_sum_cv[j]/sum);
583 : }
584 :
585 1490 : void OPESexpanded::update()
586 : {
587 1490 : if(isFirstStep_) //skip very first step, as in METAD
588 : {
589 30 : isFirstStep_=false;
590 30 : if(obs_steps_!=1) //if obs_steps_==1 go on with initialization
591 : return;
592 : }
593 1464 : if(getStep()%stride_==0)
594 : {
595 739 : if(obs_steps_>0)
596 : {
597 463 : for(unsigned j=0; j<ncv_; j++)
598 304 : obs_cvs_[counter_*ncv_+j]=getArgument(j);
599 159 : counter_++;
600 159 : if(counter_==obs_steps_)
601 : {
602 20 : log.printf("\nAction %s\n",getName().c_str());
603 20 : init_fromObs();
604 20 : log.printf("Finished initialization\n\n");
605 20 : counter_=NumWalkers_; //all preliminary observations count 1
606 20 : obs_steps_=0; //no more observation
607 : }
608 159 : return;
609 : }
610 :
611 : //update averages
612 580 : const double current_bias=getOutputQuantity(0); //the first value is always the bias
613 580 : if(NumWalkers_==1)
614 500 : updateDeltaF(current_bias);
615 : else
616 : {
617 80 : std::vector<double> cvs(ncv_);
618 240 : for(unsigned j=0; j<ncv_; j++)
619 160 : cvs[j]=getArgument(j);
620 80 : std::vector<double> all_bias(NumWalkers_);
621 80 : std::vector<double> all_cvs(NumWalkers_*ncv_);
622 80 : if(comm.Get_rank()==0)
623 : {
624 80 : multi_sim_comm.Allgather(current_bias,all_bias);
625 80 : multi_sim_comm.Allgather(cvs,all_cvs);
626 : }
627 80 : comm.Bcast(all_bias,0);
628 80 : comm.Bcast(all_cvs,0);
629 240 : for(unsigned w=0; w<NumWalkers_; w++)
630 : {
631 : //calculate ECVs
632 160 : unsigned index_wj=w*ncv_;
633 380 : for(unsigned k=0; k<pntrToECVsClass_.size(); k++)
634 : {
635 220 : pntrToECVsClass_[k]->calculateECVs(&all_cvs[index_wj]);
636 220 : index_wj+=pntrToECVsClass_[k]->getNumberOfArguments();
637 : }
638 160 : updateDeltaF(all_bias[w]);
639 : }
640 : }
641 :
642 : //write DeltaFs to file
643 580 : if((counter_/NumWalkers_-1)%print_stride_==0)
644 44 : printDeltaF();
645 :
646 : //calculate work if requested
647 580 : if(calc_work_)
648 : { //some copy and paste from calculate()
649 : //get diffMax, to avoid over/underflow
650 110 : double diffMax=-std::numeric_limits<double>::max();
651 110 : #pragma omp parallel num_threads(NumOMP_)
652 : {
653 : #pragma omp for reduction(max:diffMax)
654 : for(unsigned i=0; i<deltaF_.size(); i++)
655 : {
656 : diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
657 : if(diff_[i]>diffMax)
658 : diffMax=diff_[i];
659 : }
660 : }
661 110 : if(NumParallel_>1)
662 50 : comm.Max(diffMax);
663 : //calculate the bias
664 110 : double sum=0;
665 110 : #pragma omp parallel num_threads(NumOMP_)
666 : {
667 : #pragma omp for reduction(+:sum) nowait
668 : for(unsigned i=0; i<deltaF_.size(); i++)
669 : sum+=std::exp(diff_[i]-diffMax);
670 : }
671 110 : if(NumParallel_>1)
672 50 : comm.Sum(sum);
673 110 : const double new_bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
674 : //accumulate work
675 110 : work_+=new_bias-current_bias;
676 220 : getPntrToComponent("work")->set(work_);
677 : }
678 : }
679 :
680 : //dump state if requested
681 1305 : if( (wStateStride_>0 && getStep()%wStateStride_==0) || (wStateStride_==-1 && getCPT()) )
682 11 : dumpStateToFile();
683 : }
684 :
685 30 : void OPESexpanded::init_pntrToECVsClass()
686 : {
687 30 : std::vector<opes::ExpansionCVs*> all_pntrToECVsClass=plumed.getActionSet().select<opes::ExpansionCVs*>();
688 30 : plumed_massert(all_pntrToECVsClass.size()>0,"no Expansion CVs found");
689 67 : for(unsigned j=0; j<ncv_; j++)
690 : {
691 74 : std::string error_notECV("all the ARGs of "+getName()+" must be Expansion Collective Variables (ECV)");
692 37 : const unsigned dot_pos=getPntrToArgument(j)->getName().find(".");
693 37 : plumed_massert(dot_pos<getPntrToArgument(j)->getName().size(),error_notECV+", thus contain a dot in the name");
694 37 : unsigned foundECV_l=all_pntrToECVsClass.size();
695 44 : for(unsigned l=0; l<all_pntrToECVsClass.size(); l++)
696 : {
697 44 : if(getPntrToArgument(j)->getName().substr(0,dot_pos)==all_pntrToECVsClass[l]->getLabel())
698 : {
699 : foundECV_l=l;
700 37 : pntrToECVsClass_.push_back(all_pntrToECVsClass[l]);
701 37 : std::string missing_arg="some ECV component is missing from ARG";
702 37 : plumed_massert(j+all_pntrToECVsClass[l]->getNumberOfArguments()<=getNumberOfArguments(),missing_arg);
703 90 : for(unsigned h=0; h<all_pntrToECVsClass[l]->getNumberOfArguments(); h++)
704 : {
705 53 : std::string argName=getPntrToArgument(j+h)->getName();
706 53 : std::string expectedECVname=all_pntrToECVsClass[l]->getComponentsVector()[h];
707 53 : plumed_massert(argName==expectedECVname,missing_arg+", or is in wrong order: given ARG="+argName+" expected ARG="+expectedECVname);
708 : }
709 37 : j+=all_pntrToECVsClass[l]->getNumberOfArguments()-1;
710 : break;
711 : }
712 : }
713 37 : plumed_massert(foundECV_l<all_pntrToECVsClass.size(),error_notECV);
714 : }
715 67 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
716 44 : for(unsigned ll=l+1; ll<pntrToECVsClass_.size(); ll++)
717 7 : plumed_massert(pntrToECVsClass_[l]->getLabel()!=pntrToECVsClass_[ll]->getLabel(),"cannot use same ECV twice");
718 30 : }
719 :
720 30 : void OPESexpanded::init_linkECVs()
721 : {
722 : //TODO It should be possible to make all of this more straightforward (and probably also faster):
723 : // - get rid of index_k_, making it trivial for each ECV
724 : // - store the ECVs_ and derECVs_ vectors here as a contiguous vector, and use pointers in the ECV classes
725 : // Some caveats:
726 : // - ECVmultiThermalBaric has a nontrivial index_k_ to avoid duplicates. use duplicates instead
727 : // - can the ECVs be MPI parallel or it's too complicated?
728 30 : plumed_massert(deltaF_size_>0,"must set deltaF_size_ before calling init_linkECVs()");
729 30 : if(NumParallel_==1)
730 28 : deltaF_.resize(deltaF_size_);
731 : else
732 : {
733 2 : const unsigned extra=(rank_<(deltaF_size_%NumParallel_)?1:0);
734 2 : deltaF_.resize(deltaF_size_/NumParallel_+extra);
735 : //these are used when printing deltaF_ to file
736 2 : all_deltaF_.resize(deltaF_size_);
737 2 : all_size_.resize(NumParallel_,deltaF_size_/NumParallel_);
738 2 : disp_.resize(NumParallel_);
739 4 : for(unsigned r=0; r<NumParallel_-1; r++)
740 : {
741 2 : if(r<deltaF_size_%NumParallel_)
742 0 : all_size_[r]++;
743 2 : disp_[r+1]=disp_[r]+all_size_[r];
744 : }
745 : }
746 30 : diff_.resize(deltaF_.size());
747 30 : ECVs_.resize(ncv_);
748 30 : derECVs_.resize(ncv_);
749 30 : index_k_.resize(deltaF_.size(),std::vector<unsigned>(ncv_));
750 : unsigned index_j=0;
751 30 : unsigned sizeSkip=deltaF_size_;
752 67 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
753 : {
754 37 : std::vector< std::vector<unsigned> > l_index_k(pntrToECVsClass_[l]->getIndex_k());
755 37 : plumed_massert(deltaF_size_%l_index_k.size()==0,"buggy ECV: mismatch between getTotNumECVs() and getIndex_k().size()");
756 37 : plumed_massert(l_index_k[0].size()==pntrToECVsClass_[l]->getNumberOfArguments(),"buggy ECV: mismatch between number of ARG and underlying CVs");
757 37 : sizeSkip/=l_index_k.size();
758 90 : for(unsigned h=0; h<pntrToECVsClass_[l]->getNumberOfArguments(); h++)
759 : {
760 53 : ECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToECVs(h);
761 53 : derECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToDerECVs(h);
762 53 : if(NumParallel_==1)
763 : {
764 45589 : for(unsigned i=0; i<deltaF_size_; i++)
765 45540 : index_k_[i][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
766 : }
767 : else
768 : {
769 4 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
770 : unsigned iter=0;
771 68 : for(unsigned i=start; i<start+deltaF_.size(); i++)
772 64 : index_k_[iter++][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
773 : }
774 : }
775 37 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
776 37 : }
777 30 : plumed_massert(sizeSkip==1,"this should not happen!");
778 30 : }
779 :
780 20 : void OPESexpanded::init_fromObs() //This could probably be faster and/or require less memory...
781 : {
782 : //in case of multiple walkers gather all the statistics
783 20 : if(NumWalkers_>1)
784 : {
785 2 : std::vector<double> all_obs_cv(ncv_*obs_steps_*NumWalkers_);
786 2 : if(comm.Get_rank()==0)
787 2 : multi_sim_comm.Allgather(obs_cvs_,all_obs_cv);
788 2 : comm.Bcast(all_obs_cv,0);
789 2 : obs_cvs_=all_obs_cv; //could this lead to memory issues?
790 2 : obs_steps_*=NumWalkers_;
791 : }
792 :
793 : //initialize ECVs from observations
794 : unsigned index_j=0;
795 20 : deltaF_size_=1;
796 45 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
797 : {
798 25 : pntrToECVsClass_[l]->initECVs_observ(obs_cvs_,ncv_,index_j);
799 25 : deltaF_size_*=pntrToECVsClass_[l]->getTotNumECVs(); //ECVs from different exansions will be combined
800 25 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
801 : }
802 20 : plumed_massert(index_j==getNumberOfArguments(),"mismatch between number of linked CVs and number of ARG");
803 : //link ECVs and initialize index_k_, mapping each deltaF to a single ECVs set
804 20 : init_linkECVs();
805 :
806 : //initialize deltaF_ from obs
807 : //for the first point, t=0, the ECVs are calculated by initECVs_observ, setting also any initial guess
808 : index_j=0;
809 12379 : for(unsigned i=0; i<deltaF_.size(); i++)
810 56923 : for(unsigned j=0; j<ncv_; j++)
811 44564 : deltaF_[i]+=kbt_*ECVs_[j][index_k_[i][j]];
812 179 : for(unsigned t=1; t<obs_steps_; t++) //starts from t=1
813 : {
814 : unsigned index_j=0;
815 383 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
816 : {
817 224 : pntrToECVsClass_[l]->calculateECVs(&obs_cvs_[t*ncv_+index_j]);
818 224 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
819 : }
820 102677 : for(unsigned i=0; i<deltaF_.size(); i++)
821 : {
822 102518 : const double diff_i=(-getExpansion(i)+deltaF_[i]/kbt_-std::log(t));
823 102518 : if(diff_i>0) //save exp from overflow
824 16071 : deltaF_[i]-=kbt_*(diff_i+std::log1p(std::exp(-diff_i))+std::log1p(-1./(1.+t)));
825 : else
826 86447 : deltaF_[i]-=kbt_*(std::log1p(std::exp(diff_i))+std::log1p(-1./(1.+t)));
827 : }
828 : }
829 : obs_cvs_.clear();
830 :
831 : //set deltaF_name_
832 20 : deltaF_name_.resize(deltaF_size_,"DeltaF");
833 20 : unsigned sizeSkip=deltaF_size_;
834 45 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
835 : {
836 25 : std::vector<std::string> lambdas_l=pntrToECVsClass_[l]->getLambdas();
837 25 : plumed_massert(lambdas_l.size()==pntrToECVsClass_[l]->getTotNumECVs(),"buggy ECV: mismatch between getTotNumECVs() and getLambdas().size()");
838 25 : sizeSkip/=lambdas_l.size();
839 22457 : for(unsigned i=0; i<deltaF_size_; i++)
840 44864 : deltaF_name_[i]+="_"+lambdas_l[(i/sizeSkip)%lambdas_l.size()];
841 25 : }
842 :
843 : //print initialization to file
844 20 : log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
845 20 : printDeltaF();
846 20 : }
847 :
848 64 : void OPESexpanded::printDeltaF()
849 : {
850 64 : deltaFsOfile_.printField("time",getTime());
851 64 : deltaFsOfile_.printField("rct",rct_);
852 64 : if(NumParallel_==1)
853 : {
854 23988 : for(unsigned i=0; i<deltaF_.size(); i++)
855 23926 : deltaFsOfile_.printField(deltaF_name_[i],deltaF_[i]);
856 : }
857 : else
858 : {
859 2 : comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
860 66 : for(unsigned i=0; i<deltaF_size_; i++)
861 64 : deltaFsOfile_.printField(deltaF_name_[i],all_deltaF_[i]);
862 : }
863 64 : deltaFsOfile_.printField();
864 64 : }
865 :
866 11 : void OPESexpanded::dumpStateToFile()
867 : {
868 : //rewrite header or rewind file
869 11 : if(storeOldStates_)
870 3 : stateOfile_.clearFields();
871 8 : else if(walker_rank_==0)
872 8 : stateOfile_.rewind();
873 : //define fields
874 11 : stateOfile_.addConstantField("time");
875 11 : stateOfile_.addConstantField("counter");
876 11 : stateOfile_.addConstantField("rct");
877 : //print
878 11 : stateOfile_.printField("time",getTime());
879 11 : stateOfile_.printField("counter",counter_);
880 11 : stateOfile_.printField("rct",rct_);
881 11 : if(NumParallel_>1)
882 0 : comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
883 240 : for(unsigned i=0; i<deltaF_size_; i++)
884 : {
885 : std::size_t pos_start=7; //skip "DeltaF_"
886 687 : for(unsigned j=0; j<ncv_; j++)
887 : {
888 : plumed_dbg_massert(pos_start>6,"not enought _ in deltaF_name_"+std::to_string(i-1)+" string?");
889 458 : const std::size_t pos_end=deltaF_name_[i].find("_",pos_start);
890 916 : stateOfile_.printField(getPntrToArgument(j)->getName()," "+deltaF_name_[i].substr(pos_start,pos_end-pos_start));
891 458 : pos_start=pos_end+1;
892 : }
893 229 : if(NumParallel_==1)
894 458 : stateOfile_.printField("DeltaF",deltaF_[i]);
895 : else
896 0 : stateOfile_.printField("DeltaF",all_deltaF_[i]);
897 229 : stateOfile_.printField();
898 : }
899 : //make sure file is written even if small
900 11 : if(!storeOldStates_)
901 8 : stateOfile_.flush();
902 11 : }
903 :
904 660 : void OPESexpanded::updateDeltaF(double bias)
905 : {
906 : plumed_dbg_massert(counter_>0,"deltaF_ must be initialized");
907 660 : counter_++;
908 660 : const double arg=(bias-rct_)/kbt_-std::log(counter_-1.);
909 : double increment;
910 660 : if(arg>0) //save exp from overflow
911 28 : increment=kbt_*(arg+std::log1p(std::exp(-arg)));
912 : else
913 632 : increment=kbt_*(std::log1p(std::exp(arg)));
914 660 : #pragma omp parallel num_threads(NumOMP_)
915 : {
916 : #pragma omp for
917 : for(unsigned i=0; i<deltaF_.size(); i++)
918 : {
919 : const double diff_i=(-getExpansion(i)+(bias-rct_+deltaF_[i])/kbt_-std::log(counter_-1.));
920 : if(diff_i>0) //save exp from overflow
921 : deltaF_[i]+=increment-kbt_*(diff_i+std::log1p(std::exp(-diff_i)));
922 : else
923 : deltaF_[i]+=increment-kbt_*std::log1p(std::exp(diff_i));
924 : }
925 : }
926 660 : rct_+=increment+kbt_*std::log1p(-1./counter_);
927 660 : }
928 :
929 644469 : double OPESexpanded::getExpansion(unsigned i) const
930 : {
931 : double expansion=0;
932 3003062 : for(unsigned j=0; j<ncv_; j++)
933 2358593 : expansion+=ECVs_[j][index_k_[i][j]]; //the index_k could be trivially guessed for most ECVs, but unfourtunately not all
934 644469 : return expansion;
935 : }
936 :
937 : }
938 : }
|