ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1971
Committed: Fri Feb 28 13:25:13 2014 UTC (11 years, 2 months ago) by gezelter
File size: 35387 byte(s)
Log Message:
Fixed a few MPI issues that were triggered in mpich.

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

Properties

Name Value
svn:keywords Author Id Revision Date