ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1535
Committed: Fri Dec 31 18:31:56 2010 UTC (14 years, 4 months ago) by gezelter
File size: 27210 byte(s)
Log Message:
Well, it compiles and builds, but still has a bus error at runtime.

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    
390 gezelter 945 //allocate memory for DataStorage(circular reference, need to
391     //break it)
392 gezelter 403 info->setSnapshotManager(new SimSnapshotManager(info));
393    
394 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
395     //(only need to be set once, the global index will never change
396     //again). Local indices of atoms and rigidbodies are already set
397     //by MoleculeCreator class which actually delegates the
398     //responsibility to LocalIndexManager.
399 gezelter 403 setGlobalIndex(info);
400 chrisfen 417
401 gezelter 1287 //Although addInteractionPairs is called inside SimInfo's addMolecule
402 gezelter 945 //method, at that point atoms don't have the global index yet
403     //(their global index are all initialized to -1). Therefore we
404 gezelter 1287 //have to call addInteractionPairs explicitly here. A way to work
405 gezelter 945 //around is that we can determine the beginning global indices of
406     //atoms before they get created.
407 gezelter 403 SimInfo::MoleculeIterator mi;
408     Molecule* mol;
409     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
410 gezelter 1287 info->addInteractionPairs(mol);
411 gezelter 403 }
412    
413     if (loadInitCoords)
414 tim 1024 loadCoordinates(info, mdFileName);
415 gezelter 403
416     return info;
417     }
418 chrisfen 417
419 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
420 chrisfen 417
421 tim 749 //figure out the output file names
422 gezelter 403 std::string prefix;
423 chrisfen 417
424 gezelter 403 #ifdef IS_MPI
425 chrisfen 417
426 gezelter 403 if (worldRank == 0) {
427     #endif // is_mpi
428     Globals * simParams = info->getSimParams();
429     if (simParams->haveFinalConfig()) {
430     prefix = getPrefix(simParams->getFinalConfig());
431     } else {
432     prefix = getPrefix(mdfile);
433     }
434 chrisfen 417
435 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
436     info->setDumpFileName(prefix + ".dump");
437     info->setStatFileName(prefix + ".stat");
438 chrisfen 417 info->setRestFileName(prefix + ".zang");
439    
440 gezelter 403 #ifdef IS_MPI
441 chrisfen 417
442 gezelter 403 }
443 chrisfen 417
444 gezelter 403 #endif
445 chrisfen 417
446 gezelter 403 }
447 chrisfen 417
448 gezelter 403 #ifdef IS_MPI
449     void SimCreator::divideMolecules(SimInfo *info) {
450 tim 963 RealType numerator;
451     RealType denominator;
452     RealType precast;
453     RealType x;
454     RealType y;
455     RealType a;
456 gezelter 403 int old_atoms;
457     int add_atoms;
458     int new_atoms;
459     int nTarget;
460     int done;
461     int i;
462     int j;
463     int loops;
464     int which_proc;
465     int nProcessors;
466     std::vector<int> atomsPerProc;
467     int nGlobalMols = info->getNGlobalMolecules();
468     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
469    
470     MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
471 chrisfen 417
472 gezelter 403 if (nProcessors > nGlobalMols) {
473     sprintf(painCave.errMsg,
474     "nProcessors (%d) > nMol (%d)\n"
475     "\tThe number of processors is larger than\n"
476     "\tthe number of molecules. This will not result in a \n"
477     "\tusable division of atoms for force decomposition.\n"
478     "\tEither try a smaller number of processors, or run the\n"
479 gezelter 1390 "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols);
480 chrisfen 417
481 gezelter 403 painCave.isFatal = 1;
482     simError();
483     }
484 chrisfen 417
485 gezelter 403 int seedValue;
486     Globals * simParams = info->getSimParams();
487     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
488     if (simParams->haveSeed()) {
489     seedValue = simParams->getSeed();
490     myRandom = new SeqRandNumGen(seedValue);
491     }else {
492     myRandom = new SeqRandNumGen();
493     }
494 chrisfen 417
495    
496 gezelter 403 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
497 chrisfen 417
498 gezelter 403 //initialize atomsPerProc
499     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
500 chrisfen 417
501 gezelter 403 if (worldRank == 0) {
502     numerator = info->getNGlobalAtoms();
503     denominator = nProcessors;
504     precast = numerator / denominator;
505     nTarget = (int)(precast + 0.5);
506 chrisfen 417
507 gezelter 403 for(i = 0; i < nGlobalMols; i++) {
508     done = 0;
509     loops = 0;
510 chrisfen 417
511 gezelter 403 while (!done) {
512     loops++;
513 chrisfen 417
514 gezelter 403 // Pick a processor at random
515 chrisfen 417
516 gezelter 403 which_proc = (int) (myRandom->rand() * nProcessors);
517 chrisfen 417
518 gezelter 403 //get the molecule stamp first
519     int stampId = info->getMoleculeStampId(i);
520     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
521 chrisfen 417
522 gezelter 403 // How many atoms does this processor have so far?
523     old_atoms = atomsPerProc[which_proc];
524     add_atoms = moleculeStamp->getNAtoms();
525     new_atoms = old_atoms + add_atoms;
526 chrisfen 417
527 gezelter 403 // If we've been through this loop too many times, we need
528     // to just give up and assign the molecule to this processor
529     // and be done with it.
530 chrisfen 417
531 gezelter 403 if (loops > 100) {
532     sprintf(painCave.errMsg,
533     "I've tried 100 times to assign molecule %d to a "
534     " processor, but can't find a good spot.\n"
535     "I'm assigning it at random to processor %d.\n",
536     i, which_proc);
537 chrisfen 417
538 gezelter 403 painCave.isFatal = 0;
539     simError();
540 chrisfen 417
541 gezelter 403 molToProcMap[i] = which_proc;
542     atomsPerProc[which_proc] += add_atoms;
543 chrisfen 417
544 gezelter 403 done = 1;
545     continue;
546     }
547 chrisfen 417
548 gezelter 403 // If we can add this molecule to this processor without sending
549     // it above nTarget, then go ahead and do it:
550 chrisfen 417
551 gezelter 403 if (new_atoms <= nTarget) {
552     molToProcMap[i] = which_proc;
553     atomsPerProc[which_proc] += add_atoms;
554 chrisfen 417
555 gezelter 403 done = 1;
556     continue;
557     }
558 chrisfen 417
559 gezelter 403 // The only situation left is when new_atoms > nTarget. We
560     // want to accept this with some probability that dies off the
561     // farther we are from nTarget
562 chrisfen 417
563 gezelter 403 // roughly: x = new_atoms - nTarget
564     // Pacc(x) = exp(- a * x)
565     // where a = penalty / (average atoms per molecule)
566 chrisfen 417
567 tim 963 x = (RealType)(new_atoms - nTarget);
568 gezelter 403 y = myRandom->rand();
569 chrisfen 417
570 gezelter 403 if (y < exp(- a * x)) {
571     molToProcMap[i] = which_proc;
572     atomsPerProc[which_proc] += add_atoms;
573 chrisfen 417
574 gezelter 403 done = 1;
575     continue;
576     } else {
577     continue;
578     }
579     }
580     }
581 chrisfen 417
582 gezelter 403 delete myRandom;
583 chrisfen 417
584 gezelter 403 // Spray out this nonsense to all other processors:
585 chrisfen 417
586 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
587     } else {
588 chrisfen 417
589 gezelter 403 // Listen to your marching orders from processor 0:
590 chrisfen 417
591 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
592     }
593 chrisfen 417
594 gezelter 403 info->setMolToProcMap(molToProcMap);
595     sprintf(checkPointMsg,
596     "Successfully divided the molecules among the processors.\n");
597 gezelter 1241 errorCheckPoint();
598 gezelter 403 }
599 chrisfen 417
600 gezelter 403 #endif
601 chrisfen 417
602 gezelter 403 void SimCreator::createMolecules(SimInfo *info) {
603     MoleculeCreator molCreator;
604     int stampId;
605 chrisfen 417
606 gezelter 403 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
607 chrisfen 417
608 gezelter 403 #ifdef IS_MPI
609 chrisfen 417
610 gezelter 403 if (info->getMolToProc(i) == worldRank) {
611     #endif
612 chrisfen 417
613 gezelter 403 stampId = info->getMoleculeStampId(i);
614     Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
615     stampId, i, info->getLocalIndexManager());
616 chrisfen 417
617 gezelter 403 info->addMolecule(mol);
618 chrisfen 417
619 gezelter 403 #ifdef IS_MPI
620 chrisfen 417
621 gezelter 403 }
622 chrisfen 417
623 gezelter 403 #endif
624 chrisfen 417
625 gezelter 403 } //end for(int i=0)
626     }
627 chrisfen 417
628 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
629     SimInfo::MoleculeIterator mi;
630     Molecule::AtomIterator ai;
631     Molecule::RigidBodyIterator ri;
632     Molecule::CutoffGroupIterator ci;
633 tim 1024 Molecule::IntegrableObjectIterator ioi;
634 gezelter 403 Molecule * mol;
635     Atom * atom;
636     RigidBody * rb;
637     CutoffGroup * cg;
638     int beginAtomIndex;
639     int beginRigidBodyIndex;
640     int beginCutoffGroupIndex;
641     int nGlobalAtoms = info->getNGlobalAtoms();
642 tim 1024
643     /**@todo fixme */
644 gezelter 403 #ifndef IS_MPI
645 chrisfen 417
646 gezelter 403 beginAtomIndex = 0;
647     beginRigidBodyIndex = 0;
648     beginCutoffGroupIndex = 0;
649 chrisfen 417
650 gezelter 403 #else
651 chrisfen 417
652 gezelter 403 int nproc;
653     int myNode;
654 chrisfen 417
655 gezelter 403 myNode = worldRank;
656     MPI_Comm_size(MPI_COMM_WORLD, &nproc);
657 chrisfen 417
658 gezelter 403 std::vector < int > tmpAtomsInProc(nproc, 0);
659     std::vector < int > tmpRigidBodiesInProc(nproc, 0);
660     std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
661     std::vector < int > NumAtomsInProc(nproc, 0);
662     std::vector < int > NumRigidBodiesInProc(nproc, 0);
663     std::vector < int > NumCutoffGroupsInProc(nproc, 0);
664 chrisfen 417
665 gezelter 403 tmpAtomsInProc[myNode] = info->getNAtoms();
666     tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
667     tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
668 chrisfen 417
669 gezelter 403 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
670     MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
671     MPI_SUM, MPI_COMM_WORLD);
672     MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
673     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
674     MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
675     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
676 chrisfen 417
677 gezelter 403 beginAtomIndex = 0;
678     beginRigidBodyIndex = 0;
679     beginCutoffGroupIndex = 0;
680 chrisfen 417
681 gezelter 403 for(int i = 0; i < myNode; i++) {
682     beginAtomIndex += NumAtomsInProc[i];
683     beginRigidBodyIndex += NumRigidBodiesInProc[i];
684     beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
685     }
686 chrisfen 417
687 gezelter 403 #endif
688 chrisfen 417
689 gezelter 403 //rigidbody's index begins right after atom's
690     beginRigidBodyIndex += info->getNGlobalAtoms();
691 chrisfen 417
692 gezelter 403 for(mol = info->beginMolecule(mi); mol != NULL;
693     mol = info->nextMolecule(mi)) {
694 chrisfen 417
695 gezelter 403 //local index(index in DataStorge) of atom is important
696     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
697     atom->setGlobalIndex(beginAtomIndex++);
698     }
699 chrisfen 417
700 gezelter 403 for(rb = mol->beginRigidBody(ri); rb != NULL;
701     rb = mol->nextRigidBody(ri)) {
702     rb->setGlobalIndex(beginRigidBodyIndex++);
703     }
704 chrisfen 417
705 gezelter 403 //local index of cutoff group is trivial, it only depends on the order of travesing
706     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
707     cg = mol->nextCutoffGroup(ci)) {
708     cg->setGlobalIndex(beginCutoffGroupIndex++);
709     }
710     }
711 chrisfen 417
712 gezelter 403 //fill globalGroupMembership
713     std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
714     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
715     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
716 chrisfen 417
717 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
718     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
719     }
720 chrisfen 417
721 gezelter 403 }
722     }
723 chrisfen 417
724 gezelter 403 #ifdef IS_MPI
725     // Since the globalGroupMembership has been zero filled and we've only
726     // poked values into the atoms we know, we can do an Allreduce
727     // to get the full globalGroupMembership array (We think).
728     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
729     // docs said we could.
730 gezelter 1313 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
731 gezelter 403 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
732     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
733     info->setGlobalGroupMembership(tmpGroupMembership);
734     #else
735     info->setGlobalGroupMembership(globalGroupMembership);
736     #endif
737 chrisfen 417
738 gezelter 403 //fill molMembership
739     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
740    
741     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
742     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
743     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
744     }
745     }
746 chrisfen 417
747 gezelter 403 #ifdef IS_MPI
748 gezelter 1313 std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
749 chrisfen 417
750 gezelter 403 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
751     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
752    
753     info->setGlobalMolMembership(tmpMolMembership);
754     #else
755     info->setGlobalMolMembership(globalMolMembership);
756     #endif
757 tim 1024
758     // nIOPerMol holds the number of integrable objects per molecule
759     // here the molecules are listed by their global indices.
760    
761     std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
762     for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
763     nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
764     }
765 chrisfen 417
766 tim 1024 #ifdef IS_MPI
767     std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
768     MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
769     info->getNGlobalMolecules(), MPI_INT, MPI_SUM, MPI_COMM_WORLD);
770     #else
771     std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
772     #endif
773    
774 gezelter 1313 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
775    
776     int startingIndex = 0;
777     for (int i = 0; i < info->getNGlobalMolecules(); i++) {
778     startingIOIndexForMol[i] = startingIndex;
779     startingIndex += numIntegrableObjectsPerMol[i];
780     }
781    
782     std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
783     for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
784 tim 1024 int myGlobalIndex = mol->getGlobalIndex();
785     int globalIO = startingIOIndexForMol[myGlobalIndex];
786     for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
787     integrableObject = mol->nextIntegrableObject(ioi)) {
788 gezelter 1313 integrableObject->setGlobalIntegrableObjectIndex(globalIO);
789     IOIndexToIntegrableObject[globalIO] = integrableObject;
790     globalIO++;
791 tim 1024 }
792     }
793 gezelter 1313
794     info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
795    
796 gezelter 403 }
797 chrisfen 417
798 tim 1024 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
799 gezelter 403 Globals* simParams;
800     simParams = info->getSimParams();
801    
802 chrisfen 417
803 tim 1024 DumpReader reader(info, mdFileName);
804 gezelter 403 int nframes = reader.getNFrames();
805 chrisfen 417
806 gezelter 403 if (nframes > 0) {
807     reader.readFrame(nframes - 1);
808     } else {
809     //invalid initial coordinate file
810 chrisfen 417 sprintf(painCave.errMsg,
811     "Initial configuration file %s should at least contain one frame\n",
812 tim 1024 mdFileName.c_str());
813 gezelter 403 painCave.isFatal = 1;
814     simError();
815     }
816 chrisfen 417
817 gezelter 403 //copy the current snapshot to previous snapshot
818     info->getSnapshotManager()->advance();
819     }
820 chrisfen 417
821 gezelter 1390 } //end namespace OpenMD
822 gezelter 403
823    

Properties

Name Value
svn:keywords Author Id Revision Date