ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1983
Committed: Tue Apr 15 20:36:19 2014 UTC (11 years ago) by gezelter
File size: 34140 byte(s)
Log Message:
Fixes for ConstraintWriter in parallel, some cppcheck cleanup
starting preparation for 2.2 release

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date