ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1880
Committed: Mon Jun 17 18:28:30 2013 UTC (11 years, 10 months ago) by gezelter
File size: 32711 byte(s)
Log Message:
Preparing for official 2.1 release (clean-up)

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

Properties

Name Value
svn:keywords Author Id Revision Date