Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2011-2019 The plumed team
3 : (see the PEOPLE file at the root of the distribution for a list of names)
4 :
5 : See http://www.plumed.org for more information.
6 :
7 : This file is part of plumed, version 2.
8 :
9 : plumed is free software: you can redistribute it and/or modify
10 : it under the terms of the GNU Lesser General Public License as published by
11 : the Free Software Foundation, either version 3 of the License, or
12 : (at your option) any later version.
13 :
14 : plumed is distributed in the hope that it will be useful,
15 : but WITHOUT ANY WARRANTY; without even the implied warranty of
16 : MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 : GNU Lesser General Public License for more details.
18 :
19 : You should have received a copy of the GNU Lesser General Public License
20 : along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 : +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 : #include "PlumedMain.h"
23 : #include "ActionAtomistic.h"
24 : #include "ActionPilot.h"
25 : #include "ActionRegister.h"
26 : #include "ActionSet.h"
27 : #include "ActionWithValue.h"
28 : #include "ActionWithVirtualAtom.h"
29 : #include "Atoms.h"
30 : #include "CLToolMain.h"
31 : #include "ExchangePatterns.h"
32 : #include "GREX.h"
33 : #include "config/Config.h"
34 : #include "tools/Citations.h"
35 : #include "tools/Communicator.h"
36 : #include "tools/DLLoader.h"
37 : #include "tools/Exception.h"
38 : #include "tools/IFile.h"
39 : #include "tools/Log.h"
40 : #include "tools/OpenMP.h"
41 : #include "tools/Tools.h"
42 : #include "tools/Stopwatch.h"
43 : #include <cstdlib>
44 : #include <cstring>
45 : #include <set>
46 : #include <unordered_map>
47 :
48 : using namespace std;
49 :
50 : #include "PlumedMainEnum.inc"
51 :
52 : namespace PLMD {
53 :
54 565984 : const std::unordered_map<std::string, int> & plumedMainWordMap() {
55 567586 : static std::unordered_map<std::string, int> word_map;
56 : static bool init=false;
57 565984 : if(!init) {
58 : #include "PlumedMainMap.inc"
59 : }
60 565984 : init=true;
61 565984 : return word_map;
62 : }
63 :
64 2233 : PlumedMain::PlumedMain():
65 2233 : comm(*new Communicator),
66 2233 : multi_sim_comm(*new Communicator),
67 2233 : dlloader(*new DLLoader),
68 : cltool(NULL),
69 2233 : stopwatch(*new Stopwatch),
70 : grex(NULL),
71 : initialized(false),
72 2233 : log(*new Log),
73 2233 : citations(*new Citations),
74 : step(0),
75 : active(false),
76 : endPlumed(false),
77 2233 : atoms(*new Atoms(*this)),
78 2233 : actionSet(*new ActionSet(*this)),
79 : bias(0.0),
80 : work(0.0),
81 2233 : exchangePatterns(*new(ExchangePatterns)),
82 : exchangeStep(false),
83 : restart(false),
84 : doCheckPoint(false),
85 : stopFlag(NULL),
86 : stopNow(false),
87 : novirial(false),
88 31262 : detailedTimers(false)
89 : {
90 2233 : log.link(comm);
91 4466 : log.setLinePrefix("PLUMED: ");
92 2233 : stopwatch.start();
93 2233 : stopwatch.pause();
94 2233 : }
95 :
96 6692 : PlumedMain::~PlumedMain() {
97 2233 : stopwatch.start();
98 2233 : stopwatch.stop();
99 2233 : if(initialized) log<<stopwatch;
100 2233 : delete &exchangePatterns;
101 2233 : delete &actionSet;
102 4466 : delete &citations;
103 2233 : delete &atoms;
104 2233 : delete &log;
105 2233 : if(grex) delete grex;
106 4466 : delete &stopwatch;
107 2233 : if(cltool) delete cltool;
108 2233 : delete &dlloader;
109 2233 : delete &comm;
110 2233 : delete &multi_sim_comm;
111 4459 : }
112 :
113 : /////////////////////////////////////////////////////////////
114 : // MAIN INTERPRETER
115 :
116 : #define CHECK_INIT(ini,word) plumed_massert(ini,"cmd(\"" + word +"\") should be only used after plumed initialization")
117 : #define CHECK_NOTINIT(ini,word) plumed_massert(!(ini),"cmd(\"" + word +"\") should be only used before plumed initialization")
118 : #define CHECK_NOTNULL(val,word) plumed_massert(val,"NULL pointer received in cmd(\"" + word + "\")");
119 :
120 :
121 282992 : void PlumedMain::cmd(const std::string & word,void*val) {
122 :
123 282992 : stopwatch.start();
124 :
125 565984 : std::vector<std::string> words=Tools::getWords(word);
126 282992 : unsigned nw=words.size();
127 282992 : if(nw==0) {
128 : // do nothing
129 : } else {
130 : int iword=-1;
131 : double d;
132 282992 : const auto it=plumedMainWordMap().find(words[0]);
133 282992 : if(it!=plumedMainWordMap().end()) iword=it->second;
134 282991 : switch(iword) {
135 26838 : case cmd_setBox:
136 26838 : CHECK_INIT(initialized,word);
137 26838 : CHECK_NOTNULL(val,word);
138 26838 : atoms.setBox(val);
139 : break;
140 30656 : case cmd_setPositions:
141 30656 : CHECK_INIT(initialized,word);
142 30656 : atoms.setPositions(val);
143 : break;
144 30659 : case cmd_setMasses:
145 30659 : CHECK_INIT(initialized,word);
146 30659 : atoms.setMasses(val);
147 : break;
148 24595 : case cmd_setCharges:
149 24595 : CHECK_INIT(initialized,word);
150 24595 : atoms.setCharges(val);
151 : break;
152 0 : case cmd_setPositionsX:
153 0 : CHECK_INIT(initialized,word);
154 0 : atoms.setPositions(val,0);
155 : break;
156 0 : case cmd_setPositionsY:
157 0 : CHECK_INIT(initialized,word);
158 0 : atoms.setPositions(val,1);
159 : break;
160 0 : case cmd_setPositionsZ:
161 0 : CHECK_INIT(initialized,word);
162 0 : atoms.setPositions(val,2);
163 : break;
164 24658 : case cmd_setVirial:
165 24658 : CHECK_INIT(initialized,word);
166 24658 : CHECK_NOTNULL(val,word);
167 24658 : atoms.setVirial(val);
168 : break;
169 5988 : case cmd_setEnergy:
170 5988 : CHECK_INIT(initialized,word);
171 5988 : CHECK_NOTNULL(val,word);
172 5988 : atoms.setEnergy(val);
173 : break;
174 30656 : case cmd_setForces:
175 30656 : CHECK_INIT(initialized,word);
176 30656 : atoms.setForces(val);
177 : break;
178 0 : case cmd_setForcesX:
179 0 : CHECK_INIT(initialized,word);
180 0 : atoms.setForces(val,0);
181 : break;
182 0 : case cmd_setForcesY:
183 0 : CHECK_INIT(initialized,word);
184 0 : atoms.setForces(val,1);
185 : break;
186 0 : case cmd_setForcesZ:
187 0 : CHECK_INIT(initialized,word);
188 0 : atoms.setForces(val,2);
189 : break;
190 31072 : case cmd_calc:
191 31072 : CHECK_INIT(initialized,word);
192 31072 : calc();
193 : break;
194 10 : case cmd_prepareDependencies:
195 10 : CHECK_INIT(initialized,word);
196 10 : prepareDependencies();
197 : break;
198 10 : case cmd_shareData:
199 10 : CHECK_INIT(initialized,word);
200 10 : shareData();
201 : break;
202 260 : case cmd_prepareCalc:
203 260 : CHECK_INIT(initialized,word);
204 260 : prepareCalc();
205 : break;
206 10 : case cmd_performCalc:
207 10 : CHECK_INIT(initialized,word);
208 10 : performCalc();
209 : break;
210 260 : case cmd_performCalcNoUpdate:
211 260 : CHECK_INIT(initialized,word);
212 260 : performCalcNoUpdate();
213 : break;
214 10 : case cmd_update:
215 10 : CHECK_INIT(initialized,word);
216 10 : update();
217 : break;
218 6074 : case cmd_setStep:
219 6074 : CHECK_INIT(initialized,word);
220 6086 : CHECK_NOTNULL(val,word);
221 6071 : step=(*static_cast<int*>(val));
222 6071 : atoms.startStep();
223 : break;
224 25261 : case cmd_setStepLong:
225 25261 : CHECK_INIT(initialized,word);
226 25261 : CHECK_NOTNULL(val,word);
227 25261 : step=(*static_cast<long int*>(val));
228 25261 : atoms.startStep();
229 : break;
230 : // words used less frequently:
231 1056 : case cmd_setAtomsNlocal:
232 1056 : CHECK_INIT(initialized,word);
233 1056 : CHECK_NOTNULL(val,word);
234 1056 : atoms.setAtomsNlocal(*static_cast<int*>(val));
235 : break;
236 920 : case cmd_setAtomsGatindex:
237 920 : CHECK_INIT(initialized,word);
238 920 : atoms.setAtomsGatindex(static_cast<int*>(val),false);
239 : break;
240 0 : case cmd_setAtomsFGatindex:
241 0 : CHECK_INIT(initialized,word);
242 0 : atoms.setAtomsGatindex(static_cast<int*>(val),true);
243 : break;
244 136 : case cmd_setAtomsContiguous:
245 136 : CHECK_INIT(initialized,word);
246 136 : CHECK_NOTNULL(val,word);
247 136 : atoms.setAtomsContiguous(*static_cast<int*>(val));
248 : break;
249 15 : case cmd_createFullList:
250 15 : CHECK_INIT(initialized,word);
251 15 : CHECK_NOTNULL(val,word);
252 15 : atoms.createFullList(static_cast<int*>(val));
253 : break;
254 15 : case cmd_getFullList:
255 15 : CHECK_INIT(initialized,word);
256 15 : CHECK_NOTNULL(val,word);
257 15 : atoms.getFullList(static_cast<int**>(val));
258 : break;
259 15 : case cmd_clearFullList:
260 15 : CHECK_INIT(initialized,word);
261 15 : atoms.clearFullList();
262 : break;
263 0 : case cmd_read:
264 0 : CHECK_INIT(initialized,word);
265 0 : if(val)readInputFile(static_cast<char*>(val));
266 0 : else readInputFile("plumed.dat");
267 : break;
268 113 : case cmd_readInputLine:
269 113 : CHECK_INIT(initialized,word);
270 113 : CHECK_NOTNULL(val,word);
271 226 : readInputLine(static_cast<char*>(val));
272 74 : break;
273 0 : case cmd_clear:
274 0 : CHECK_INIT(initialized,word);
275 0 : actionSet.clearDelete();
276 : break;
277 5 : case cmd_getApiVersion:
278 5 : CHECK_NOTNULL(val,word);
279 5 : *(static_cast<int*>(val))=5;
280 5 : break;
281 : // commands which can be used only before initialization:
282 635 : case cmd_init:
283 635 : CHECK_NOTINIT(initialized,word);
284 635 : init();
285 : break;
286 556 : case cmd_setRealPrecision:
287 556 : CHECK_NOTINIT(initialized,word);
288 556 : CHECK_NOTNULL(val,word);
289 556 : atoms.setRealPrecision(*static_cast<int*>(val));
290 : break;
291 512 : case cmd_setMDLengthUnits:
292 512 : CHECK_NOTINIT(initialized,word);
293 512 : CHECK_NOTNULL(val,word);
294 512 : atoms.MD2double(val,d);
295 512 : atoms.setMDLengthUnits(d);
296 : break;
297 512 : case cmd_setMDChargeUnits:
298 512 : CHECK_NOTINIT(initialized,word);
299 512 : CHECK_NOTNULL(val,word);
300 512 : atoms.MD2double(val,d);
301 512 : atoms.setMDChargeUnits(d);
302 : break;
303 512 : case cmd_setMDMassUnits:
304 512 : CHECK_NOTINIT(initialized,word);
305 512 : CHECK_NOTNULL(val,word);
306 512 : atoms.MD2double(val,d);
307 512 : atoms.setMDMassUnits(d);
308 : break;
309 37 : case cmd_setMDEnergyUnits:
310 37 : CHECK_NOTINIT(initialized,word);
311 37 : CHECK_NOTNULL(val,word);
312 37 : atoms.MD2double(val,d);
313 37 : atoms.setMDEnergyUnits(d);
314 : break;
315 0 : case cmd_setMDTimeUnits:
316 0 : CHECK_NOTINIT(initialized,word);
317 0 : CHECK_NOTNULL(val,word);
318 0 : atoms.MD2double(val,d);
319 0 : atoms.setMDTimeUnits(d);
320 : break;
321 0 : case cmd_setNaturalUnits:
322 : // set the boltzman constant for MD in natural units (kb=1)
323 : // only needed in LJ codes if the MD is passing temperatures to plumed (so, not yet...)
324 : // use as cmd("setNaturalUnits")
325 0 : CHECK_NOTINIT(initialized,word);
326 0 : atoms.setMDNaturalUnits(true);
327 : break;
328 44 : case cmd_setNoVirial:
329 44 : CHECK_NOTINIT(initialized,word);
330 44 : novirial=true;
331 44 : break;
332 559 : case cmd_setPlumedDat:
333 559 : CHECK_NOTINIT(initialized,word);
334 559 : CHECK_NOTNULL(val,word);
335 559 : plumedDat=static_cast<char*>(val);
336 : break;
337 245 : case cmd_setMPIComm:
338 245 : CHECK_NOTINIT(initialized,word);
339 245 : comm.Set_comm(val);
340 245 : atoms.setDomainDecomposition(comm);
341 : break;
342 0 : case cmd_setMPIFComm:
343 0 : CHECK_NOTINIT(initialized,word);
344 0 : comm.Set_fcomm(val);
345 0 : atoms.setDomainDecomposition(comm);
346 : break;
347 0 : case cmd_setMPImultiSimComm:
348 0 : CHECK_NOTINIT(initialized,word);
349 0 : multi_sim_comm.Set_comm(val);
350 : break;
351 635 : case cmd_setNatoms:
352 635 : CHECK_NOTINIT(initialized,word);
353 635 : CHECK_NOTNULL(val,word);
354 635 : atoms.setNatoms(*static_cast<int*>(val));
355 : break;
356 556 : case cmd_setTimestep:
357 556 : CHECK_NOTINIT(initialized,word);
358 556 : CHECK_NOTNULL(val,word);
359 556 : atoms.setTimeStep(val);
360 : break;
361 : /* ADDED WITH API==2 */
362 43 : case cmd_setKbT:
363 43 : CHECK_NOTINIT(initialized,word);
364 43 : CHECK_NOTNULL(val,word);
365 43 : atoms.setKbT(val);
366 : break;
367 : /* ADDED WITH API==3 */
368 0 : case cmd_setRestart:
369 0 : CHECK_NOTINIT(initialized,word);
370 0 : CHECK_NOTNULL(val,word);
371 0 : if(*static_cast<int*>(val)!=0) restart=true;
372 : break;
373 : /* ADDED WITH API==4 */
374 0 : case cmd_doCheckPoint:
375 0 : CHECK_INIT(initialized,word);
376 0 : CHECK_NOTNULL(val,word);
377 0 : doCheckPoint = false;
378 0 : if(*static_cast<int*>(val)!=0) doCheckPoint = true;
379 : break;
380 : /* STOP API */
381 555 : case cmd_setMDEngine:
382 555 : CHECK_NOTINIT(initialized,word);
383 555 : CHECK_NOTNULL(val,word);
384 555 : MDEngine=static_cast<char*>(val);
385 : break;
386 549 : case cmd_setLog:
387 549 : CHECK_NOTINIT(initialized,word);
388 549 : log.link(static_cast<FILE*>(val));
389 : break;
390 43 : case cmd_setLogFile:
391 43 : CHECK_NOTINIT(initialized,word);
392 43 : CHECK_NOTNULL(val,word);
393 86 : log.open(static_cast<char*>(val));
394 43 : break;
395 : // other commands that should be used after initialization:
396 30971 : case cmd_setStopFlag:
397 30971 : CHECK_INIT(initialized,word);
398 30971 : CHECK_NOTNULL(val,word);
399 30971 : stopFlag=static_cast<int*>(val);
400 30971 : break;
401 0 : case cmd_getExchangesFlag:
402 0 : CHECK_INIT(initialized,word);
403 0 : CHECK_NOTNULL(val,word);
404 0 : exchangePatterns.getFlag((*static_cast<int*>(val)));
405 : break;
406 0 : case cmd_setExchangesSeed:
407 0 : CHECK_INIT(initialized,word);
408 0 : CHECK_NOTNULL(val,word);
409 0 : exchangePatterns.setSeed((*static_cast<int*>(val)));
410 : break;
411 0 : case cmd_setNumberOfReplicas:
412 0 : CHECK_INIT(initialized,word);
413 0 : CHECK_NOTNULL(val,word);
414 0 : exchangePatterns.setNofR((*static_cast<int*>(val)));
415 : break;
416 0 : case cmd_getExchangesList:
417 0 : CHECK_INIT(initialized,word);
418 0 : CHECK_NOTNULL(val,word);
419 0 : exchangePatterns.getList((static_cast<int*>(val)));
420 : break;
421 517 : case cmd_runFinalJobs:
422 517 : CHECK_INIT(initialized,word);
423 517 : runJobsAtEndOfCalculation();
424 : break;
425 0 : case cmd_isEnergyNeeded:
426 0 : CHECK_INIT(initialized,word);
427 0 : CHECK_NOTNULL(val,word);
428 0 : if(atoms.isEnergyNeeded()) *(static_cast<int*>(val))=1;
429 0 : else *(static_cast<int*>(val))=0;
430 : break;
431 275 : case cmd_getBias:
432 275 : CHECK_INIT(initialized,word);
433 275 : CHECK_NOTNULL(val,word);
434 322 : atoms.double2MD(getBias()/(atoms.getMDUnits().getEnergy()/atoms.getUnits().getEnergy()),val);
435 275 : break;
436 0 : case cmd_checkAction:
437 0 : CHECK_NOTNULL(val,word);
438 0 : plumed_assert(nw==2);
439 47 : *(static_cast<int*>(val))=(actionRegister().check(words[1]) ? 1:0);
440 0 : break;
441 797 : case cmd_GREX:
442 797 : if(!grex) grex=new GREX(*this);
443 797 : plumed_massert(grex,"error allocating grex");
444 : {
445 : std::string kk=words[1];
446 2098 : for(unsigned i=2; i<words.size(); i++) kk+=" "+words[i];
447 2391 : grex->cmd(kk.c_str(),val);
448 : }
449 797 : break;
450 5146 : case cmd_CLTool:
451 5146 : CHECK_NOTINIT(initialized,word);
452 5146 : if(!cltool) cltool=new CLToolMain;
453 : {
454 : std::string kk=words[1];
455 10292 : for(unsigned i=2; i<words.size(); i++) kk+=" "+words[i];
456 15438 : cltool->cmd(kk.c_str(),val);
457 : }
458 5146 : break;
459 1 : default:
460 5 : plumed_merror("cannot interpret cmd(\"" + word + "\"). check plumed developers manual to see the available commands.");
461 : break;
462 : }
463 : }
464 282945 : stopwatch.pause();
465 282945 : }
466 :
467 : ////////////////////////////////////////////////////////////////////////
468 :
469 635 : void PlumedMain::init() {
470 : // check that initialization just happens once
471 635 : initialized=true;
472 635 : atoms.init();
473 635 : if(!log.isOpen()) log.link(stdout);
474 635 : log<<"PLUMED is starting\n";
475 1905 : log<<"Version: "<<config::getVersionLong()<<" (git: "<<config::getVersionGit()<<") compiled on " __DATE__ " at " __TIME__ "\n";
476 635 : log<<"Please cite this paper when using PLUMED ";
477 1905 : log<<cite("Tribello, Bonomi, Branduardi, Camilloni, and Bussi, Comput. Phys. Commun. 185, 604 (2014)");
478 635 : log<<"\n";
479 635 : log<<"For further information see the PLUMED web page at http://www.plumed.org\n";
480 1270 : log<<"Root: "<<config::getPlumedRoot()<<"\n";
481 1905 : log<<"For installed feature, see "<<config::getPlumedRoot() + "/src/config/config.txt\n";
482 1270 : log.printf("Molecular dynamics engine: %s\n",MDEngine.c_str());
483 635 : log.printf("Precision of reals: %d\n",atoms.getRealPrecision());
484 635 : log.printf("Running over %d %s\n",comm.Get_size(),(comm.Get_size()>1?"nodes":"node"));
485 635 : log<<"Number of threads: "<<OpenMP::getNumThreads()<<"\n";
486 635 : log<<"Cache line size: "<<OpenMP::getCachelineSize()<<"\n";
487 635 : log.printf("Number of atoms: %d\n",atoms.getNatoms());
488 635 : if(grex) log.printf("GROMACS-like replica exchange is on\n");
489 1270 : log.printf("File suffix: %s\n",getSuffix().c_str());
490 635 : if(plumedDat.length()>0) {
491 1118 : readInputFile(plumedDat);
492 559 : plumedDat="";
493 : }
494 635 : atoms.updateUnits();
495 635 : log.printf("Timestep: %f\n",atoms.getTimeStep());
496 635 : if(atoms.getKbT()>0.0)
497 43 : log.printf("KbT: %f\n",atoms.getKbT());
498 : else {
499 592 : log.printf("KbT has not been set by the MD engine\n");
500 592 : log.printf("It should be set by hand where needed\n");
501 : }
502 635 : log<<"Relevant bibliography:\n";
503 635 : log<<citations;
504 635 : log<<"Please read and cite where appropriate!\n";
505 635 : log<<"Finished setup\n";
506 635 : }
507 :
508 572 : void PlumedMain::readInputFile(std::string str) {
509 572 : plumed_assert(initialized);
510 1144 : log.printf("FILE: %s\n",str.c_str());
511 1144 : IFile ifile;
512 572 : ifile.link(*this);
513 572 : ifile.open(str);
514 572 : ifile.allowNoEOL();
515 572 : std::vector<std::string> words;
516 5678 : while(Tools::getParsedLine(ifile,words) && !endPlumed) readInputWords(words);
517 572 : endPlumed=false;
518 1144 : log.printf("END FILE: %s\n",str.c_str());
519 572 : log.flush();
520 :
521 1144 : pilots=actionSet.select<ActionPilot*>();
522 572 : }
523 :
524 156 : void PlumedMain::readInputLine(const std::string & str) {
525 156 : plumed_assert(initialized);
526 156 : if(str.empty()) return;
527 312 : std::vector<std::string> words=Tools::getWords(str);
528 156 : citations.clear();
529 156 : readInputWords(words);
530 117 : if(!citations.empty()) {
531 4 : log<<"Relevant bibliography:\n";
532 4 : log<<citations;
533 4 : log<<"Please read and cite where appropriate!\n";
534 : }
535 : }
536 :
537 5283 : void PlumedMain::readInputWords(const std::vector<std::string> & words) {
538 5283 : plumed_assert(initialized);
539 5283 : if(words.empty())return;
540 5283 : else if(words[0]=="_SET_SUFFIX") {
541 3 : plumed_assert(words.size()==2);
542 : setSuffix(words[1]);
543 : } else {
544 10560 : std::vector<std::string> interpreted(words);
545 5280 : Tools::interpretLabel(interpreted);
546 10562 : Action* action=actionRegister().create(ActionOptions(*this,interpreted));
547 5243 : if(!action) {
548 : std::string msg;
549 : msg ="ERROR\nI cannot understand line:";
550 30 : for(unsigned i=0; i<interpreted.size(); ++i) msg+=" "+interpreted[i];
551 : msg+="\nMaybe a missing space or a typo?";
552 2 : log << msg;
553 2 : log.flush();
554 6 : plumed_merror(msg);
555 : };
556 5241 : action->checkRead();
557 5241 : actionSet.push_back(action);
558 : };
559 :
560 10488 : pilots=actionSet.select<ActionPilot*>();
561 : }
562 :
563 : ////////////////////////////////////////////////////////////////////////
564 :
565 0 : void PlumedMain::exit(int c) {
566 0 : comm.Abort(c);
567 0 : }
568 :
569 5278 : Log& PlumedMain::getLog() {
570 5278 : return log;
571 : }
572 :
573 31072 : void PlumedMain::calc() {
574 31072 : prepareCalc();
575 31072 : performCalc();
576 31072 : }
577 :
578 31332 : void PlumedMain::prepareCalc() {
579 31332 : prepareDependencies();
580 31332 : shareData();
581 31332 : }
582 :
583 : //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
584 : // here we have the main steps in "calc()"
585 : // they can be called individually, but the standard thing is to
586 : // traverse them in this order:
587 31510 : void PlumedMain::prepareDependencies() {
588 :
589 63020 : stopwatch.start("1 Prepare dependencies");
590 :
591 : // activate all the actions which are on step
592 : // activation is recursive and enables also the dependencies
593 : // before doing that, the prepare() method is called to see if there is some
594 : // new/changed dependency (up to now, only useful for dependences on virtual atoms,
595 : // which can be dynamically changed).
596 :
597 : // First switch off all actions
598 299717 : for(const auto & p : actionSet) {
599 236697 : p->deactivate();
600 : }
601 :
602 : // for optimization, an "active" flag remains false if no action at all is active
603 31510 : active=false;
604 310475 : for(unsigned i=0; i<pilots.size(); ++i) {
605 82485 : if(pilots[i]->onStep()) {
606 73746 : pilots[i]->activate();
607 73746 : active=true;
608 : }
609 : };
610 :
611 : // also, if one of them is the total energy, tell to atoms that energy should be collected
612 299717 : for(const auto & p : actionSet) {
613 473394 : if(p->isActive()) {
614 190362 : if(p->checkNeedsGradients()) p->setOption("GRADIENTS");
615 : }
616 : }
617 :
618 63020 : stopwatch.stop("1 Prepare dependencies");
619 31510 : }
620 :
621 31342 : void PlumedMain::shareData() {
622 : // atom positions are shared (but only if there is something to do)
623 31342 : if(!active)return;
624 59708 : stopwatch.start("2 Sharing data");
625 29854 : if(atoms.getNatoms()>0) atoms.share();
626 59708 : stopwatch.stop("2 Sharing data");
627 : }
628 :
629 260 : void PlumedMain::performCalcNoUpdate() {
630 260 : waitData();
631 260 : justCalculate();
632 260 : backwardPropagate();
633 260 : }
634 :
635 31082 : void PlumedMain::performCalc() {
636 31082 : waitData();
637 31082 : justCalculate();
638 31082 : backwardPropagate();
639 31082 : update();
640 31082 : }
641 :
642 31426 : void PlumedMain::waitData() {
643 31426 : if(!active)return;
644 59876 : stopwatch.start("3 Waiting for data");
645 29938 : if(atoms.getNatoms()>0) atoms.wait();
646 59876 : stopwatch.stop("3 Waiting for data");
647 : }
648 :
649 31426 : void PlumedMain::justCalculate() {
650 31426 : if(!active)return;
651 59876 : stopwatch.start("4 Calculating (forward loop)");
652 29938 : bias=0.0;
653 29938 : work=0.0;
654 :
655 : int iaction=0;
656 : // calculate the active actions in order (assuming *backward* dependence)
657 287384 : for(const auto & p : actionSet) {
658 455016 : if(p->isActive()) {
659 : std::string actionNumberLabel;
660 189634 : if(detailedTimers) {
661 0 : Tools::convert(iaction,actionNumberLabel);
662 0 : actionNumberLabel="4A "+actionNumberLabel+" "+p->getLabel();
663 0 : stopwatch.start(actionNumberLabel);
664 : }
665 189634 : ActionWithValue*av=dynamic_cast<ActionWithValue*>(p);
666 189634 : ActionAtomistic*aa=dynamic_cast<ActionAtomistic*>(p);
667 : {
668 189634 : if(av) av->clearInputForces();
669 189634 : if(av) av->clearDerivatives();
670 : }
671 : {
672 189634 : if(aa) aa->clearOutputForces();
673 282404 : if(aa) if(aa->isActive()) aa->retrieveAtoms();
674 : }
675 189634 : if(p->checkNumericalDerivatives()) p->calculateNumericalDerivatives();
676 189092 : else p->calculate();
677 : // This retrieves components called bias
678 337152 : if(av) bias+=av->getOutputQuantity("bias");
679 337152 : if(av) work+=av->getOutputQuantity("work");
680 189634 : if(av)av->setGradientsIfNeeded();
681 189634 : ActionWithVirtualAtom*avv=dynamic_cast<ActionWithVirtualAtom*>(p);
682 189634 : if(avv)avv->setGradientsIfNeeded();
683 189634 : if(detailedTimers) stopwatch.stop(actionNumberLabel);
684 : }
685 227508 : iaction++;
686 : }
687 59876 : stopwatch.stop("4 Calculating (forward loop)");
688 : }
689 :
690 0 : void PlumedMain::justApply() {
691 0 : backwardPropagate();
692 0 : update();
693 0 : }
694 :
695 31342 : void PlumedMain::backwardPropagate() {
696 31342 : if(!active)return;
697 : int iaction=0;
698 59708 : stopwatch.start("5 Applying (backward loop)");
699 : // apply them in reverse order
700 542934 : for(auto pp=actionSet.rbegin(); pp!=actionSet.rend(); ++pp) {
701 : const auto & p(*pp);
702 453372 : if(p->isActive()) {
703 :
704 : std::string actionNumberLabel;
705 188992 : if(detailedTimers) {
706 0 : Tools::convert(iaction,actionNumberLabel);
707 0 : actionNumberLabel="5A "+actionNumberLabel+" "+p->getLabel();
708 0 : stopwatch.start(actionNumberLabel);
709 : }
710 :
711 188992 : p->apply();
712 188992 : ActionAtomistic*a=dynamic_cast<ActionAtomistic*>(p);
713 : // still ActionAtomistic has a special treatment, since they may need to add forces on atoms
714 188992 : if(a) a->applyForces();
715 :
716 188992 : if(detailedTimers) stopwatch.stop(actionNumberLabel);
717 : }
718 226686 : iaction++;
719 : }
720 :
721 : // this is updating the MD copy of the forces
722 29854 : if(detailedTimers) stopwatch.start("5B Update forces");
723 29854 : if(atoms.getNatoms()>0) atoms.updateForces();
724 29854 : if(detailedTimers) stopwatch.stop("5B Update forces");
725 59708 : stopwatch.stop("5 Applying (backward loop)");
726 : }
727 :
728 31092 : void PlumedMain::update() {
729 31092 : if(!active)return;
730 :
731 59208 : stopwatch.start("6 Update");
732 : // update step (for statistics, etc)
733 59208 : updateFlags.push(true);
734 284024 : for(const auto & p : actionSet) {
735 224816 : p->beforeUpdate();
736 636254 : if(p->isActive() && p->checkUpdate() && updateFlagsTop()) p->update();
737 : }
738 59212 : while(!updateFlags.empty()) updateFlags.pop();
739 : if(!updateFlags.empty()) plumed_merror("non matching changes in the update flags");
740 : // Check that no action has told the calculation to stop
741 29604 : if(stopNow) {
742 36 : if(stopFlag) (*stopFlag)=1;
743 0 : else plumed_merror("your md code cannot handle plumed stop events - add a call to plumed.comm(stopFlag,stopCondition)");
744 : }
745 :
746 : // flush by default every 10000 steps
747 : // hopefully will not affect performance
748 : // also if receive checkpointing signal
749 29604 : if(step%10000==0||doCheckPoint) {
750 651 : fflush();
751 651 : log.flush();
752 6937 : for(const auto & p : actionSet) p->fflush();
753 : }
754 59208 : stopwatch.stop("6 Update");
755 : }
756 :
757 2 : void PlumedMain::load(const std::string& ss) {
758 2 : if(DLLoader::installed()) {
759 : std::string s=ss;
760 : size_t n=s.find_last_of(".");
761 3 : std::string extension="";
762 : std::string base=s;
763 6 : if(n!=std::string::npos && n<s.length()-1) extension=s.substr(n+1);
764 6 : if(n!=std::string::npos && n<s.length()) base=s.substr(0,n);
765 2 : if(extension=="cpp") {
766 : // full path command, including environment setup
767 : // this will work even if plumed is not in the execution path or if it has been
768 : // installed with a name different from "plumed"
769 13 : std::string cmd=config::getEnvCommand()+" \""+config::getPlumedRoot()+"\"/scripts/mklib.sh "+s;
770 2 : log<<"Executing: "<<cmd;
771 2 : if(comm.Get_size()>0) log<<" (only on master node)";
772 2 : log<<"\n";
773 4 : if(comm.Get_rank()==0) system(cmd.c_str());
774 2 : comm.Barrier();
775 4 : base="./"+base;
776 : }
777 9 : s=base+"."+config::getSoExt();
778 2 : void *p=dlloader.load(s);
779 2 : if(!p) {
780 4 : const std::string error_msg="I cannot load library " + ss + " " + dlloader.error();
781 1 : log<<"ERROR\n";
782 1 : log<<error_msg<<"\n";
783 3 : plumed_merror(error_msg);
784 : }
785 3 : log<<"Loading shared library "<<s.c_str()<<"\n";
786 1 : log<<"Here is the new list of available actions\n";
787 1 : log<<actionRegister();
788 0 : } else plumed_merror("loading not enabled, please recompile with -D__PLUMED_HAS_DLOPEN");
789 1 : }
790 :
791 848 : double PlumedMain::getBias() const {
792 848 : return bias;
793 : }
794 :
795 450 : double PlumedMain::getWork() const {
796 450 : return work;
797 : }
798 :
799 39 : FILE* PlumedMain::fopen(const char *path, const char *mode) {
800 39 : std::string mmode(mode);
801 39 : std::string ppath(path);
802 : std::string suffix(getSuffix());
803 39 : std::string ppathsuf=ppath+suffix;
804 39 : FILE*fp=std::fopen(const_cast<char*>(ppathsuf.c_str()),const_cast<char*>(mmode.c_str()));
805 39 : if(!fp) fp=std::fopen(const_cast<char*>(ppath.c_str()),const_cast<char*>(mmode.c_str()));
806 39 : plumed_massert(fp,"file " + ppath + " cannot be found");
807 39 : return fp;
808 : }
809 :
810 39 : int PlumedMain::fclose(FILE*fp) {
811 39 : return std::fclose(fp);
812 : }
813 :
814 1924 : std::string PlumedMain::cite(const std::string&item) {
815 1924 : return citations.cite(item);
816 : }
817 :
818 1104 : void PlumedMain::fflush() {
819 3498 : for(const auto & p : files) {
820 2394 : p->flush();
821 : }
822 1104 : }
823 :
824 3338 : void PlumedMain::insertFile(FileBase&f) {
825 6676 : files.insert(&f);
826 3338 : }
827 :
828 3553 : void PlumedMain::eraseFile(FileBase&f) {
829 7106 : files.erase(&f);
830 3553 : }
831 :
832 35 : void PlumedMain::stop() {
833 35 : stopNow=true;
834 35 : }
835 :
836 517 : void PlumedMain::runJobsAtEndOfCalculation() {
837 6010 : for(const auto & p : actionSet) {
838 4976 : p->runFinalJobs();
839 : }
840 517 : }
841 :
842 4839 : }
843 :
844 : //////////////////////////////////////////////////////////////////////////////////////////////////////////////////
|