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