ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1938
Committed: Thu Oct 31 15:32:17 2013 UTC (11 years, 6 months ago) by gezelter
File size: 32713 byte(s)
Log Message:
Some MPI include re-ordering to work with the Intel MPI implementation.

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

Properties

Name Value
svn:keywords Author Id Revision Date