ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1577
Committed: Wed Jun 8 20:26:56 2011 UTC (13 years, 10 months ago) by gezelter
File size: 27186 byte(s)
Log Message:
bug fixes

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

Properties

Name Value
svn:keywords Author Id Revision Date