OpenMD 3.0
Molecular Dynamics in the Open
Loading...
Searching...
No Matches
SimCreator.cpp
Go to the documentation of this file.
1/*
2 * Copyright (c) 2004-present, The University of Notre Dame. All rights
3 * reserved.
4 *
5 * Redistribution and use in source and binary forms, with or without
6 * modification, are permitted provided that the following conditions are met:
7 *
8 * 1. Redistributions of source code must retain the above copyright notice,
9 * this list of conditions and the following disclaimer.
10 *
11 * 2. Redistributions in binary form must reproduce the above copyright notice,
12 * this list of conditions and the following disclaimer in the documentation
13 * and/or other materials provided with the distribution.
14 *
15 * 3. Neither the name of the copyright holder nor the names of its
16 * contributors may be used to endorse or promote products derived from
17 * this software without specific prior written permission.
18 *
19 * THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
20 * AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
21 * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
22 * ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
23 * LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
24 * CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF
25 * SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS
26 * INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN
27 * CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE)
28 * ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
29 * POSSIBILITY OF SUCH DAMAGE.
30 *
31 * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
32 * research, please cite the appropriate papers when you publish your
33 * work. Good starting points are:
34 *
35 * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
36 * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
37 * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 234107 (2008).
38 * [4] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
39 * [5] Kuang & Gezelter, Mol. Phys., 110, 691-701 (2012).
40 * [6] Lamichhane, Gezelter & Newman, J. Chem. Phys. 141, 134109 (2014).
41 * [7] Lamichhane, Newman & Gezelter, J. Chem. Phys. 141, 134110 (2014).
42 * [8] Bhattarai, Newman & Gezelter, Phys. Rev. B 99, 094106 (2019).
43 */
44
45/**
46 * @file SimCreator.cpp
47 * @author tlin
48 * @date 11/03/2004
49 * @version 1.0
50 */
51
52#include "brains/SimCreator.hpp"
53
54#include <exception>
55#include <iostream>
56#include <sstream>
57#include <string>
58
59#ifdef IS_MPI
60#include <mpi.h>
61#endif
62
63#include "antlr/ANTLRException.hpp"
64#include "antlr/CharStreamException.hpp"
65#include "antlr/MismatchedCharException.hpp"
66#include "antlr/MismatchedTokenException.hpp"
67#include "antlr/NoViableAltException.hpp"
68#include "antlr/NoViableAltForCharException.hpp"
69#include "antlr/RecognitionException.hpp"
70#include "antlr/TokenStreamException.hpp"
71#include "antlr/TokenStreamIOException.hpp"
72#include "antlr/TokenStreamRecognitionException.hpp"
73#include "brains/ForceField.hpp"
76#include "io/DumpReader.hpp"
77#include "mdParser/MDLexer.hpp"
78#include "mdParser/MDParser.hpp"
79#include "mdParser/MDTreeParser.hpp"
80#include "mdParser/SimplePreprocessor.hpp"
81#include "types/DirectionalAdapter.hpp"
82#include "types/EAMAdapter.hpp"
83#include "types/FixedChargeAdapter.hpp"
84#include "types/FluctuatingChargeAdapter.hpp"
85#include "types/MultipoleAdapter.hpp"
86#include "types/PolarizableAdapter.hpp"
87#include "types/SuttonChenAdapter.hpp"
88#include "utils/RandNumGen.hpp"
89#include "utils/Revision.hpp"
90#include "utils/Trim.hpp"
91#include "utils/simError.h"
92
93namespace OpenMD {
94
95 Globals* SimCreator::parseFile(std::istream& rawMetaDataStream,
96 const std::string& filename, int mdFileVersion,
97 int startOfMetaDataBlock) {
98 Globals* simParams = NULL;
99 try {
100 // Create a preprocessor that preprocesses md file into an ostringstream
101 std::stringstream ppStream;
102#ifdef IS_MPI
103 int streamSize;
104 const int primaryNode = 0;
105
106 if (worldRank == primaryNode) {
107 MPI_Bcast(&mdFileVersion, 1, MPI_INT, primaryNode, MPI_COMM_WORLD);
108#endif
109 SimplePreprocessor preprocessor;
110 preprocessor.preprocess(rawMetaDataStream, filename,
111 startOfMetaDataBlock, ppStream);
112
113#ifdef IS_MPI
114 // broadcasting the stream size
115 streamSize = ppStream.str().size() + 1;
116 MPI_Bcast(&streamSize, 1, MPI_INT, primaryNode, MPI_COMM_WORLD);
117 MPI_Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())),
118 streamSize, MPI_CHAR, primaryNode, MPI_COMM_WORLD);
119 } else {
120 MPI_Bcast(&mdFileVersion, 1, MPI_INT, primaryNode, MPI_COMM_WORLD);
121
122 // get stream size
123 MPI_Bcast(&streamSize, 1, MPI_INT, primaryNode, MPI_COMM_WORLD);
124 char* buf = new char[streamSize];
125 assert(buf);
126
127 // receive file content
128 MPI_Bcast(buf, streamSize, MPI_CHAR, primaryNode, MPI_COMM_WORLD);
129
130 ppStream.str(buf);
131 delete[] buf;
132 }
133#endif
134 // Create a scanner that reads from the input stream
135 MDLexer lexer(ppStream);
136 lexer.setFilename(filename);
137 lexer.initDeferredLineCount();
138
139 // Create a parser that reads from the scanner
140 MDParser parser(lexer);
141 parser.setFilename(filename);
142
143 // Create an observer that synchorizes file name change
144 FilenameObserver observer;
145 observer.setLexer(&lexer);
146 observer.setParser(&parser);
147 lexer.setObserver(&observer);
148 antlr::ASTFactory factory;
149 parser.initializeASTFactory(factory);
150 parser.setASTFactory(&factory);
151 parser.mdfile();
152 // Create a tree parser that reads information into Globals
153 MDTreeParser treeParser;
154 treeParser.initializeASTFactory(factory);
155 treeParser.setASTFactory(&factory);
156 simParams = treeParser.walkTree(parser.getAST());
157 }
158
160 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
161 "Mismatched Character: %s in file %s at line %d, column %d\n",
162 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
163 e.getColumn());
164 painCave.isFatal = 1;
165 simError();
167 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
168 "Mismatched Token: %s in file %s at line %d, column %d\n",
169 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
170 e.getColumn());
171 painCave.isFatal = 1;
172 simError();
174 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
175 "No Viable Alternative for Character: %s in file %s at line %d, "
176 "column %d\n",
177 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
178 e.getColumn());
179 painCave.isFatal = 1;
180 simError();
181 } catch (antlr::NoViableAltException& e) {
182 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
183 "No Viable Alternative: %s in file %s at line %d, column %d\n",
184 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
185 e.getColumn());
186 painCave.isFatal = 1;
187 simError();
188 }
189
191 snprintf(
192 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
193 "Token Stream Recognition: %s in file %s at line %d, column %d\n",
194 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
195 e.getColumn());
196 painCave.isFatal = 1;
197 simError();
198 }
199
201 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
202 "Token Stream IO exception: %s\n", e.getMessage().c_str());
203 painCave.isFatal = 1;
204 simError();
205 }
206
207 catch (antlr::TokenStreamException& e) {
208 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
209 "Token Stream exception: %s\n", e.getMessage().c_str());
210 painCave.isFatal = 1;
211 simError();
212 } catch (antlr::RecognitionException& e) {
213 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
214 "Recognition exception: %s in file %s at line %d, column %d\n",
215 e.getMessage().c_str(), e.getFilename().c_str(), e.getLine(),
216 e.getColumn());
217 painCave.isFatal = 1;
218 simError();
219 } catch (antlr::CharStreamException& e) {
220 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
221 "Character Stream exception: %s\n", e.getMessage().c_str());
222 painCave.isFatal = 1;
223 simError();
224 } catch (OpenMDException& e) {
225 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH, "%s\n", e.what());
226 painCave.isFatal = 1;
227 simError();
228 } catch (std::exception& e) {
229 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
230 "parser exception: %s\n", e.what());
231 painCave.isFatal = 1;
232 simError();
233 }
234
235 simParams->setMDfileVersion(mdFileVersion);
236 return simParams;
237 }
238
239 SimInfo* SimCreator::createSim(const std::string& mdFileName,
240 bool loadInitCoords) {
241 const int bufferSize = 65535;
242 char buffer[bufferSize];
243 int lineNo = 0;
244 std::string mdRawData;
245 int metaDataBlockStart = -1;
246 int metaDataBlockEnd = -1;
247 streamoff mdOffset {};
248 int mdFileVersion(2);
249
250 // Create a string for embedding the version information in the MetaData
251 std::string version;
252 version.assign("## Last run using OpenMD version: ");
253 version.append(OPENMD_VERSION_MAJOR);
254 version.append(".");
255 version.append(OPENMD_VERSION_MINOR);
256 version.append(",");
257
258 std::string rev(revision, strnlen(revision, 40));
259 rev.append(40 - rev.length(), ' ');
260
261 version.append(" revision: ");
262 // If there's no GIT revision, just call this the RELEASE revision.
263 if (!rev.empty()) {
264 version.append(rev);
265 } else {
266 version.append("RELEASE");
267 }
268
269#ifdef IS_MPI
270 const int primaryNode = 0;
271 if (worldRank == primaryNode) {
272#endif
273
274 std::ifstream mdFile_;
275 mdFile_.open(mdFileName.c_str(), ifstream::in | ifstream::binary);
276
277 if (mdFile_.fail()) {
278 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
279 "SimCreator: Cannot open file: %s\n", mdFileName.c_str());
280 painCave.isFatal = 1;
281 simError();
282 }
283
284 mdFile_.getline(buffer, bufferSize);
285 ++lineNo;
286 std::string line = Utils::trimLeftCopy(buffer);
287 std::size_t i = CaseInsensitiveFind(line, "<OpenMD");
288 if (i == string::npos) {
289 // try the older file strings to see if that works:
290 i = CaseInsensitiveFind(line, "<OOPSE");
291 }
292
293 if (i == string::npos) {
294 // still no luck!
295 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
296 "SimCreator: File: %s is not a valid OpenMD file!\n",
297 mdFileName.c_str());
298 painCave.isFatal = 1;
299 simError();
300 }
301
302 // found the correct opening string, now try to get the file
303 // format version number.
304
305 StringTokenizer tokenizer(line, "=<> \t\n\r");
306 std::string fileType = tokenizer.nextToken();
307 toUpper(fileType);
308
309 mdFileVersion = 0;
310
311 if (fileType == "OPENMD") {
312 while (tokenizer.hasMoreTokens()) {
313 std::string token(tokenizer.nextToken());
314 toUpper(token);
315 if (token == "VERSION") {
316 mdFileVersion = tokenizer.nextTokenAsInt();
317 break;
318 }
319 }
320 }
321
322 bool startFound(false);
323 bool endFound(false);
324 // scan through the input stream and find MetaData tag
325 while (mdFile_.getline(buffer, bufferSize)) {
326 ++lineNo;
327
328 std::string line = Utils::trimLeftCopy(buffer);
329 if (metaDataBlockStart == -1) {
330 std::size_t i = CaseInsensitiveFind(line, "<MetaData>");
331 if (i != string::npos) {
332 metaDataBlockStart = lineNo;
333 startFound = true;
334 mdOffset = mdFile_.tellg();
335 }
336 } else {
337 std::size_t i = CaseInsensitiveFind(line, "</MetaData>");
338 if (i != string::npos) {
339 metaDataBlockEnd = lineNo;
340 endFound = true;
341 }
342 }
343 if (startFound && endFound) break;
344 }
345
346 if (metaDataBlockStart == -1) {
347 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
348 "SimCreator: File: %s did not contain a <MetaData> tag!\n",
349 mdFileName.c_str());
350 painCave.isFatal = 1;
351 simError();
352 }
353 if (metaDataBlockEnd == -1) {
354 snprintf(
355 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
356 "SimCreator: File: %s did not contain a closed MetaData block!\n",
357 mdFileName.c_str());
358 painCave.isFatal = 1;
359 simError();
360 }
361
362 mdFile_.clear();
363 mdFile_.seekg(0);
364 mdFile_.seekg(mdOffset);
365
366 mdRawData.clear();
367
368 bool foundVersion = false;
369
370 for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
371 mdFile_.getline(buffer, bufferSize);
372 std::string line = Utils::trimLeftCopy(buffer);
373 std::size_t j =
374 CaseInsensitiveFind(line, "## Last run using OpenMD Version");
375 if (j != string::npos) {
376 foundVersion = true;
377 mdRawData += version;
378 } else {
379 mdRawData += buffer;
380 }
381 mdRawData += "\n";
382 }
383
384 if (!foundVersion) mdRawData += version + "\n";
385
386 mdFile_.close();
387
388#ifdef IS_MPI
389 }
390#endif
391
392 std::stringstream rawMetaDataStream(mdRawData);
393
394 // parse meta-data file
395 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
396 metaDataBlockStart + 1);
397
398 // create the force field
399 ForceField* ff = new ForceField(simParams->getForceField());
400
401 if (ff == NULL) {
402 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
403 "ForceField Factory can not create %s force field\n",
404 simParams->getForceField().c_str());
405 painCave.isFatal = 1;
406 simError();
407 }
408
409 if (simParams->haveForceFieldFileName()) {
410 ff->setForceFieldFileName(simParams->getForceFieldFileName());
411 }
412
413 std::string forcefieldFileName;
414 forcefieldFileName = ff->getForceFieldFileName();
415
416 if (simParams->haveForceFieldVariant()) {
417 // If the force field has variant, the variant force field name will be
418 // Base.variant.frc. For exampel EAM.u6.frc
419
420 std::string variant = simParams->getForceFieldVariant();
421
422 std::string::size_type pos = forcefieldFileName.rfind(".frc");
423 variant = "." + variant;
424 if (pos != std::string::npos) {
425 forcefieldFileName.insert(pos, variant);
426 } else {
427 // If the default force field file name does not containt .frc suffix,
428 // just append the .variant
429 forcefieldFileName.append(variant);
430 }
431 }
432
433 ff->parse(forcefieldFileName);
434 // create SimInfo
435 SimInfo* info = new SimInfo(ff, simParams);
436
437 info->setRawMetaData(mdRawData);
438
439 // gather parameters (SimCreator only retrieves part of the
440 // parameters)
441 gatherParameters(info, mdFileName);
442
443 // divide the molecules and determine the global index of molecules
444#ifdef IS_MPI
445 divideMolecules(info);
446#endif
447
448 // create the molecules
449 createMolecules(info);
450
451 // find the storage layout
452
453 computeStorageLayouts(info);
454
455 int asl = info->getAtomStorageLayout();
456 int rbsl = info->getRigidBodyStorageLayout();
457 int cgsl = info->getCutoffGroupStorageLayout();
458
459 // allocate memory for DataStorage(circular reference, need to
460 // break it)
461 info->setSnapshotManager(new SimSnapshotManager(info, asl, rbsl, cgsl));
462
463 // set the global index of atoms, rigidbodies and cutoffgroups
464 //(only need to be set once, the global index will never change
465 // again). Local indices of atoms and rigidbodies are already set
466 // by MoleculeCreator class which actually delegates the
467 // responsibility to LocalIndexManager.
468 setGlobalIndex(info);
469
470 // Although addInteractionPairs is called inside SimInfo's addMolecule
471 // method, at that point atoms don't have the global index yet
472 //(their global index are all initialized to -1). Therefore we
473 // have to call addInteractionPairs explicitly here. A way to work
474 // around is that we can determine the beginning global indices of
475 // atoms before they get created.
476 SimInfo::MoleculeIterator mi;
477 Molecule* mol;
478 for (mol = info->beginMolecule(mi); mol != NULL;
479 mol = info->nextMolecule(mi)) {
480 info->addInteractionPairs(mol);
481 }
482
483 if (loadInitCoords) loadCoordinates(info, mdFileName);
484 return info;
485 }
486
487 void SimCreator::gatherParameters(SimInfo* info, const std::string& mdfile) {
488 // figure out the output file names
489 std::string prefix;
490
491#ifdef IS_MPI
492
493 if (worldRank == 0) {
494#endif // is_mpi
495 Globals* simParams = info->getSimParams();
496 if (simParams->haveFinalConfig()) {
497 prefix = getPrefix(simParams->getFinalConfig());
498 } else {
499 prefix = getPrefix(mdfile);
500 }
501
502 info->setFinalConfigFileName(prefix + ".eor");
503 info->setDumpFileName(prefix + ".dump");
504 info->setStatFileName(prefix + ".stat");
505 info->setReportFileName(prefix + ".report");
506 info->setRestFileName(prefix + ".zang");
507
508#ifdef IS_MPI
509 }
510#endif
511 }
512
513#ifdef IS_MPI
514 void SimCreator::divideMolecules(SimInfo* info) {
515 RealType a;
516 int nProcessors;
517 std::vector<int> atomsPerProc;
518 int nGlobalMols = info->getNGlobalMolecules();
519 std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error
520 // condition:
521
522 MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
523
524 if (nProcessors > nGlobalMols) {
525 snprintf(painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
526 "nProcessors (%d) > nMol (%d)\n"
527 "\tThe number of processors is larger than\n"
528 "\tthe number of molecules. This will not result in a \n"
529 "\tusable division of atoms for force decomposition.\n"
530 "\tEither try a smaller number of processors, or run the\n"
531 "\tsingle-processor version of OpenMD.\n",
532 nProcessors, nGlobalMols);
533
534 painCave.isFatal = 1;
535 simError();
536 }
537
538 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
539
540 // initialize atomsPerProc
541 atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
542
543 if (worldRank == 0) {
544 Utils::RandNumGenPtr myRandom = info->getRandomNumberGenerator();
545
546 std::uniform_int_distribution<> processorDistribution {0,
547 nProcessors - 1};
548 std::uniform_real_distribution<RealType> yDistribution {0, 1};
549
550 RealType numerator = info->getNGlobalAtoms();
551 RealType denominator = nProcessors;
552 RealType precast = numerator / denominator;
553 int nTarget = (int)(precast + 0.5);
554 int which_proc {0};
555
556 for (int i = 0; i < nGlobalMols; i++) {
557 // get the molecule stamp first
558 int stampId = info->getMoleculeStampId(i);
559 MoleculeStamp* moleculeStamp = info->getMoleculeStamp(stampId);
560 int add_atoms = moleculeStamp->getNAtoms();
561
562 if (nProcessors > 1) {
563 int done = 0;
564 int loops = 0;
565
566 while (!done) {
567 loops++;
568
569 // Pick a processor at random
570 which_proc = processorDistribution(*myRandom);
571
572 // How many atoms does this processor have so far?
573 int old_atoms = atomsPerProc[which_proc];
574 int new_atoms = old_atoms + add_atoms;
575
576 // If we've been through this loop too many times, we need
577 // to just give up and assign the molecule to this processor
578 // and be done with it.
579
580 if (loops > 100) {
581 snprintf(
582 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
583 "There have been 100 attempts to assign molecule %d to an\n"
584 "\tunderworked processor, but there's no good place to\n"
585 "\tleave it. OpenMD is assigning it at random to processor "
586 "%d.\n",
587 i, which_proc);
588
589 painCave.isFatal = 0;
590 painCave.severity = OPENMD_INFO;
591 simError();
592
593 molToProcMap[i] = which_proc;
594 atomsPerProc[which_proc] += add_atoms;
595
596 done = 1;
597 continue;
598 }
599
600 // If we can add this molecule to this processor without sending
601 // it above nTarget, then go ahead and do it:
602
603 if (new_atoms <= nTarget) {
604 molToProcMap[i] = which_proc;
605 atomsPerProc[which_proc] += add_atoms;
606
607 done = 1;
608 continue;
609 }
610
611 // The only situation left is when new_atoms > nTarget. We
612 // want to accept this with some probability that dies off the
613 // farther we are from nTarget
614
615 // roughly: x = new_atoms - nTarget
616 // Pacc(x) = exp(- a * x)
617 // where a = penalty / (average atoms per molecule)
618
619 RealType x = (RealType)(new_atoms - nTarget);
620 RealType y = yDistribution(*myRandom);
621
622 if (y < exp(-a * x)) {
623 molToProcMap[i] = which_proc;
624 atomsPerProc[which_proc] += add_atoms;
625
626 done = 1;
627 continue;
628 } else {
629 continue;
630 }
631 }
632 } else {
633 which_proc = 0;
634 molToProcMap[i] = which_proc;
635 atomsPerProc[which_proc] += add_atoms;
636 }
637 }
638
639 // Spray out this nonsense to all other processors:
640 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
641
642 } else {
643 // Listen to your marching orders from processor 0:
644 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
645 }
646
647 info->setMolToProcMap(molToProcMap);
648 snprintf(checkPointMsg, MAX_SIM_ERROR_MSG_LENGTH,
649 "Successfully divided the molecules among the processors.\n");
650 errorCheckPoint();
651 }
652
653#endif
654
655 void SimCreator::createMolecules(SimInfo* info) {
656 MoleculeCreator molCreator;
657 int stampId;
658
659 for (int i = 0; i < info->getNGlobalMolecules(); i++) {
660#ifdef IS_MPI
661 if (info->getMolToProc(i) == worldRank) {
662#endif
663
664 stampId = info->getMoleculeStampId(i);
665 Molecule* mol = molCreator.createMolecule(
666 info->getForceField(), info->getMoleculeStamp(stampId), stampId, i,
667 info->getLocalIndexManager());
668
669 info->addMolecule(mol);
670
671#ifdef IS_MPI
672 }
673#endif
674 }
675 }
676
677 void SimCreator::computeStorageLayouts(SimInfo* info) {
678 Globals* simParams = info->getSimParams();
679 int nRigidBodies = info->getNGlobalRigidBodies();
680 AtomTypeSet atomTypes = info->getSimulatedAtomTypes();
681 AtomTypeSet::iterator i;
682 bool hasDirectionalAtoms = false;
683 bool hasFixedCharge = false;
684 bool hasDipoles = false;
685 bool hasQuadrupoles = false;
686 bool hasPolarizable = false;
687 bool hasFluctuatingCharge = false;
688 bool hasMetallic = false;
689 int atomStorageLayout = 0;
690 int rigidBodyStorageLayout = 0;
691 int cutoffGroupStorageLayout = 0;
692
693 atomStorageLayout |= DataStorage::dslPosition;
694 atomStorageLayout |= DataStorage::dslVelocity;
695 atomStorageLayout |= DataStorage::dslForce;
696 cutoffGroupStorageLayout |= DataStorage::dslPosition;
697
698 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
701 EAMAdapter ea = EAMAdapter((*i));
706
707 if (da.isDirectional()) { hasDirectionalAtoms = true; }
708 if (ma.isDipole()) { hasDipoles = true; }
709 if (ma.isQuadrupole()) { hasQuadrupoles = true; }
710 if (ea.isEAM() || sca.isSuttonChen()) { hasMetallic = true; }
711 if (fca.isFixedCharge()) { hasFixedCharge = true; }
712 if (fqa.isFluctuatingCharge()) { hasFluctuatingCharge = true; }
713 if (pa.isPolarizable()) { hasPolarizable = true; }
714 }
715
716 if (nRigidBodies > 0) {
717 rigidBodyStorageLayout |= DataStorage::dslPosition;
718 rigidBodyStorageLayout |= DataStorage::dslVelocity;
719 rigidBodyStorageLayout |= DataStorage::dslForce;
720 rigidBodyStorageLayout |= DataStorage::dslAmat;
721 rigidBodyStorageLayout |= DataStorage::dslAngularMomentum;
722 rigidBodyStorageLayout |= DataStorage::dslTorque;
723 }
724 if (hasDirectionalAtoms) {
725 atomStorageLayout |= DataStorage::dslAmat;
726 if (atomStorageLayout & DataStorage::dslVelocity) {
727 atomStorageLayout |= DataStorage::dslAngularMomentum;
728 }
729 if (atomStorageLayout & DataStorage::dslForce) {
730 atomStorageLayout |= DataStorage::dslTorque;
731 }
732 }
733 if (hasDipoles) { atomStorageLayout |= DataStorage::dslDipole; }
734 if (hasQuadrupoles) { atomStorageLayout |= DataStorage::dslQuadrupole; }
735 if (hasFixedCharge || hasFluctuatingCharge) {
736 atomStorageLayout |= DataStorage::dslSkippedCharge;
737 }
738 if (hasMetallic) {
739 atomStorageLayout |= DataStorage::dslDensity;
740 atomStorageLayout |= DataStorage::dslFunctional;
741 atomStorageLayout |= DataStorage::dslFunctionalDerivative;
742 }
743 if (hasPolarizable) { atomStorageLayout |= DataStorage::dslElectricField; }
744 if (hasFluctuatingCharge) {
745 atomStorageLayout |= DataStorage::dslFlucQPosition;
746 if (atomStorageLayout & DataStorage::dslVelocity) {
747 atomStorageLayout |= DataStorage::dslFlucQVelocity;
748 }
749 if (atomStorageLayout & DataStorage::dslForce) {
750 atomStorageLayout |= DataStorage::dslFlucQForce;
751 }
752 }
753
754 // if the user has asked for them, make sure we've got the memory for the
755 // objects defined.
756
757 if (simParams->getOutputParticlePotential()) {
758 atomStorageLayout |= DataStorage::dslParticlePot;
759 }
760
761 if (simParams->havePrintHeatFlux()) {
762 if (simParams->getPrintHeatFlux()) {
763 atomStorageLayout |= DataStorage::dslParticlePot;
764 }
765 }
766
767 if (simParams->getOutputElectricField() | simParams->haveElectricField() |
768 simParams->haveUniformField() |
769 simParams->haveUniformGradientStrength() |
770 simParams->haveUniformGradientDirection1() |
771 simParams->haveUniformGradientDirection2()) {
772 atomStorageLayout |= DataStorage::dslElectricField;
773 rigidBodyStorageLayout |= DataStorage::dslElectricField;
774 }
775
776 if (simParams->getRNEMDParameters()->haveUseRNEMD()) {
777 if (simParams->getRNEMDParameters()->getUseRNEMD()) {
778 if (simParams->getRNEMDParameters()->requiresElectricField()) {
779 atomStorageLayout |= DataStorage::dslElectricField;
780 rigidBodyStorageLayout |= DataStorage::dslElectricField;
781 }
782 }
783 }
784
785 if (simParams->getOutputSitePotential()) {
786 atomStorageLayout |= DataStorage::dslSitePotential;
787 rigidBodyStorageLayout |= DataStorage::dslSitePotential;
788 }
789
790 if (simParams->getOutputFluctuatingCharges()) {
791 atomStorageLayout |= DataStorage::dslFlucQPosition;
792 atomStorageLayout |= DataStorage::dslFlucQVelocity;
793 atomStorageLayout |= DataStorage::dslFlucQForce;
794 }
795
796 info->setAtomStorageLayout(atomStorageLayout);
797 info->setRigidBodyStorageLayout(rigidBodyStorageLayout);
798 info->setCutoffGroupStorageLayout(cutoffGroupStorageLayout);
799
800 return;
801 }
802
803 void SimCreator::setGlobalIndex(SimInfo* info) {
804 SimInfo::MoleculeIterator mi;
805 Molecule::AtomIterator ai;
806 Molecule::RigidBodyIterator ri;
807 Molecule::CutoffGroupIterator ci;
808 Molecule::BondIterator boi;
809 Molecule::BendIterator bei;
810 Molecule::TorsionIterator ti;
811 Molecule::InversionIterator ii;
812 Molecule::IntegrableObjectIterator ioi;
813 Molecule* mol;
814 Atom* atom;
815 RigidBody* rb;
816 CutoffGroup* cg;
817 Bond* bond;
818 Bend* bend;
819 Torsion* torsion;
820 Inversion* inversion;
821 int beginAtomIndex;
822 int beginRigidBodyIndex;
823 int beginCutoffGroupIndex;
824 int beginBondIndex;
825 int beginBendIndex;
826 int beginTorsionIndex;
827 int beginInversionIndex;
828#ifdef IS_MPI
829 int nGlobalAtoms = info->getNGlobalAtoms();
830 int nGlobalRigidBodies = info->getNGlobalRigidBodies();
831#endif
832
833 beginAtomIndex = 0;
834 // The rigid body indices begin immediately after the atom indices:
835 beginRigidBodyIndex = info->getNGlobalAtoms();
836 beginCutoffGroupIndex = 0;
837 beginBondIndex = 0;
838 beginBendIndex = 0;
839 beginTorsionIndex = 0;
840 beginInversionIndex = 0;
841
842 for (int i = 0; i < info->getNGlobalMolecules(); i++) {
843#ifdef IS_MPI
844 if (info->getMolToProc(i) == worldRank) {
845#endif
846 // stuff to do if I own this molecule
847 mol = info->getMoleculeByGlobalIndex(i);
848
849 // The local index(index in DataStorge) of the atom is important:
850 for (atom = mol->beginAtom(ai); atom != NULL;
851 atom = mol->nextAtom(ai)) {
852 atom->setGlobalIndex(beginAtomIndex++);
853 }
854
855 for (rb = mol->beginRigidBody(ri); rb != NULL;
856 rb = mol->nextRigidBody(ri)) {
857 rb->setGlobalIndex(beginRigidBodyIndex++);
858 }
859
860 // The local index of other objects only depends on the order
861 // of traversal:
862 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
863 cg = mol->nextCutoffGroup(ci)) {
864 cg->setGlobalIndex(beginCutoffGroupIndex++);
865 }
866 for (bond = mol->beginBond(boi); bond != NULL;
867 bond = mol->nextBond(boi)) {
868 bond->setGlobalIndex(beginBondIndex++);
869 }
870 for (bend = mol->beginBend(bei); bend != NULL;
871 bend = mol->nextBend(bei)) {
872 bend->setGlobalIndex(beginBendIndex++);
873 }
874 for (torsion = mol->beginTorsion(ti); torsion != NULL;
875 torsion = mol->nextTorsion(ti)) {
876 torsion->setGlobalIndex(beginTorsionIndex++);
877 }
878 for (inversion = mol->beginInversion(ii); inversion != NULL;
879 inversion = mol->nextInversion(ii)) {
880 inversion->setGlobalIndex(beginInversionIndex++);
881 }
882
883#ifdef IS_MPI
884 } else {
885 // stuff to do if I don't own this molecule
886
887 int stampId = info->getMoleculeStampId(i);
888 MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
889
890 beginAtomIndex += stamp->getNAtoms();
891 beginRigidBodyIndex += stamp->getNRigidBodies();
892 beginCutoffGroupIndex +=
893 stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
894 beginBondIndex += stamp->getNBonds();
895 beginBendIndex += stamp->getNBends();
896 beginTorsionIndex += stamp->getNTorsions();
897 beginInversionIndex += stamp->getNInversions();
898 }
899#endif
900
901 } // end for(int i=0)
902
903 // fill globalGroupMembership
904 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
905 for (mol = info->beginMolecule(mi); mol != NULL;
906 mol = info->nextMolecule(mi)) {
907 for (cg = mol->beginCutoffGroup(ci); cg != NULL;
908 cg = mol->nextCutoffGroup(ci)) {
909 for (atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
910 globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
911 }
912 }
913 }
914
915#ifdef IS_MPI
916 // Since the globalGroupMembership has been zero filled and we've only
917 // poked values into the atoms we know, we can do an Allreduce
918 // to get the full globalGroupMembership array (We think).
919 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
920 // docs said we could.
921 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
922 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0],
923 nGlobalAtoms, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
924
925 info->setGlobalGroupMembership(tmpGroupMembership);
926#else
927 info->setGlobalGroupMembership(globalGroupMembership);
928#endif
929
930 // fill molMembership
931 std::vector<int> globalMolMembership(
932 info->getNGlobalAtoms() + info->getNGlobalRigidBodies(), 0);
933
934 for (mol = info->beginMolecule(mi); mol != NULL;
935 mol = info->nextMolecule(mi)) {
936 for (atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
937 globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
938 }
939 for (rb = mol->beginRigidBody(ri); rb != NULL;
940 rb = mol->nextRigidBody(ri)) {
941 globalMolMembership[rb->getGlobalIndex()] = mol->getGlobalIndex();
942 }
943 }
944
945#ifdef IS_MPI
946 std::vector<int> tmpMolMembership(
947 info->getNGlobalAtoms() + info->getNGlobalRigidBodies(), 0);
948 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
949 nGlobalAtoms + nGlobalRigidBodies, MPI_INT, MPI_SUM,
950 MPI_COMM_WORLD);
951
952 info->setGlobalMolMembership(tmpMolMembership);
953#else
954 info->setGlobalMolMembership(globalMolMembership);
955#endif
956
957 // nIOPerMol holds the number of integrable objects per molecule
958 // here the molecules are listed by their global indices.
959
960 std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
961 for (mol = info->beginMolecule(mi); mol != NULL;
962 mol = info->nextMolecule(mi)) {
963 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
964 }
965
966#ifdef IS_MPI
967 std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
968 MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
969 info->getNGlobalMolecules(), MPI_INT, MPI_SUM,
970 MPI_COMM_WORLD);
971#else
972 std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
973#endif
974
975 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
976
977 int startingIndex = 0;
978 for (int i = 0; i < info->getNGlobalMolecules(); i++) {
979 startingIOIndexForMol[i] = startingIndex;
980 startingIndex += numIntegrableObjectsPerMol[i];
981 }
982
983 std::vector<StuntDouble*> IOIndexToIntegrableObject(
985 for (mol = info->beginMolecule(mi); mol != NULL;
986 mol = info->nextMolecule(mi)) {
987 int myGlobalIndex = mol->getGlobalIndex();
988 int globalIO = startingIOIndexForMol[myGlobalIndex];
989 for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
990 sd = mol->nextIntegrableObject(ioi)) {
991 sd->setGlobalIntegrableObjectIndex(globalIO);
992 IOIndexToIntegrableObject[globalIO] = sd;
993 globalIO++;
994 }
995 }
996
997 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
998 }
999
1000 void SimCreator::loadCoordinates(SimInfo* info,
1001 const std::string& mdFileName) {
1002 DumpReader reader(info, mdFileName);
1003 int nframes = reader.getNFrames();
1004
1005 if (nframes > 0) {
1006 reader.readFrame(nframes - 1);
1007 } else {
1008 // invalid initial coordinate file
1009 snprintf(
1010 painCave.errMsg, MAX_SIM_ERROR_MSG_LENGTH,
1011 "Initial configuration file %s should at least contain one frame\n",
1012 mdFileName.c_str());
1013 painCave.isFatal = 1;
1014 simError();
1015 }
1016 // copy the current snapshot to previous snapshot
1017 info->getSnapshotManager()->advance();
1018 }
1019
1020} // namespace OpenMD
size_t getNIntegrableObjects()
Returns the total number of integrable objects in this molecule.
Definition Molecule.hpp:179
void setGlobalIndex(int index)
Sets the global index of this molecule.
Definition Molecule.hpp:126
int getGlobalIndex()
Returns the global index of this molecule.
Definition Molecule.hpp:107
SimInfo * createSim(const std::string &mdFileName, bool loadInitCoords=true)
Setup Simulation.
One of the heavy-weight classes of OpenMD, SimInfo maintains objects and variables relating to the cu...
Definition SimInfo.hpp:93
Molecule * getMoleculeByGlobalIndex(int index)
Finds a molecule with a specified global index.
Definition SimInfo.hpp:300
int getNGlobalIntegrableObjects()
Returns the total number of integrable objects (total number of rigid bodies plus the total number of...
Definition SimInfo.hpp:139
void setGlobalMolMembership(const std::vector< int > &gmm)
Sets GlobalMolMembership.
Definition SimInfo.hpp:370
void setAtomStorageLayout(int asl)
Sets the storage layouts (computed by SimCreator)
Definition SimInfo.hpp:256
int getMolToProc(int globalIndex)
Finds the processor where a molecule resides.
Definition SimInfo.hpp:667
ForceField * getForceField()
Returns the force field.
Definition SimInfo.hpp:266
Molecule * beginMolecule(MoleculeIterator &i)
Returns the first molecule in this SimInfo and intialize the iterator.
Definition SimInfo.cpp:240
MoleculeStamp * getMoleculeStamp(int id)
Returns the molecule stamp.
Definition SimInfo.hpp:290
int getNGlobalRigidBodies()
Returns the total number of integrable objects (total number of rigid bodies plus the total number of...
Definition SimInfo.hpp:146
int getAtomStorageLayout()
Returns the storage layouts (computed by SimCreator)
Definition SimInfo.hpp:251
int getNGlobalAtoms()
Returns the total number of atoms in the system.
Definition SimInfo.hpp:129
int getNGlobalMolecules()
Returns the total number of molecules in the system.
Definition SimInfo.hpp:126
void setSnapshotManager(SnapshotManager *sman)
Sets the snapshot manager.
Definition SimInfo.cpp:971
void setGlobalGroupMembership(const std::vector< int > &ggm)
Sets GlobalGroupMembership.
Definition SimInfo.hpp:362
bool addMolecule(Molecule *mol)
Adds a molecule.
Definition SimInfo.cpp:187
void addInteractionPairs(Molecule *mol)
add all special interaction pairs (including excluded interactions) in a molecule into the appropriat...
Definition SimInfo.cpp:376
Molecule * nextMolecule(MoleculeIterator &i)
Returns the next avaliable Molecule based on the iterator.
Definition SimInfo.cpp:245
SnapshotManager * getSnapshotManager()
Returns the snapshot manager.
Definition SimInfo.hpp:248
void setMolToProcMap(const std::vector< int > &molToProcMap)
Set MolToProcMap array.
Definition SimInfo.hpp:675
AtomTypeSet getSimulatedAtomTypes()
Returns the set of atom types present in this simulation.
Definition SimInfo.cpp:712
LocalIndexManager * getLocalIndexManager()
Returns the local index manager.
Definition SimInfo.hpp:282
"brains/SimSnapshotManager.hpp"
The string tokenizer class allows an application to break a string into tokens The set of delimiters ...
"Don't move, or you're dead! Stand up! Captain, we've got them!"
virtual std::string getMessage() const
Return error message without additional info (if present)
AST Super Factory shared by TreeParser and Parser.
std::string getMessage() const
Returns a clean error message (no line number/column information)
std::string getMessage() const
Returns a clean error message (no line number/column information)
std::string getMessage() const
Returns a clean error message (no line number/column information)
std::string getMessage() const
Returns a clean error message (no line number/column information)
virtual std::string getFilename() const
Return file where mishap occurred.
Baseclass for exceptions thrown by classes implementing the TokenStream interface.
Exception thrown from generated lexers when there's no default error handler specified.
This basic Periodic Table class was originally taken from the data.cpp file in OpenBabel.
std::string getPrefix(const std::string &str)