ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1613
Committed: Thu Aug 18 20:18:19 2011 UTC (13 years, 8 months ago) by gezelter
File size: 27279 byte(s)
Log Message:
Fixed a parallel bug in computing exclude lists.
Added file versioning information in MD files.
Still tracking down cutoff group bugs.

File Contents

# User Rev Content
1 gezelter 403 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
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     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39     * [4] Vardeman & Gezelter, in progress (2009).
40 gezelter 403 */
41    
42     /**
43     * @file SimCreator.cpp
44     * @author tlin
45     * @date 11/03/2004
46     * @time 13:51am
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     #include "UseTheForce/ForceFieldFactory.hpp"
59     #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 403 #ifdef IS_MPI
79     #include "math/ParallelRandNumGen.hpp"
80     #endif
81    
82 gezelter 1390 namespace OpenMD {
83 chrisfen 417
84 gezelter 1613 Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int mdFileVersion, int startOfMetaDataBlock ){
85 tim 1024 Globals* simParams = NULL;
86     try {
87 tim 770
88 tim 1024 // Create a preprocessor that preprocesses md file into an ostringstream
89     std::stringstream ppStream;
90 tim 770 #ifdef IS_MPI
91 tim 1024 int streamSize;
92     const int masterNode = 0;
93     int commStatus;
94     if (worldRank == masterNode) {
95 gezelter 1613 commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
96     #endif
97 tim 1024 SimplePreprocessor preprocessor;
98     preprocessor.preprocess(rawMetaDataStream, filename, startOfMetaDataBlock, ppStream);
99 tim 770
100     #ifdef IS_MPI
101 tim 1024 //brocasting the stream size
102     streamSize = ppStream.str().size() +1;
103     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
104 tim 770
105 tim 1024 commStatus = MPI_Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())), streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
106 tim 770
107    
108 tim 1024 } else {
109 gezelter 1613
110     commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
111    
112 tim 1024 //get stream size
113     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
114 gezelter 1313
115 tim 1024 char* buf = new char[streamSize];
116     assert(buf);
117 tim 770
118 tim 1024 //receive file content
119     commStatus = MPI_Bcast(buf, streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
120 tim 770
121 tim 1024 ppStream.str(buf);
122 gezelter 1313 delete [] buf;
123 tim 770
124 tim 1024 }
125 tim 770 #endif
126 tim 1024 // Create a scanner that reads from the input stream
127     MDLexer lexer(ppStream);
128     lexer.setFilename(filename);
129     lexer.initDeferredLineCount();
130 chrisfen 417
131 tim 1024 // Create a parser that reads from the scanner
132     MDParser parser(lexer);
133     parser.setFilename(filename);
134 tim 770
135 tim 1024 // Create an observer that synchorizes file name change
136     FilenameObserver observer;
137     observer.setLexer(&lexer);
138     observer.setParser(&parser);
139     lexer.setObserver(&observer);
140 chrisfen 417
141 tim 1024 antlr::ASTFactory factory;
142     parser.initializeASTFactory(factory);
143     parser.setASTFactory(&factory);
144     parser.mdfile();
145 tim 770
146 tim 1024 // Create a tree parser that reads information into Globals
147     MDTreeParser treeParser;
148     treeParser.initializeASTFactory(factory);
149     treeParser.setASTFactory(&factory);
150     simParams = treeParser.walkTree(parser.getAST());
151     }
152 tim 845
153 tim 816
154 tim 1024 catch(antlr::MismatchedCharException& e) {
155     sprintf(painCave.errMsg,
156     "parser exception: %s %s:%d:%d\n",
157     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
158     painCave.isFatal = 1;
159     simError();
160     }
161     catch(antlr::MismatchedTokenException &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::NoViableAltForCharException &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::NoViableAltException &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 tim 845
183 tim 1024 catch(antlr::TokenStreamRecognitionException& e) {
184     sprintf(painCave.errMsg,
185     "parser exception: %s %s:%d:%d\n",
186     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
187     painCave.isFatal = 1;
188     simError();
189     }
190 tim 845
191 tim 1024 catch(antlr::TokenStreamIOException& e) {
192     sprintf(painCave.errMsg,
193     "parser exception: %s\n",
194     e.getMessage().c_str());
195     painCave.isFatal = 1;
196     simError();
197     }
198 tim 845
199 tim 1024 catch(antlr::TokenStreamException& e) {
200     sprintf(painCave.errMsg,
201     "parser exception: %s\n",
202     e.getMessage().c_str());
203     painCave.isFatal = 1;
204     simError();
205     }
206     catch (antlr::RecognitionException& e) {
207     sprintf(painCave.errMsg,
208     "parser exception: %s %s:%d:%d\n",
209     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
210     painCave.isFatal = 1;
211     simError();
212     }
213     catch (antlr::CharStreamException& e) {
214     sprintf(painCave.errMsg,
215     "parser exception: %s\n",
216     e.getMessage().c_str());
217     painCave.isFatal = 1;
218     simError();
219     }
220 gezelter 1390 catch (OpenMDException& e) {
221 tim 1024 sprintf(painCave.errMsg,
222     "%s\n",
223     e.getMessage().c_str());
224     painCave.isFatal = 1;
225     simError();
226     }
227     catch (std::exception& e) {
228     sprintf(painCave.errMsg,
229     "parser exception: %s\n",
230     e.what());
231     painCave.isFatal = 1;
232     simError();
233     }
234 tim 770
235 gezelter 1613 simParams->setMDfileVersion(mdFileVersion);
236 tim 1024 return simParams;
237 gezelter 403 }
238 chrisfen 417
239 chrisfen 514 SimInfo* SimCreator::createSim(const std::string & mdFileName,
240     bool loadInitCoords) {
241 gezelter 1469
242 tim 1024 const int bufferSize = 65535;
243     char buffer[bufferSize];
244     int lineNo = 0;
245     std::string mdRawData;
246     int metaDataBlockStart = -1;
247     int metaDataBlockEnd = -1;
248     int i;
249     int mdOffset;
250 gezelter 1613 int mdFileVersion;
251 tim 1024
252     #ifdef IS_MPI
253     const int masterNode = 0;
254     if (worldRank == masterNode) {
255     #endif
256    
257     std::ifstream mdFile_(mdFileName.c_str());
258    
259     if (mdFile_.fail()) {
260     sprintf(painCave.errMsg,
261     "SimCreator: Cannot open file: %s\n",
262     mdFileName.c_str());
263     painCave.isFatal = 1;
264     simError();
265     }
266    
267     mdFile_.getline(buffer, bufferSize);
268     ++lineNo;
269     std::string line = trimLeftCopy(buffer);
270 gezelter 1390 i = CaseInsensitiveFind(line, "<OpenMD");
271 gezelter 1287 if (static_cast<size_t>(i) == string::npos) {
272 gezelter 1390 // try the older file strings to see if that works:
273     i = CaseInsensitiveFind(line, "<OOPSE");
274     }
275    
276     if (static_cast<size_t>(i) == string::npos) {
277     // still no luck!
278 tim 1024 sprintf(painCave.errMsg,
279 gezelter 1390 "SimCreator: File: %s is not a valid OpenMD file!\n",
280 tim 1024 mdFileName.c_str());
281     painCave.isFatal = 1;
282     simError();
283     }
284 gezelter 1613
285     // found the correct opening string, now try to get the file
286     // format version number.
287 tim 1024
288 gezelter 1613 StringTokenizer tokenizer(line, "=<> \t\n\r");
289     std::string fileType = tokenizer.nextToken();
290     toUpper(fileType);
291    
292     mdFileVersion = 0;
293    
294     if (fileType == "OPENMD") {
295     while (tokenizer.hasMoreTokens()) {
296     std::string token(tokenizer.nextToken());
297     toUpper(token);
298     if (token == "VERSION") {
299     mdFileVersion = tokenizer.nextTokenAsInt();
300     break;
301     }
302     }
303     }
304    
305 tim 1024 //scan through the input stream and find MetaData tag
306     while(mdFile_.getline(buffer, bufferSize)) {
307     ++lineNo;
308    
309     std::string line = trimLeftCopy(buffer);
310     if (metaDataBlockStart == -1) {
311     i = CaseInsensitiveFind(line, "<MetaData>");
312     if (i != string::npos) {
313     metaDataBlockStart = lineNo;
314     mdOffset = mdFile_.tellg();
315     }
316     } else {
317     i = CaseInsensitiveFind(line, "</MetaData>");
318     if (i != string::npos) {
319     metaDataBlockEnd = lineNo;
320     }
321     }
322     }
323    
324     if (metaDataBlockStart == -1) {
325     sprintf(painCave.errMsg,
326     "SimCreator: File: %s did not contain a <MetaData> tag!\n",
327     mdFileName.c_str());
328     painCave.isFatal = 1;
329     simError();
330     }
331     if (metaDataBlockEnd == -1) {
332     sprintf(painCave.errMsg,
333     "SimCreator: File: %s did not contain a closed MetaData block!\n",
334     mdFileName.c_str());
335     painCave.isFatal = 1;
336     simError();
337     }
338    
339     mdFile_.clear();
340     mdFile_.seekg(0);
341     mdFile_.seekg(mdOffset);
342    
343     mdRawData.clear();
344    
345     for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
346     mdFile_.getline(buffer, bufferSize);
347     mdRawData += buffer;
348     mdRawData += "\n";
349     }
350    
351     mdFile_.close();
352    
353     #ifdef IS_MPI
354     }
355     #endif
356    
357     std::stringstream rawMetaDataStream(mdRawData);
358    
359 gezelter 403 //parse meta-data file
360 gezelter 1613 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
361     metaDataBlockStart + 1);
362 chrisfen 417
363 gezelter 403 //create the force field
364 gezelter 1277 ForceField * ff = ForceFieldFactory::getInstance()->createForceField(simParams->getForceField());
365    
366 gezelter 403 if (ff == NULL) {
367 chrisfen 514 sprintf(painCave.errMsg,
368     "ForceField Factory can not create %s force field\n",
369 tim 665 simParams->getForceField().c_str());
370 gezelter 403 painCave.isFatal = 1;
371     simError();
372     }
373 chrisfen 417
374 gezelter 403 if (simParams->haveForceFieldFileName()) {
375     ff->setForceFieldFileName(simParams->getForceFieldFileName());
376     }
377    
378     std::string forcefieldFileName;
379     forcefieldFileName = ff->getForceFieldFileName();
380 chrisfen 417
381 gezelter 403 if (simParams->haveForceFieldVariant()) {
382     //If the force field has variant, the variant force field name will be
383     //Base.variant.frc. For exampel EAM.u6.frc
384 chrisfen 417
385 gezelter 403 std::string variant = simParams->getForceFieldVariant();
386 chrisfen 417
387 gezelter 403 std::string::size_type pos = forcefieldFileName.rfind(".frc");
388     variant = "." + variant;
389     if (pos != std::string::npos) {
390     forcefieldFileName.insert(pos, variant);
391     } else {
392     //If the default force field file name does not containt .frc suffix, just append the .variant
393     forcefieldFileName.append(variant);
394     }
395     }
396    
397     ff->parse(forcefieldFileName);
398     //create SimInfo
399 tim 770 SimInfo * info = new SimInfo(ff, simParams);
400 tim 1024
401     info->setRawMetaData(mdRawData);
402 tim 490
403 gezelter 945 //gather parameters (SimCreator only retrieves part of the
404     //parameters)
405 gezelter 403 gatherParameters(info, mdFileName);
406 chrisfen 417
407 gezelter 403 //divide the molecules and determine the global index of molecules
408     #ifdef IS_MPI
409     divideMolecules(info);
410     #endif
411 chrisfen 417
412 gezelter 403 //create the molecules
413     createMolecules(info);
414 chrisfen 417
415 gezelter 945 //allocate memory for DataStorage(circular reference, need to
416     //break it)
417 gezelter 403 info->setSnapshotManager(new SimSnapshotManager(info));
418    
419 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
420     //(only need to be set once, the global index will never change
421     //again). Local indices of atoms and rigidbodies are already set
422     //by MoleculeCreator class which actually delegates the
423     //responsibility to LocalIndexManager.
424 gezelter 403 setGlobalIndex(info);
425 chrisfen 417
426 gezelter 1287 //Although addInteractionPairs is called inside SimInfo's addMolecule
427 gezelter 945 //method, at that point atoms don't have the global index yet
428     //(their global index are all initialized to -1). Therefore we
429 gezelter 1287 //have to call addInteractionPairs explicitly here. A way to work
430 gezelter 945 //around is that we can determine the beginning global indices of
431     //atoms before they get created.
432 gezelter 403 SimInfo::MoleculeIterator mi;
433     Molecule* mol;
434     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
435 gezelter 1287 info->addInteractionPairs(mol);
436 gezelter 403 }
437    
438     if (loadInitCoords)
439 tim 1024 loadCoordinates(info, mdFileName);
440 gezelter 403 return info;
441     }
442 chrisfen 417
443 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
444 chrisfen 417
445 tim 749 //figure out the output file names
446 gezelter 403 std::string prefix;
447 chrisfen 417
448 gezelter 403 #ifdef IS_MPI
449 chrisfen 417
450 gezelter 403 if (worldRank == 0) {
451     #endif // is_mpi
452     Globals * simParams = info->getSimParams();
453     if (simParams->haveFinalConfig()) {
454     prefix = getPrefix(simParams->getFinalConfig());
455     } else {
456     prefix = getPrefix(mdfile);
457     }
458 chrisfen 417
459 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
460     info->setDumpFileName(prefix + ".dump");
461     info->setStatFileName(prefix + ".stat");
462 chrisfen 417 info->setRestFileName(prefix + ".zang");
463    
464 gezelter 403 #ifdef IS_MPI
465 chrisfen 417
466 gezelter 403 }
467 chrisfen 417
468 gezelter 403 #endif
469 chrisfen 417
470 gezelter 403 }
471 chrisfen 417
472 gezelter 403 #ifdef IS_MPI
473     void SimCreator::divideMolecules(SimInfo *info) {
474 tim 963 RealType numerator;
475     RealType denominator;
476     RealType precast;
477     RealType x;
478     RealType y;
479     RealType a;
480 gezelter 403 int old_atoms;
481     int add_atoms;
482     int new_atoms;
483     int nTarget;
484     int done;
485     int i;
486     int j;
487     int loops;
488     int which_proc;
489     int nProcessors;
490     std::vector<int> atomsPerProc;
491     int nGlobalMols = info->getNGlobalMolecules();
492     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
493    
494     MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
495 chrisfen 417
496 gezelter 403 if (nProcessors > nGlobalMols) {
497     sprintf(painCave.errMsg,
498     "nProcessors (%d) > nMol (%d)\n"
499     "\tThe number of processors is larger than\n"
500     "\tthe number of molecules. This will not result in a \n"
501     "\tusable division of atoms for force decomposition.\n"
502     "\tEither try a smaller number of processors, or run the\n"
503 gezelter 1390 "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols);
504 chrisfen 417
505 gezelter 403 painCave.isFatal = 1;
506     simError();
507     }
508 chrisfen 417
509 gezelter 403 int seedValue;
510     Globals * simParams = info->getSimParams();
511     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
512     if (simParams->haveSeed()) {
513     seedValue = simParams->getSeed();
514     myRandom = new SeqRandNumGen(seedValue);
515     }else {
516     myRandom = new SeqRandNumGen();
517     }
518 chrisfen 417
519    
520 gezelter 403 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
521 chrisfen 417
522 gezelter 403 //initialize atomsPerProc
523     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
524 chrisfen 417
525 gezelter 403 if (worldRank == 0) {
526     numerator = info->getNGlobalAtoms();
527     denominator = nProcessors;
528     precast = numerator / denominator;
529     nTarget = (int)(precast + 0.5);
530 chrisfen 417
531 gezelter 403 for(i = 0; i < nGlobalMols; i++) {
532     done = 0;
533     loops = 0;
534 chrisfen 417
535 gezelter 403 while (!done) {
536     loops++;
537 chrisfen 417
538 gezelter 403 // Pick a processor at random
539 chrisfen 417
540 gezelter 403 which_proc = (int) (myRandom->rand() * nProcessors);
541 chrisfen 417
542 gezelter 403 //get the molecule stamp first
543     int stampId = info->getMoleculeStampId(i);
544     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
545 chrisfen 417
546 gezelter 403 // How many atoms does this processor have so far?
547     old_atoms = atomsPerProc[which_proc];
548     add_atoms = moleculeStamp->getNAtoms();
549     new_atoms = old_atoms + add_atoms;
550 chrisfen 417
551 gezelter 403 // If we've been through this loop too many times, we need
552     // to just give up and assign the molecule to this processor
553     // and be done with it.
554 chrisfen 417
555 gezelter 403 if (loops > 100) {
556     sprintf(painCave.errMsg,
557     "I've tried 100 times to assign molecule %d to a "
558     " processor, but can't find a good spot.\n"
559     "I'm assigning it at random to processor %d.\n",
560     i, which_proc);
561 chrisfen 417
562 gezelter 403 painCave.isFatal = 0;
563     simError();
564 chrisfen 417
565 gezelter 403 molToProcMap[i] = which_proc;
566     atomsPerProc[which_proc] += add_atoms;
567 chrisfen 417
568 gezelter 403 done = 1;
569     continue;
570     }
571 chrisfen 417
572 gezelter 403 // If we can add this molecule to this processor without sending
573     // it above nTarget, then go ahead and do it:
574 chrisfen 417
575 gezelter 403 if (new_atoms <= nTarget) {
576     molToProcMap[i] = which_proc;
577     atomsPerProc[which_proc] += add_atoms;
578 chrisfen 417
579 gezelter 403 done = 1;
580     continue;
581     }
582 chrisfen 417
583 gezelter 403 // The only situation left is when new_atoms > nTarget. We
584     // want to accept this with some probability that dies off the
585     // farther we are from nTarget
586 chrisfen 417
587 gezelter 403 // roughly: x = new_atoms - nTarget
588     // Pacc(x) = exp(- a * x)
589     // where a = penalty / (average atoms per molecule)
590 chrisfen 417
591 tim 963 x = (RealType)(new_atoms - nTarget);
592 gezelter 403 y = myRandom->rand();
593 chrisfen 417
594 gezelter 403 if (y < exp(- a * x)) {
595     molToProcMap[i] = which_proc;
596     atomsPerProc[which_proc] += add_atoms;
597 chrisfen 417
598 gezelter 403 done = 1;
599     continue;
600     } else {
601     continue;
602     }
603     }
604     }
605 chrisfen 417
606 gezelter 403 delete myRandom;
607 chrisfen 417
608 gezelter 403 // Spray out this nonsense to all other processors:
609 chrisfen 417
610 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
611     } else {
612 chrisfen 417
613 gezelter 403 // Listen to your marching orders from processor 0:
614 chrisfen 417
615 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
616     }
617 chrisfen 417
618 gezelter 403 info->setMolToProcMap(molToProcMap);
619     sprintf(checkPointMsg,
620     "Successfully divided the molecules among the processors.\n");
621 gezelter 1241 errorCheckPoint();
622 gezelter 403 }
623 chrisfen 417
624 gezelter 403 #endif
625 chrisfen 417
626 gezelter 403 void SimCreator::createMolecules(SimInfo *info) {
627     MoleculeCreator molCreator;
628     int stampId;
629 chrisfen 417
630 gezelter 403 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
631 chrisfen 417
632 gezelter 403 #ifdef IS_MPI
633 chrisfen 417
634 gezelter 403 if (info->getMolToProc(i) == worldRank) {
635     #endif
636 chrisfen 417
637 gezelter 403 stampId = info->getMoleculeStampId(i);
638 gezelter 1600 Molecule * mol = molCreator.createMolecule(info->getForceField(),
639     info->getMoleculeStamp(stampId),
640     stampId, i,
641     info->getLocalIndexManager());
642 chrisfen 417
643 gezelter 403 info->addMolecule(mol);
644 chrisfen 417
645 gezelter 403 #ifdef IS_MPI
646 chrisfen 417
647 gezelter 403 }
648 chrisfen 417
649 gezelter 403 #endif
650 chrisfen 417
651 gezelter 403 } //end for(int i=0)
652     }
653 chrisfen 417
654 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
655     SimInfo::MoleculeIterator mi;
656     Molecule::AtomIterator ai;
657     Molecule::RigidBodyIterator ri;
658     Molecule::CutoffGroupIterator ci;
659 tim 1024 Molecule::IntegrableObjectIterator ioi;
660 gezelter 403 Molecule * mol;
661     Atom * atom;
662     RigidBody * rb;
663     CutoffGroup * cg;
664     int beginAtomIndex;
665     int beginRigidBodyIndex;
666     int beginCutoffGroupIndex;
667     int nGlobalAtoms = info->getNGlobalAtoms();
668 chrisfen 417
669 gezelter 403 beginAtomIndex = 0;
670     beginRigidBodyIndex = 0;
671     beginCutoffGroupIndex = 0;
672 gezelter 1600
673     for(int i = 0; i < info->getNGlobalMolecules(); i++) {
674 chrisfen 417
675 gezelter 1600 #ifdef IS_MPI
676     if (info->getMolToProc(i) == worldRank) {
677     #endif
678     // stuff to do if I own this molecule
679     mol = info->getMoleculeByGlobalIndex(i);
680    
681     //local index(index in DataStorge) of atom is important
682     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
683     atom->setGlobalIndex(beginAtomIndex++);
684     }
685    
686     for(rb = mol->beginRigidBody(ri); rb != NULL;
687     rb = mol->nextRigidBody(ri)) {
688     rb->setGlobalIndex(beginRigidBodyIndex++);
689     }
690    
691     //local index of cutoff group is trivial, it only depends on
692     //the order of travesing
693     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
694     cg = mol->nextCutoffGroup(ci)) {
695     cg->setGlobalIndex(beginCutoffGroupIndex++);
696     }
697    
698     #ifdef IS_MPI
699     } else {
700    
701     // stuff to do if I don't own this molecule
702    
703     int stampId = info->getMoleculeStampId(i);
704     MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
705    
706     beginAtomIndex += stamp->getNAtoms();
707     beginRigidBodyIndex += stamp->getNRigidBodies();
708     beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
709 gezelter 403 }
710 gezelter 1600 #endif
711    
712     } //end for(int i=0)
713    
714 gezelter 403 //fill globalGroupMembership
715 gezelter 1600 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
716 gezelter 403 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
717     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
718 chrisfen 417
719 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
720     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
721     }
722 chrisfen 417
723 gezelter 403 }
724     }
725 gezelter 1577
726 gezelter 403 #ifdef IS_MPI
727     // Since the globalGroupMembership has been zero filled and we've only
728     // poked values into the atoms we know, we can do an Allreduce
729     // to get the full globalGroupMembership array (We think).
730     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
731     // docs said we could.
732 gezelter 1313 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
733 gezelter 403 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
734 gezelter 1600 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
735 gezelter 403 info->setGlobalGroupMembership(tmpGroupMembership);
736     #else
737     info->setGlobalGroupMembership(globalGroupMembership);
738     #endif
739 chrisfen 417
740 gezelter 403 //fill molMembership
741     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
742    
743     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
744     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
745     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
746     }
747     }
748 chrisfen 417
749 gezelter 403 #ifdef IS_MPI
750 gezelter 1313 std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
751 chrisfen 417
752 gezelter 403 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
753     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
754    
755     info->setGlobalMolMembership(tmpMolMembership);
756     #else
757     info->setGlobalMolMembership(globalMolMembership);
758     #endif
759 tim 1024
760     // nIOPerMol holds the number of integrable objects per molecule
761     // here the molecules are listed by their global indices.
762    
763     std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
764     for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
765     nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
766     }
767 chrisfen 417
768 tim 1024 #ifdef IS_MPI
769     std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
770     MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
771     info->getNGlobalMolecules(), MPI_INT, MPI_SUM, MPI_COMM_WORLD);
772     #else
773     std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
774     #endif
775    
776 gezelter 1313 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
777    
778     int startingIndex = 0;
779     for (int i = 0; i < info->getNGlobalMolecules(); i++) {
780     startingIOIndexForMol[i] = startingIndex;
781     startingIndex += numIntegrableObjectsPerMol[i];
782     }
783    
784     std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
785     for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
786 tim 1024 int myGlobalIndex = mol->getGlobalIndex();
787     int globalIO = startingIOIndexForMol[myGlobalIndex];
788     for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
789     integrableObject = mol->nextIntegrableObject(ioi)) {
790 gezelter 1313 integrableObject->setGlobalIntegrableObjectIndex(globalIO);
791     IOIndexToIntegrableObject[globalIO] = integrableObject;
792     globalIO++;
793 tim 1024 }
794     }
795 gezelter 1597
796 gezelter 1313 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
797    
798 gezelter 403 }
799 chrisfen 417
800 tim 1024 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
801 gezelter 403 Globals* simParams;
802 gezelter 1540
803 gezelter 403 simParams = info->getSimParams();
804    
805 tim 1024 DumpReader reader(info, mdFileName);
806 gezelter 403 int nframes = reader.getNFrames();
807 gezelter 1540
808 gezelter 403 if (nframes > 0) {
809     reader.readFrame(nframes - 1);
810     } else {
811     //invalid initial coordinate file
812 chrisfen 417 sprintf(painCave.errMsg,
813     "Initial configuration file %s should at least contain one frame\n",
814 tim 1024 mdFileName.c_str());
815 gezelter 403 painCave.isFatal = 1;
816     simError();
817     }
818     //copy the current snapshot to previous snapshot
819     info->getSnapshotManager()->advance();
820     }
821 chrisfen 417
822 gezelter 1390 } //end namespace OpenMD
823 gezelter 403
824    

Properties

Name Value
svn:keywords Author Id Revision Date