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 32 : Bias::registerKeywords(keys);
183 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");
184 32 : keys.add("compulsory","PACE","how often the bias is updated");
185 32 : keys.add("compulsory","OBSERVATION_STEPS","100","number of unbiased initial PACE steps to collect statistics for initialization");
186 : //DeltaFs and state files
187 32 : 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)");
188 32 : keys.add("compulsory","PRINT_STRIDE","100","stride for printing to DELTAFS file, in units of PACE");
189 32 : keys.add("optional","FMT","specify format for DELTAFS file");
190 32 : keys.add("optional","STATE_RFILE","read from this file the Delta F estimates and all the info needed to RESTART the simulation");
191 32 : keys.add("optional","STATE_WFILE","write to this file the Delta F estimates and all the info needed to RESTART the simulation");
192 32 : 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 32 : keys.addFlag("STORE_STATES",false,"append to STATE_WFILE instead of ovewriting it each time");
194 : //miscellaneous
195 32 : keys.addFlag("CALC_WORK",false,"calculate the total accumulated work done by the bias since last restart");
196 32 : keys.addFlag("WALKERS_MPI",false,"switch on MPI version of multiple walkers");
197 32 : keys.addFlag("SERIAL",false,"perform calculations in serial");
198 32 : keys.use("RESTART");
199 32 : keys.use("UPDATE_FROM");
200 32 : keys.use("UPDATE_UNTIL");
201 :
202 : //output components
203 64 : keys.addOutputComponent("work","CALC_WORK","scalar","total accumulated work done by the bias");
204 32 : }
205 :
206 30 : OPESexpanded::OPESexpanded(const ActionOptions&ao)
207 : : PLUMED_BIAS_INIT(ao)
208 30 : , isFirstStep_(true)
209 30 : , counter_(0)
210 30 : , ncv_(getNumberOfArguments())
211 30 : , deltaF_size_(0)
212 30 : , rct_(0)
213 30 : , work_(0) {
214 : //set pace
215 30 : parse("PACE",stride_);
216 30 : parse("OBSERVATION_STEPS",obs_steps_);
217 30 : plumed_massert(obs_steps_!=0,"minimum is OBSERVATION_STEPS=1");
218 30 : obs_cvs_.resize(obs_steps_*ncv_);
219 :
220 : //deltaFs file
221 : std::string deltaFsFileName;
222 30 : parse("FILE",deltaFsFileName);
223 60 : parse("PRINT_STRIDE",print_stride_);
224 : std::string fmt;
225 60 : parse("FMT",fmt);
226 : //output checkpoint of current state
227 : std::string restartFileName;
228 60 : parse("STATE_RFILE",restartFileName);
229 : std::string stateFileName;
230 30 : parse("STATE_WFILE",stateFileName);
231 30 : wStateStride_=0;
232 30 : parse("STATE_WSTRIDE",wStateStride_);
233 30 : storeOldStates_=false;
234 30 : parseFlag("STORE_STATES",storeOldStates_);
235 30 : if(wStateStride_!=0 || storeOldStates_) {
236 5 : plumed_massert(stateFileName.length()>0,"filename for storing simulation status not specified, use STATE_WFILE");
237 : }
238 30 : if(wStateStride_>0) {
239 5 : plumed_massert(wStateStride_>=(int)stride_,"STATE_WSTRIDE is in units of MD steps, thus should be a multiple of PACE");
240 : }
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 :
245 : //work flag
246 30 : parseFlag("CALC_WORK",calc_work_);
247 :
248 : //multiple walkers //external MW for cp2k not supported, but anyway cp2k cannot put bias on energy!
249 30 : bool walkers_mpi=false;
250 30 : parseFlag("WALKERS_MPI",walkers_mpi);
251 30 : if(walkers_mpi) {
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 4 : NumWalkers_=multi_sim_comm.Get_size();
258 4 : walker_rank_=multi_sim_comm.Get_rank();
259 : }
260 4 : comm.Bcast(NumWalkers_,0); //if each walker has more than one processor update them all
261 4 : comm.Bcast(walker_rank_,0);
262 : } else {
263 26 : NumWalkers_=1;
264 26 : walker_rank_=0;
265 : }
266 :
267 : //parallelization stuff
268 30 : NumOMP_=OpenMP::getNumThreads();
269 30 : NumParallel_=comm.Get_size();
270 30 : rank_=comm.Get_rank();
271 30 : bool serial=false;
272 30 : parseFlag("SERIAL",serial);
273 30 : if(serial) {
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 :
289 : //restart if needed
290 30 : if(getRestart()) {
291 : bool stateRestart=true;
292 10 : if(restartFileName.length()==0) {
293 : stateRestart=false;
294 : restartFileName=deltaFsFileName;
295 : }
296 10 : IFile ifile;
297 10 : ifile.link(*this);
298 10 : if(ifile.FileExist(restartFileName)) {
299 10 : log.printf(" RESTART - make sure all ECVs used are the same as before\n");
300 10 : log.printf(" restarting from: %s\n",restartFileName.c_str());
301 10 : ifile.open(restartFileName);
302 10 : if(stateRestart) { //get all info
303 2 : log.printf(" it should be a STATE file (not a DELTAFS file)\n");
304 : double time; //not used
305 2 : ifile.scanField("time",time);
306 2 : ifile.scanField("counter",counter_);
307 4 : ifile.scanField("rct",rct_);
308 : std::string tmp_lambda;
309 66 : while(ifile.scanField(getPntrToArgument(0)->getName(),tmp_lambda)) {
310 64 : std::string subs="DeltaF_"+tmp_lambda;
311 128 : for(unsigned jj=1; jj<ncv_; jj++) {
312 : tmp_lambda.clear();
313 64 : ifile.scanField(getPntrToArgument(jj)->getName(),tmp_lambda);
314 128 : subs+="_"+tmp_lambda;
315 : }
316 64 : deltaF_name_.push_back(subs);
317 : double tmp_deltaF;
318 64 : ifile.scanField("DeltaF",tmp_deltaF);
319 64 : deltaF_.push_back(tmp_deltaF);
320 64 : ifile.scanField();
321 : tmp_lambda.clear();
322 : }
323 2 : log.printf(" successfully read %lu DeltaF values\n",deltaF_name_.size());
324 2 : if(NumParallel_>1) {
325 2 : all_deltaF_=deltaF_;
326 : }
327 : } else { //get just deltaFs names
328 8 : ifile.scanFieldList(deltaF_name_);
329 8 : plumed_massert(deltaF_name_.size()>=4,"RESTART - fewer than expected FIELDS found in '"+deltaFsFileName+"' file");
330 8 : plumed_massert(deltaF_name_[deltaF_name_.size()-1]=="print_stride","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
331 8 : plumed_massert(deltaF_name_[0]=="time","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
332 8 : plumed_massert(deltaF_name_[1]=="rct","RESTART - coult not find expected FIELDS in '"+deltaFsFileName+"' file");
333 : deltaF_name_.pop_back();
334 : deltaF_name_.erase(deltaF_name_.begin(),deltaF_name_.begin()+2);
335 : std::size_t pos=5; //each name starts with "DeltaF"
336 22 : for(unsigned j=0; j<ncv_; j++) {
337 14 : pos=deltaF_name_[0].find("_",pos+1); //checking only first one, hopefully is enough
338 : }
339 8 : plumed_massert(pos<deltaF_name_[0].length(),"RESTART - fewer '_' than expected in DeltaF fields: did you remove any CV?");
340 8 : pos=deltaF_name_[0].find("_",pos+1);
341 8 : plumed_massert(pos>deltaF_name_[0].length(),"RESTART - more '_' than expected in DeltaF fields: did you add new CV?");
342 : }
343 : //get lambdas, init ECVs and Link them
344 10 : deltaF_size_=deltaF_name_.size();
345 525 : auto getLambdaName=[](const std::string& name,const unsigned start,const unsigned dim) {
346 : std::size_t pos_start=5; //each name starts with "DeltaF"
347 1068 : for(unsigned j=0; j<=start; j++) {
348 543 : pos_start=name.find("_",pos_start+1);
349 : }
350 : std::size_t pos_end=pos_start;
351 1527 : for(unsigned j=0; j<dim; j++) {
352 1002 : pos_end=name.find("_",pos_end+1);
353 : }
354 525 : pos_start++; //do not include heading "_"
355 525 : return name.substr(pos_start,pos_end-pos_start);
356 : };
357 10 : unsigned index_j=ncv_;
358 : unsigned sizeSkip=1;
359 22 : for(int l=pntrToECVsClass_.size()-1; l>=0; l--) {
360 12 : const unsigned dim_l=pntrToECVsClass_[l]->getNumberOfArguments();
361 12 : index_j-=dim_l;
362 12 : std::vector<std::string> lambdas_l(1);
363 12 : lambdas_l[0]=getLambdaName(deltaF_name_[0],index_j,dim_l);
364 523 : for(unsigned i=sizeSkip; i<deltaF_size_; i+=sizeSkip) {
365 513 : std::string tmp_lambda=getLambdaName(deltaF_name_[i],index_j,dim_l);
366 513 : if(tmp_lambda==lambdas_l[0]) {
367 : break;
368 : }
369 511 : lambdas_l.push_back(tmp_lambda);
370 : }
371 12 : pntrToECVsClass_[l]->initECVs_restart(lambdas_l);
372 12 : sizeSkip*=lambdas_l.size();
373 12 : }
374 10 : plumed_massert(sizeSkip==deltaF_size_,"RESTART - this should not happen");
375 10 : init_linkECVs(); //link ECVs and initializes index_k_
376 10 : log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
377 10 : obs_steps_=0; //avoid initializing again
378 10 : if(stateRestart) {
379 2 : if(NumParallel_>1) {
380 2 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
381 : unsigned iter=0;
382 34 : for(unsigned i=start; i<start+deltaF_.size(); i++) {
383 32 : deltaF_[iter++]=all_deltaF_[i];
384 : }
385 : }
386 : } else { //read each step
387 8 : counter_=1;
388 : unsigned count_lines=0;
389 8 : ifile.allowIgnoredFields(); //this allows for multiple restart, but without checking for consistency between them!
390 : double time;
391 48 : while(ifile.scanField("time",time)) { //only number of lines and last line is important
392 : unsigned restart_stride;
393 16 : ifile.scanField("print_stride",restart_stride);
394 16 : ifile.scanField("rct",rct_);
395 16 : if(NumParallel_==1) {
396 1014 : for(unsigned i=0; i<deltaF_size_; i++) {
397 998 : ifile.scanField(deltaF_name_[i],deltaF_[i]);
398 : }
399 : } else {
400 0 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
401 : unsigned iter=0;
402 0 : for(unsigned i=start; i<start+deltaF_.size(); i++) {
403 0 : ifile.scanField(deltaF_name_[i],deltaF_[iter++]);
404 : }
405 : }
406 16 : ifile.scanField();
407 16 : if(count_lines>0) {
408 8 : counter_+=restart_stride;
409 : }
410 16 : count_lines++;
411 : }
412 8 : counter_*=NumWalkers_;
413 8 : log.printf(" successfully read %u lines, up to t=%g\n",count_lines,time);
414 : }
415 10 : ifile.reset(false);
416 10 : ifile.close();
417 : } else { //same behaviour as METAD
418 0 : plumed_merror("RESTART requested, but file '"+restartFileName+"' was not found!\n Set RESTART=NO or provide a restart file");
419 : }
420 10 : if(NumWalkers_>1) { //make sure that all walkers are doing the same thing
421 2 : std::vector<unsigned long long> all_counter(NumWalkers_);
422 2 : if(comm.Get_rank()==0) {
423 2 : multi_sim_comm.Allgather(counter_,all_counter);
424 : }
425 2 : comm.Bcast(all_counter,0);
426 : bool same_number_of_steps=true;
427 4 : for(unsigned w=1; w<NumWalkers_; w++)
428 2 : if(all_counter[0]!=all_counter[w]) {
429 : same_number_of_steps=false;
430 : }
431 2 : plumed_massert(same_number_of_steps,"RESTART - not all walkers are reading the same file!");
432 : }
433 30 : } else if(restartFileName.length()>0) {
434 0 : log.printf(" +++ WARNING +++ the provided STATE_RFILE will be ignored, since RESTART was not requested\n");
435 : }
436 :
437 : //sync all walkers to avoid opening files before reding is over (see also METAD)
438 30 : comm.Barrier();
439 30 : if(comm.Get_rank()==0 && walkers_mpi) {
440 4 : multi_sim_comm.Barrier();
441 : }
442 :
443 : //setup DeltaFs file
444 30 : deltaFsOfile_.link(*this);
445 30 : if(NumWalkers_>1) {
446 4 : if(walker_rank_>0) {
447 : deltaFsFileName="/dev/null"; //only first walker writes on file
448 : }
449 8 : deltaFsOfile_.enforceSuffix("");
450 : }
451 30 : deltaFsOfile_.open(deltaFsFileName);
452 30 : if(fmt.length()>0) {
453 60 : deltaFsOfile_.fmtField(" "+fmt);
454 : }
455 : deltaFsOfile_.setHeavyFlush(); //do I need it?
456 30 : deltaFsOfile_.addConstantField("print_stride");
457 30 : deltaFsOfile_.printField("print_stride",print_stride_);
458 :
459 : //open file for storing state
460 30 : if(wStateStride_!=0) {
461 6 : stateOfile_.link(*this);
462 6 : if(NumWalkers_>1) {
463 0 : if(walker_rank_>0) {
464 : stateFileName="/dev/null"; //only first walker writes on file
465 : }
466 0 : stateOfile_.enforceSuffix("");
467 : }
468 6 : stateOfile_.open(stateFileName);
469 6 : if(fmt.length()>0) {
470 12 : stateOfile_.fmtField(" "+fmt);
471 : }
472 : }
473 :
474 : //add output components
475 30 : if(calc_work_) {
476 12 : addComponent("work");
477 12 : componentIsNotPeriodic("work");
478 : }
479 :
480 : //printing some info
481 30 : log.printf(" updating the bias with PACE = %u\n",stride_);
482 30 : log.printf(" initial unbiased OBSERVATION_STEPS = %u (in units of PACE)\n",obs_steps_);
483 30 : if(wStateStride_>0) {
484 5 : log.printf(" state checkpoints are written on file %s every %d MD steps\n",stateFileName.c_str(),wStateStride_);
485 : }
486 30 : if(wStateStride_==-1) {
487 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());
488 : }
489 30 : if(walkers_mpi) {
490 4 : log.printf(" -- WALKERS_MPI: if multiple replicas are present, they will share the same bias via MPI\n");
491 : }
492 30 : if(NumWalkers_>1) {
493 4 : log.printf(" using multiple walkers\n");
494 4 : log.printf(" number of walkers: %u\n",NumWalkers_);
495 4 : log.printf(" walker rank: %u\n",walker_rank_);
496 : }
497 30 : int mw_warning=0;
498 30 : if(!walkers_mpi && comm.Get_rank()==0 && multi_sim_comm.Get_size()>(int)NumWalkers_) {
499 0 : mw_warning=1;
500 : }
501 30 : comm.Bcast(mw_warning,0);
502 30 : if(mw_warning) { //log.printf messes up with comm, so never use it without Bcast!
503 0 : log.printf(" +++ WARNING +++ multiple replicas will NOT communicate unless the flag WALKERS_MPI is used\n");
504 : }
505 30 : if(NumParallel_>1) {
506 2 : log.printf(" using multiple MPI processes per simulation: %u\n",NumParallel_);
507 : }
508 30 : if(NumOMP_>1) {
509 25 : log.printf(" using multiple OpenMP threads per simulation: %u\n",NumOMP_);
510 : }
511 30 : if(serial) {
512 5 : log.printf(" -- SERIAL: no loop parallelization, despite %d MPI processes and %u OpenMP threads available\n",comm.Get_size(),OpenMP::getNumThreads());
513 : }
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 1490 : if(deltaF_size_==0) { //no bias before initialization
521 325 : return;
522 : }
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 : diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
531 : if(diff_[i]>diffMax) {
532 : diffMax=diff_[i];
533 : }
534 : }
535 : }
536 1165 : if(NumParallel_>1) {
537 102 : comm.Max(diffMax);
538 : }
539 :
540 : //calculate the bias and the forces
541 1165 : double sum=0;
542 1165 : std::vector<double> der_sum_cv(ncv_,0);
543 1165 : if(NumOMP_==1) {
544 2730 : for(unsigned i=0; i<deltaF_.size(); i++) {
545 2520 : double add_i=std::exp(diff_[i]-diffMax);
546 2520 : sum+=add_i;
547 : //set derivatives
548 6960 : for(unsigned j=0; j<ncv_; j++) {
549 4440 : der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
550 : }
551 : }
552 : } else {
553 955 : #pragma omp parallel num_threads(NumOMP_)
554 : {
555 : std::vector<double> omp_der_sum_cv(ncv_,0);
556 : #pragma omp for reduction(+:sum) nowait
557 : for(unsigned i=0; i<deltaF_.size(); i++) {
558 : double add_i=std::exp(diff_[i]-diffMax);
559 : sum+=add_i;
560 : //set derivatives
561 : for(unsigned j=0; j<ncv_; j++) {
562 : omp_der_sum_cv[j]-=derECVs_[j][index_k_[i][j]]*add_i;
563 : }
564 : }
565 : #pragma omp critical
566 : for(unsigned j=0; j<ncv_; j++) {
567 : der_sum_cv[j]+=omp_der_sum_cv[j];
568 : }
569 : }
570 : }
571 1165 : if(NumParallel_>1) {
572 : //each MPI process has part of the full deltaF_ vector, so must Sum
573 102 : comm.Sum(sum);
574 102 : comm.Sum(der_sum_cv);
575 : }
576 :
577 : //set bias and forces
578 1165 : const double bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
579 1165 : setBias(bias);
580 3163 : for(unsigned j=0; j<ncv_; j++) {
581 1998 : setOutputForce(j,kbt_*der_sum_cv[j]/sum);
582 : }
583 : }
584 :
585 1490 : void OPESexpanded::update() {
586 1490 : if(isFirstStep_) { //skip very first step, as in METAD
587 30 : isFirstStep_=false;
588 30 : if(obs_steps_!=1) { //if obs_steps_==1 go on with initialization
589 : return;
590 : }
591 : }
592 1464 : if(getStep()%stride_==0) {
593 739 : if(obs_steps_>0) {
594 463 : for(unsigned j=0; j<ncv_; j++) {
595 304 : obs_cvs_[counter_*ncv_+j]=getArgument(j);
596 : }
597 159 : counter_++;
598 159 : if(counter_==obs_steps_) {
599 20 : log.printf("\nAction %s\n",getName().c_str());
600 20 : init_fromObs();
601 20 : log.printf("Finished initialization\n\n");
602 20 : counter_=NumWalkers_; //all preliminary observations count 1
603 20 : obs_steps_=0; //no more observation
604 : }
605 159 : return;
606 : }
607 :
608 : //update averages
609 580 : const double current_bias=getOutputQuantity(0); //the first value is always the bias
610 580 : if(NumWalkers_==1) {
611 500 : updateDeltaF(current_bias);
612 : } else {
613 80 : std::vector<double> cvs(ncv_);
614 240 : for(unsigned j=0; j<ncv_; j++) {
615 160 : cvs[j]=getArgument(j);
616 : }
617 80 : std::vector<double> all_bias(NumWalkers_);
618 80 : std::vector<double> all_cvs(NumWalkers_*ncv_);
619 80 : if(comm.Get_rank()==0) {
620 80 : multi_sim_comm.Allgather(current_bias,all_bias);
621 80 : multi_sim_comm.Allgather(cvs,all_cvs);
622 : }
623 80 : comm.Bcast(all_bias,0);
624 80 : comm.Bcast(all_cvs,0);
625 240 : for(unsigned w=0; w<NumWalkers_; w++) {
626 : //calculate ECVs
627 160 : unsigned index_wj=w*ncv_;
628 380 : for(unsigned k=0; k<pntrToECVsClass_.size(); k++) {
629 220 : pntrToECVsClass_[k]->calculateECVs(&all_cvs[index_wj]);
630 220 : index_wj+=pntrToECVsClass_[k]->getNumberOfArguments();
631 : }
632 160 : updateDeltaF(all_bias[w]);
633 : }
634 : }
635 :
636 : //write DeltaFs to file
637 580 : if((counter_/NumWalkers_-1)%print_stride_==0) {
638 44 : printDeltaF();
639 : }
640 :
641 : //calculate work if requested
642 580 : if(calc_work_) {
643 : //some copy and paste from calculate()
644 : //get diffMax, to avoid over/underflow
645 110 : double diffMax=-std::numeric_limits<double>::max();
646 110 : #pragma omp parallel num_threads(NumOMP_)
647 : {
648 : #pragma omp for reduction(max:diffMax)
649 : for(unsigned i=0; i<deltaF_.size(); i++) {
650 : diff_[i]=(-getExpansion(i)+deltaF_[i]/kbt_);
651 : if(diff_[i]>diffMax) {
652 : diffMax=diff_[i];
653 : }
654 : }
655 : }
656 110 : if(NumParallel_>1) {
657 50 : comm.Max(diffMax);
658 : }
659 : //calculate the bias
660 110 : double sum=0;
661 110 : #pragma omp parallel num_threads(NumOMP_)
662 : {
663 : #pragma omp for reduction(+:sum) nowait
664 : for(unsigned i=0; i<deltaF_.size(); i++) {
665 : sum+=std::exp(diff_[i]-diffMax);
666 : }
667 : }
668 110 : if(NumParallel_>1) {
669 50 : comm.Sum(sum);
670 : }
671 110 : const double new_bias=-kbt_*(diffMax+std::log(sum/deltaF_size_));
672 : //accumulate work
673 110 : work_+=new_bias-current_bias;
674 220 : getPntrToComponent("work")->set(work_);
675 : }
676 : }
677 :
678 : //dump state if requested
679 1305 : if( (wStateStride_>0 && getStep()%wStateStride_==0) || (wStateStride_==-1 && getCPT()) ) {
680 11 : dumpStateToFile();
681 : }
682 : }
683 :
684 30 : void OPESexpanded::init_pntrToECVsClass() {
685 30 : std::vector<opes::ExpansionCVs*> all_pntrToECVsClass=plumed.getActionSet().select<opes::ExpansionCVs*>();
686 30 : plumed_massert(all_pntrToECVsClass.size()>0,"no Expansion CVs found");
687 67 : for(unsigned j=0; j<ncv_; j++) {
688 74 : std::string error_notECV("all the ARGs of "+getName()+" must be Expansion Collective Variables (ECV)");
689 37 : const unsigned dot_pos=getPntrToArgument(j)->getName().find(".");
690 37 : plumed_massert(dot_pos<getPntrToArgument(j)->getName().size(),error_notECV+", thus contain a dot in the name");
691 37 : unsigned foundECV_l=all_pntrToECVsClass.size();
692 44 : for(unsigned l=0; l<all_pntrToECVsClass.size(); l++) {
693 44 : if(getPntrToArgument(j)->getName().substr(0,dot_pos)==all_pntrToECVsClass[l]->getLabel()) {
694 : foundECV_l=l;
695 37 : pntrToECVsClass_.push_back(all_pntrToECVsClass[l]);
696 37 : std::string missing_arg="some ECV component is missing from ARG";
697 37 : plumed_massert(j+all_pntrToECVsClass[l]->getNumberOfArguments()<=getNumberOfArguments(),missing_arg);
698 90 : for(unsigned h=0; h<all_pntrToECVsClass[l]->getNumberOfArguments(); h++) {
699 53 : std::string argName=getPntrToArgument(j+h)->getName();
700 53 : std::string expectedECVname=all_pntrToECVsClass[l]->getComponentsVector()[h];
701 53 : plumed_massert(argName==expectedECVname,missing_arg+", or is in wrong order: given ARG="+argName+" expected ARG="+expectedECVname);
702 : }
703 37 : j+=all_pntrToECVsClass[l]->getNumberOfArguments()-1;
704 : break;
705 : }
706 : }
707 37 : plumed_massert(foundECV_l<all_pntrToECVsClass.size(),error_notECV);
708 : }
709 67 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++)
710 44 : for(unsigned ll=l+1; ll<pntrToECVsClass_.size(); ll++) {
711 7 : plumed_massert(pntrToECVsClass_[l]->getLabel()!=pntrToECVsClass_[ll]->getLabel(),"cannot use same ECV twice");
712 : }
713 30 : }
714 :
715 30 : void OPESexpanded::init_linkECVs() {
716 : //TODO It should be possible to make all of this more straightforward (and probably also faster):
717 : // - get rid of index_k_, making it trivial for each ECV
718 : // - store the ECVs_ and derECVs_ vectors here as a contiguous vector, and use pointers in the ECV classes
719 : // Some caveats:
720 : // - ECVmultiThermalBaric has a nontrivial index_k_ to avoid duplicates. use duplicates instead
721 : // - can the ECVs be MPI parallel or it's too complicated?
722 30 : plumed_massert(deltaF_size_>0,"must set deltaF_size_ before calling init_linkECVs()");
723 30 : if(NumParallel_==1) {
724 28 : deltaF_.resize(deltaF_size_);
725 : } else {
726 2 : const unsigned extra=(rank_<(deltaF_size_%NumParallel_)?1:0);
727 2 : deltaF_.resize(deltaF_size_/NumParallel_+extra);
728 : //these are used when printing deltaF_ to file
729 2 : all_deltaF_.resize(deltaF_size_);
730 2 : all_size_.resize(NumParallel_,deltaF_size_/NumParallel_);
731 2 : disp_.resize(NumParallel_);
732 4 : for(unsigned r=0; r<NumParallel_-1; r++) {
733 2 : if(r<deltaF_size_%NumParallel_) {
734 0 : all_size_[r]++;
735 : }
736 2 : disp_[r+1]=disp_[r]+all_size_[r];
737 : }
738 : }
739 30 : diff_.resize(deltaF_.size());
740 30 : ECVs_.resize(ncv_);
741 30 : derECVs_.resize(ncv_);
742 30 : index_k_.resize(deltaF_.size(),std::vector<unsigned>(ncv_));
743 : unsigned index_j=0;
744 30 : unsigned sizeSkip=deltaF_size_;
745 67 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++) {
746 37 : std::vector< std::vector<unsigned> > l_index_k(pntrToECVsClass_[l]->getIndex_k());
747 37 : plumed_massert(deltaF_size_%l_index_k.size()==0,"buggy ECV: mismatch between getTotNumECVs() and getIndex_k().size()");
748 37 : plumed_massert(l_index_k[0].size()==pntrToECVsClass_[l]->getNumberOfArguments(),"buggy ECV: mismatch between number of ARG and underlying CVs");
749 37 : sizeSkip/=l_index_k.size();
750 90 : for(unsigned h=0; h<pntrToECVsClass_[l]->getNumberOfArguments(); h++) {
751 53 : ECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToECVs(h);
752 53 : derECVs_[index_j+h]=pntrToECVsClass_[l]->getPntrToDerECVs(h);
753 53 : if(NumParallel_==1) {
754 45589 : for(unsigned i=0; i<deltaF_size_; i++) {
755 45540 : index_k_[i][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
756 : }
757 : } else {
758 4 : const unsigned start=(deltaF_size_/NumParallel_)*rank_+std::min(rank_,deltaF_size_%NumParallel_);
759 : unsigned iter=0;
760 68 : for(unsigned i=start; i<start+deltaF_.size(); i++) {
761 64 : index_k_[iter++][index_j+h]=l_index_k[(i/sizeSkip)%l_index_k.size()][h];
762 : }
763 : }
764 : }
765 37 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
766 37 : }
767 30 : plumed_massert(sizeSkip==1,"this should not happen!");
768 30 : }
769 :
770 20 : void OPESexpanded::init_fromObs() { //This could probably be faster and/or require less memory...
771 : //in case of multiple walkers gather all the statistics
772 20 : if(NumWalkers_>1) {
773 2 : std::vector<double> all_obs_cv(ncv_*obs_steps_*NumWalkers_);
774 2 : if(comm.Get_rank()==0) {
775 2 : multi_sim_comm.Allgather(obs_cvs_,all_obs_cv);
776 : }
777 2 : comm.Bcast(all_obs_cv,0);
778 2 : obs_cvs_=all_obs_cv; //could this lead to memory issues?
779 2 : obs_steps_*=NumWalkers_;
780 : }
781 :
782 : //initialize ECVs from observations
783 : unsigned index_j=0;
784 20 : deltaF_size_=1;
785 45 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++) {
786 25 : pntrToECVsClass_[l]->initECVs_observ(obs_cvs_,ncv_,index_j);
787 25 : deltaF_size_*=pntrToECVsClass_[l]->getTotNumECVs(); //ECVs from different exansions will be combined
788 25 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
789 : }
790 20 : plumed_massert(index_j==getNumberOfArguments(),"mismatch between number of linked CVs and number of ARG");
791 : //link ECVs and initialize index_k_, mapping each deltaF to a single ECVs set
792 20 : init_linkECVs();
793 :
794 : //initialize deltaF_ from obs
795 : //for the first point, t=0, the ECVs are calculated by initECVs_observ, setting also any initial guess
796 : index_j=0;
797 12379 : for(unsigned i=0; i<deltaF_.size(); i++)
798 56923 : for(unsigned j=0; j<ncv_; j++) {
799 44564 : deltaF_[i]+=kbt_*ECVs_[j][index_k_[i][j]];
800 : }
801 179 : for(unsigned t=1; t<obs_steps_; t++) { //starts from t=1
802 : unsigned index_j=0;
803 383 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++) {
804 224 : pntrToECVsClass_[l]->calculateECVs(&obs_cvs_[t*ncv_+index_j]);
805 224 : index_j+=pntrToECVsClass_[l]->getNumberOfArguments();
806 : }
807 102677 : for(unsigned i=0; i<deltaF_.size(); i++) {
808 102518 : const double diff_i=(-getExpansion(i)+deltaF_[i]/kbt_-std::log(t));
809 102518 : if(diff_i>0) { //save exp from overflow
810 16071 : deltaF_[i]-=kbt_*(diff_i+std::log1p(std::exp(-diff_i))+std::log1p(-1./(1.+t)));
811 : } else {
812 86447 : deltaF_[i]-=kbt_*(std::log1p(std::exp(diff_i))+std::log1p(-1./(1.+t)));
813 : }
814 : }
815 : }
816 : obs_cvs_.clear();
817 :
818 : //set deltaF_name_
819 20 : deltaF_name_.resize(deltaF_size_,"DeltaF");
820 20 : unsigned sizeSkip=deltaF_size_;
821 45 : for(unsigned l=0; l<pntrToECVsClass_.size(); l++) {
822 25 : std::vector<std::string> lambdas_l=pntrToECVsClass_[l]->getLambdas();
823 25 : plumed_massert(lambdas_l.size()==pntrToECVsClass_[l]->getTotNumECVs(),"buggy ECV: mismatch between getTotNumECVs() and getLambdas().size()");
824 25 : sizeSkip/=lambdas_l.size();
825 22457 : for(unsigned i=0; i<deltaF_size_; i++) {
826 44864 : deltaF_name_[i]+="_"+lambdas_l[(i/sizeSkip)%lambdas_l.size()];
827 : }
828 25 : }
829 :
830 : //print initialization to file
831 20 : log.printf(" ->%4u DeltaFs in total\n",deltaF_size_);
832 20 : printDeltaF();
833 20 : }
834 :
835 64 : void OPESexpanded::printDeltaF() {
836 64 : deltaFsOfile_.printField("time",getTime());
837 64 : deltaFsOfile_.printField("rct",rct_);
838 64 : if(NumParallel_==1) {
839 23988 : for(unsigned i=0; i<deltaF_.size(); i++) {
840 23926 : deltaFsOfile_.printField(deltaF_name_[i],deltaF_[i]);
841 : }
842 : } else {
843 2 : comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
844 66 : for(unsigned i=0; i<deltaF_size_; i++) {
845 64 : deltaFsOfile_.printField(deltaF_name_[i],all_deltaF_[i]);
846 : }
847 : }
848 64 : deltaFsOfile_.printField();
849 64 : }
850 :
851 11 : void OPESexpanded::dumpStateToFile() {
852 : //rewrite header or rewind file
853 11 : if(storeOldStates_) {
854 3 : stateOfile_.clearFields();
855 8 : } else if(walker_rank_==0) {
856 8 : stateOfile_.rewind();
857 : }
858 : //define fields
859 11 : stateOfile_.addConstantField("time");
860 11 : stateOfile_.addConstantField("counter");
861 11 : stateOfile_.addConstantField("rct");
862 : //print
863 11 : stateOfile_.printField("time",getTime());
864 11 : stateOfile_.printField("counter",counter_);
865 11 : stateOfile_.printField("rct",rct_);
866 11 : if(NumParallel_>1) {
867 0 : comm.Allgatherv(deltaF_,all_deltaF_,&all_size_[0],&disp_[0]); //can we avoid using this big vector?
868 : }
869 240 : for(unsigned i=0; i<deltaF_size_; i++) {
870 : std::size_t pos_start=7; //skip "DeltaF_"
871 687 : for(unsigned j=0; j<ncv_; j++) {
872 : plumed_dbg_massert(pos_start>6,"not enought _ in deltaF_name_"+std::to_string(i-1)+" string?");
873 458 : const std::size_t pos_end=deltaF_name_[i].find("_",pos_start);
874 916 : stateOfile_.printField(getPntrToArgument(j)->getName()," "+deltaF_name_[i].substr(pos_start,pos_end-pos_start));
875 458 : pos_start=pos_end+1;
876 : }
877 229 : if(NumParallel_==1) {
878 458 : stateOfile_.printField("DeltaF",deltaF_[i]);
879 : } else {
880 0 : stateOfile_.printField("DeltaF",all_deltaF_[i]);
881 : }
882 229 : stateOfile_.printField();
883 : }
884 : //make sure file is written even if small
885 11 : if(!storeOldStates_) {
886 8 : stateOfile_.flush();
887 : }
888 11 : }
889 :
890 660 : void OPESexpanded::updateDeltaF(double bias) {
891 : plumed_dbg_massert(counter_>0,"deltaF_ must be initialized");
892 660 : counter_++;
893 660 : const double arg=(bias-rct_)/kbt_-std::log(counter_-1.);
894 : double increment;
895 660 : if(arg>0) { //save exp from overflow
896 28 : increment=kbt_*(arg+std::log1p(std::exp(-arg)));
897 : } else {
898 632 : increment=kbt_*(std::log1p(std::exp(arg)));
899 : }
900 660 : #pragma omp parallel num_threads(NumOMP_)
901 : {
902 : #pragma omp for
903 : for(unsigned i=0; i<deltaF_.size(); i++) {
904 : const double diff_i=(-getExpansion(i)+(bias-rct_+deltaF_[i])/kbt_-std::log(counter_-1.));
905 : if(diff_i>0) { //save exp from overflow
906 : deltaF_[i]+=increment-kbt_*(diff_i+std::log1p(std::exp(-diff_i)));
907 : } else {
908 : deltaF_[i]+=increment-kbt_*std::log1p(std::exp(diff_i));
909 : }
910 : }
911 : }
912 660 : rct_+=increment+kbt_*std::log1p(-1./counter_);
913 660 : }
914 :
915 644469 : double OPESexpanded::getExpansion(unsigned i) const {
916 : double expansion=0;
917 3003062 : for(unsigned j=0; j<ncv_; j++) {
918 2358593 : expansion+=ECVs_[j][index_k_[i][j]]; //the index_k could be trivially guessed for most ECVs, but unfourtunately not all
919 : }
920 644469 : return expansion;
921 : }
922 :
923 : }
924 : }
|