ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1025
Committed: Wed Aug 30 20:33:44 2006 UTC (18 years, 8 months ago) by gezelter
Original Path: trunk/src/brains/SimCreator.cpp
File size: 26930 byte(s)
Log Message:
fixed some bugs in DumpWriter, eliminated some extra printing of
debugging information

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     * 1. Acknowledgement of the program authors must be made in any
10     * publication of scientific results based in part on use of the
11     * program. An acceptable form of acknowledgement is citation of
12     * the article in which the program was described (Matthew
13     * A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher
14     * J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented
15     * Parallel Simulation Engine for Molecular Dynamics,"
16     * J. Comput. Chem. 26, pp. 252-271 (2005))
17     *
18     * 2. Redistributions of source code must retain the above copyright
19     * notice, this list of conditions and the following disclaimer.
20     *
21     * 3. Redistributions in binary form must reproduce the above copyright
22     * notice, this list of conditions and the following disclaimer in the
23     * documentation and/or other materials provided with the
24     * distribution.
25     *
26     * This software is provided "AS IS," without a warranty of any
27     * kind. All express or implied conditions, representations and
28     * warranties, including any implied warranty of merchantability,
29     * fitness for a particular purpose or non-infringement, are hereby
30     * excluded. The University of Notre Dame and its licensors shall not
31     * be liable for any damages suffered by licensee as a result of
32     * using, modifying or distributing the software or its
33     * derivatives. In no event will the University of Notre Dame or its
34     * licensors be liable for any lost revenue, profit or data, or for
35     * direct, indirect, special, consequential, incidental or punitive
36     * damages, however caused and regardless of the theory of liability,
37     * arising out of the use of or inability to use software, even if the
38     * University of Notre Dame has been advised of the possibility of
39     * such damages.
40     */
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     namespace oopse {
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 tim 770
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     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 tim 770
149 tim 1024 }
150 tim 845
151 tim 816
152 tim 1024 catch(antlr::MismatchedCharException& e) {
153     sprintf(painCave.errMsg,
154     "parser exception: %s %s:%d:%d\n",
155     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
156     painCave.isFatal = 1;
157     simError();
158     }
159     catch(antlr::MismatchedTokenException &e) {
160     sprintf(painCave.errMsg,
161     "parser exception: %s %s:%d:%d\n",
162     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
163     painCave.isFatal = 1;
164     simError();
165     }
166     catch(antlr::NoViableAltForCharException &e) {
167     sprintf(painCave.errMsg,
168     "parser exception: %s %s:%d:%d\n",
169     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
170     painCave.isFatal = 1;
171     simError();
172     }
173     catch(antlr::NoViableAltException &e) {
174     sprintf(painCave.errMsg,
175     "parser exception: %s %s:%d:%d\n",
176     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
177     painCave.isFatal = 1;
178     simError();
179     }
180 tim 845
181 tim 1024 catch(antlr::TokenStreamRecognitionException& e) {
182     sprintf(painCave.errMsg,
183     "parser exception: %s %s:%d:%d\n",
184     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
185     painCave.isFatal = 1;
186     simError();
187     }
188 tim 845
189 tim 1024 catch(antlr::TokenStreamIOException& e) {
190     sprintf(painCave.errMsg,
191     "parser exception: %s\n",
192     e.getMessage().c_str());
193     painCave.isFatal = 1;
194     simError();
195     }
196 tim 845
197 tim 1024 catch(antlr::TokenStreamException& e) {
198     sprintf(painCave.errMsg,
199     "parser exception: %s\n",
200     e.getMessage().c_str());
201     painCave.isFatal = 1;
202     simError();
203     }
204     catch (antlr::RecognitionException& e) {
205     sprintf(painCave.errMsg,
206     "parser exception: %s %s:%d:%d\n",
207     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
208     painCave.isFatal = 1;
209     simError();
210     }
211     catch (antlr::CharStreamException& e) {
212     sprintf(painCave.errMsg,
213     "parser exception: %s\n",
214     e.getMessage().c_str());
215     painCave.isFatal = 1;
216     simError();
217     }
218     catch (OOPSEException& e) {
219     sprintf(painCave.errMsg,
220     "%s\n",
221     e.getMessage().c_str());
222     painCave.isFatal = 1;
223     simError();
224     }
225     catch (std::exception& e) {
226     sprintf(painCave.errMsg,
227     "parser exception: %s\n",
228     e.what());
229     painCave.isFatal = 1;
230     simError();
231     }
232 tim 770
233 tim 1024 return simParams;
234 gezelter 403 }
235 chrisfen 417
236 chrisfen 514 SimInfo* SimCreator::createSim(const std::string & mdFileName,
237     bool loadInitCoords) {
238 tim 770
239 tim 1024 const int bufferSize = 65535;
240     char buffer[bufferSize];
241     int lineNo = 0;
242     std::string mdRawData;
243     int metaDataBlockStart = -1;
244     int metaDataBlockEnd = -1;
245     int i;
246     int mdOffset;
247    
248     #ifdef IS_MPI
249     const int masterNode = 0;
250     if (worldRank == masterNode) {
251     #endif
252    
253     std::ifstream mdFile_(mdFileName.c_str());
254    
255     if (mdFile_.fail()) {
256     sprintf(painCave.errMsg,
257     "SimCreator: Cannot open file: %s\n",
258     mdFileName.c_str());
259     painCave.isFatal = 1;
260     simError();
261     }
262    
263     mdFile_.getline(buffer, bufferSize);
264     ++lineNo;
265     std::string line = trimLeftCopy(buffer);
266     i = CaseInsensitiveFind(line, "<OOPSE");
267     if (i == string::npos) {
268     sprintf(painCave.errMsg,
269     "SimCreator: File: %s is not an OOPSE file!\n",
270     mdFileName.c_str());
271     painCave.isFatal = 1;
272     simError();
273     }
274    
275     //scan through the input stream and find MetaData tag
276     while(mdFile_.getline(buffer, bufferSize)) {
277     ++lineNo;
278    
279     std::string line = trimLeftCopy(buffer);
280     if (metaDataBlockStart == -1) {
281     i = CaseInsensitiveFind(line, "<MetaData>");
282     if (i != string::npos) {
283     metaDataBlockStart = lineNo;
284     mdOffset = mdFile_.tellg();
285     }
286     } else {
287     i = CaseInsensitiveFind(line, "</MetaData>");
288     if (i != string::npos) {
289     metaDataBlockEnd = lineNo;
290     }
291     }
292     }
293    
294     if (metaDataBlockStart == -1) {
295     sprintf(painCave.errMsg,
296     "SimCreator: File: %s did not contain a <MetaData> tag!\n",
297     mdFileName.c_str());
298     painCave.isFatal = 1;
299     simError();
300     }
301     if (metaDataBlockEnd == -1) {
302     sprintf(painCave.errMsg,
303     "SimCreator: File: %s did not contain a closed MetaData block!\n",
304     mdFileName.c_str());
305     painCave.isFatal = 1;
306     simError();
307     }
308    
309     mdFile_.clear();
310     mdFile_.seekg(0);
311     mdFile_.seekg(mdOffset);
312    
313     mdRawData.clear();
314    
315     for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
316     mdFile_.getline(buffer, bufferSize);
317     mdRawData += buffer;
318     mdRawData += "\n";
319     }
320    
321     mdFile_.close();
322    
323     #ifdef IS_MPI
324     }
325     #endif
326    
327     std::stringstream rawMetaDataStream(mdRawData);
328    
329 gezelter 403 //parse meta-data file
330 tim 1024 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, metaDataBlockStart+1);
331 chrisfen 417
332 gezelter 403 //create the force field
333 chrisfen 514 ForceField * ff = ForceFieldFactory::getInstance()
334     ->createForceField(simParams->getForceField());
335 gezelter 403
336     if (ff == NULL) {
337 chrisfen 514 sprintf(painCave.errMsg,
338     "ForceField Factory can not create %s force field\n",
339 tim 665 simParams->getForceField().c_str());
340 gezelter 403 painCave.isFatal = 1;
341     simError();
342     }
343 chrisfen 417
344 gezelter 403 if (simParams->haveForceFieldFileName()) {
345     ff->setForceFieldFileName(simParams->getForceFieldFileName());
346     }
347    
348     std::string forcefieldFileName;
349     forcefieldFileName = ff->getForceFieldFileName();
350 chrisfen 417
351 gezelter 403 if (simParams->haveForceFieldVariant()) {
352     //If the force field has variant, the variant force field name will be
353     //Base.variant.frc. For exampel EAM.u6.frc
354 chrisfen 417
355 gezelter 403 std::string variant = simParams->getForceFieldVariant();
356 chrisfen 417
357 gezelter 403 std::string::size_type pos = forcefieldFileName.rfind(".frc");
358     variant = "." + variant;
359     if (pos != std::string::npos) {
360     forcefieldFileName.insert(pos, variant);
361     } else {
362     //If the default force field file name does not containt .frc suffix, just append the .variant
363     forcefieldFileName.append(variant);
364     }
365     }
366    
367     ff->parse(forcefieldFileName);
368 chuckv 823 ff->setFortranForceOptions();
369 gezelter 403 //create SimInfo
370 tim 770 SimInfo * info = new SimInfo(ff, simParams);
371 tim 1024
372     info->setRawMetaData(mdRawData);
373 tim 490
374 gezelter 945 //gather parameters (SimCreator only retrieves part of the
375     //parameters)
376 gezelter 403 gatherParameters(info, mdFileName);
377 chrisfen 417
378 gezelter 403 //divide the molecules and determine the global index of molecules
379     #ifdef IS_MPI
380     divideMolecules(info);
381     #endif
382 chrisfen 417
383 gezelter 403 //create the molecules
384     createMolecules(info);
385 chrisfen 417
386    
387 gezelter 945 //allocate memory for DataStorage(circular reference, need to
388     //break it)
389 gezelter 403 info->setSnapshotManager(new SimSnapshotManager(info));
390    
391 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
392     //(only need to be set once, the global index will never change
393     //again). Local indices of atoms and rigidbodies are already set
394     //by MoleculeCreator class which actually delegates the
395     //responsibility to LocalIndexManager.
396 gezelter 403 setGlobalIndex(info);
397 chrisfen 417
398 gezelter 945 //Although addExcludePairs is called inside SimInfo's addMolecule
399     //method, at that point atoms don't have the global index yet
400     //(their global index are all initialized to -1). Therefore we
401     //have to call addExcludePairs explicitly here. A way to work
402     //around is that we can determine the beginning global indices of
403     //atoms before they get created.
404 gezelter 403 SimInfo::MoleculeIterator mi;
405     Molecule* mol;
406     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
407     info->addExcludePairs(mol);
408     }
409    
410     if (loadInitCoords)
411 tim 1024 loadCoordinates(info, mdFileName);
412 gezelter 403
413     return info;
414     }
415 chrisfen 417
416 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
417 chrisfen 417
418 tim 749 //figure out the output file names
419 gezelter 403 std::string prefix;
420 chrisfen 417
421 gezelter 403 #ifdef IS_MPI
422 chrisfen 417
423 gezelter 403 if (worldRank == 0) {
424     #endif // is_mpi
425     Globals * simParams = info->getSimParams();
426     if (simParams->haveFinalConfig()) {
427     prefix = getPrefix(simParams->getFinalConfig());
428     } else {
429     prefix = getPrefix(mdfile);
430     }
431 chrisfen 417
432 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
433     info->setDumpFileName(prefix + ".dump");
434     info->setStatFileName(prefix + ".stat");
435 chrisfen 417 info->setRestFileName(prefix + ".zang");
436    
437 gezelter 403 #ifdef IS_MPI
438 chrisfen 417
439 gezelter 403 }
440 chrisfen 417
441 gezelter 403 #endif
442 chrisfen 417
443 gezelter 403 }
444 chrisfen 417
445 gezelter 403 #ifdef IS_MPI
446     void SimCreator::divideMolecules(SimInfo *info) {
447 tim 963 RealType numerator;
448     RealType denominator;
449     RealType precast;
450     RealType x;
451     RealType y;
452     RealType a;
453 gezelter 403 int old_atoms;
454     int add_atoms;
455     int new_atoms;
456     int nTarget;
457     int done;
458     int i;
459     int j;
460     int loops;
461     int which_proc;
462     int nProcessors;
463     std::vector<int> atomsPerProc;
464     int nGlobalMols = info->getNGlobalMolecules();
465     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
466    
467     MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
468 chrisfen 417
469 gezelter 403 if (nProcessors > nGlobalMols) {
470     sprintf(painCave.errMsg,
471     "nProcessors (%d) > nMol (%d)\n"
472     "\tThe number of processors is larger than\n"
473     "\tthe number of molecules. This will not result in a \n"
474     "\tusable division of atoms for force decomposition.\n"
475     "\tEither try a smaller number of processors, or run the\n"
476     "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
477 chrisfen 417
478 gezelter 403 painCave.isFatal = 1;
479     simError();
480     }
481 chrisfen 417
482 gezelter 403 int seedValue;
483     Globals * simParams = info->getSimParams();
484     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
485     if (simParams->haveSeed()) {
486     seedValue = simParams->getSeed();
487     myRandom = new SeqRandNumGen(seedValue);
488     }else {
489     myRandom = new SeqRandNumGen();
490     }
491 chrisfen 417
492    
493 gezelter 403 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
494 chrisfen 417
495 gezelter 403 //initialize atomsPerProc
496     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
497 chrisfen 417
498 gezelter 403 if (worldRank == 0) {
499     numerator = info->getNGlobalAtoms();
500     denominator = nProcessors;
501     precast = numerator / denominator;
502     nTarget = (int)(precast + 0.5);
503 chrisfen 417
504 gezelter 403 for(i = 0; i < nGlobalMols; i++) {
505     done = 0;
506     loops = 0;
507 chrisfen 417
508 gezelter 403 while (!done) {
509     loops++;
510 chrisfen 417
511 gezelter 403 // Pick a processor at random
512 chrisfen 417
513 gezelter 403 which_proc = (int) (myRandom->rand() * nProcessors);
514 chrisfen 417
515 gezelter 403 //get the molecule stamp first
516     int stampId = info->getMoleculeStampId(i);
517     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
518 chrisfen 417
519 gezelter 403 // How many atoms does this processor have so far?
520     old_atoms = atomsPerProc[which_proc];
521     add_atoms = moleculeStamp->getNAtoms();
522     new_atoms = old_atoms + add_atoms;
523 chrisfen 417
524 gezelter 403 // If we've been through this loop too many times, we need
525     // to just give up and assign the molecule to this processor
526     // and be done with it.
527 chrisfen 417
528 gezelter 403 if (loops > 100) {
529     sprintf(painCave.errMsg,
530     "I've tried 100 times to assign molecule %d to a "
531     " processor, but can't find a good spot.\n"
532     "I'm assigning it at random to processor %d.\n",
533     i, which_proc);
534 chrisfen 417
535 gezelter 403 painCave.isFatal = 0;
536     simError();
537 chrisfen 417
538 gezelter 403 molToProcMap[i] = which_proc;
539     atomsPerProc[which_proc] += add_atoms;
540 chrisfen 417
541 gezelter 403 done = 1;
542     continue;
543     }
544 chrisfen 417
545 gezelter 403 // If we can add this molecule to this processor without sending
546     // it above nTarget, then go ahead and do it:
547 chrisfen 417
548 gezelter 403 if (new_atoms <= nTarget) {
549     molToProcMap[i] = which_proc;
550     atomsPerProc[which_proc] += add_atoms;
551 chrisfen 417
552 gezelter 403 done = 1;
553     continue;
554     }
555 chrisfen 417
556 gezelter 403 // The only situation left is when new_atoms > nTarget. We
557     // want to accept this with some probability that dies off the
558     // farther we are from nTarget
559 chrisfen 417
560 gezelter 403 // roughly: x = new_atoms - nTarget
561     // Pacc(x) = exp(- a * x)
562     // where a = penalty / (average atoms per molecule)
563 chrisfen 417
564 tim 963 x = (RealType)(new_atoms - nTarget);
565 gezelter 403 y = myRandom->rand();
566 chrisfen 417
567 gezelter 403 if (y < exp(- a * x)) {
568     molToProcMap[i] = which_proc;
569     atomsPerProc[which_proc] += add_atoms;
570 chrisfen 417
571 gezelter 403 done = 1;
572     continue;
573     } else {
574     continue;
575     }
576     }
577     }
578 chrisfen 417
579 gezelter 403 delete myRandom;
580 chrisfen 417
581 gezelter 403 // Spray out this nonsense to all other processors:
582 chrisfen 417
583 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
584     } else {
585 chrisfen 417
586 gezelter 403 // Listen to your marching orders from processor 0:
587 chrisfen 417
588 gezelter 403 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
589     }
590 chrisfen 417
591 gezelter 403 info->setMolToProcMap(molToProcMap);
592     sprintf(checkPointMsg,
593     "Successfully divided the molecules among the processors.\n");
594     MPIcheckPoint();
595     }
596 chrisfen 417
597 gezelter 403 #endif
598 chrisfen 417
599 gezelter 403 void SimCreator::createMolecules(SimInfo *info) {
600     MoleculeCreator molCreator;
601     int stampId;
602 chrisfen 417
603 gezelter 403 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
604 chrisfen 417
605 gezelter 403 #ifdef IS_MPI
606 chrisfen 417
607 gezelter 403 if (info->getMolToProc(i) == worldRank) {
608     #endif
609 chrisfen 417
610 gezelter 403 stampId = info->getMoleculeStampId(i);
611     Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
612     stampId, i, info->getLocalIndexManager());
613 chrisfen 417
614 gezelter 403 info->addMolecule(mol);
615 chrisfen 417
616 gezelter 403 #ifdef IS_MPI
617 chrisfen 417
618 gezelter 403 }
619 chrisfen 417
620 gezelter 403 #endif
621 chrisfen 417
622 gezelter 403 } //end for(int i=0)
623     }
624 chrisfen 417
625 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
626     SimInfo::MoleculeIterator mi;
627     Molecule::AtomIterator ai;
628     Molecule::RigidBodyIterator ri;
629     Molecule::CutoffGroupIterator ci;
630 tim 1024 Molecule::IntegrableObjectIterator ioi;
631 gezelter 403 Molecule * mol;
632     Atom * atom;
633     RigidBody * rb;
634     CutoffGroup * cg;
635     int beginAtomIndex;
636     int beginRigidBodyIndex;
637     int beginCutoffGroupIndex;
638     int nGlobalAtoms = info->getNGlobalAtoms();
639 tim 1024
640     /**@todo fixme */
641 gezelter 403 #ifndef IS_MPI
642 chrisfen 417
643 gezelter 403 beginAtomIndex = 0;
644     beginRigidBodyIndex = 0;
645     beginCutoffGroupIndex = 0;
646 chrisfen 417
647 gezelter 403 #else
648 chrisfen 417
649 gezelter 403 int nproc;
650     int myNode;
651 chrisfen 417
652 gezelter 403 myNode = worldRank;
653     MPI_Comm_size(MPI_COMM_WORLD, &nproc);
654 chrisfen 417
655 gezelter 403 std::vector < int > tmpAtomsInProc(nproc, 0);
656     std::vector < int > tmpRigidBodiesInProc(nproc, 0);
657     std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
658     std::vector < int > NumAtomsInProc(nproc, 0);
659     std::vector < int > NumRigidBodiesInProc(nproc, 0);
660     std::vector < int > NumCutoffGroupsInProc(nproc, 0);
661 chrisfen 417
662 gezelter 403 tmpAtomsInProc[myNode] = info->getNAtoms();
663     tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
664     tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
665 chrisfen 417
666 gezelter 403 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
667     MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
668     MPI_SUM, MPI_COMM_WORLD);
669     MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
670     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
671     MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
672     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
673 chrisfen 417
674 gezelter 403 beginAtomIndex = 0;
675     beginRigidBodyIndex = 0;
676     beginCutoffGroupIndex = 0;
677 chrisfen 417
678 gezelter 403 for(int i = 0; i < myNode; i++) {
679     beginAtomIndex += NumAtomsInProc[i];
680     beginRigidBodyIndex += NumRigidBodiesInProc[i];
681     beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
682     }
683 chrisfen 417
684 gezelter 403 #endif
685 chrisfen 417
686 gezelter 403 //rigidbody's index begins right after atom's
687     beginRigidBodyIndex += info->getNGlobalAtoms();
688 chrisfen 417
689 gezelter 403 for(mol = info->beginMolecule(mi); mol != NULL;
690     mol = info->nextMolecule(mi)) {
691 chrisfen 417
692 gezelter 403 //local index(index in DataStorge) of atom is important
693     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
694     atom->setGlobalIndex(beginAtomIndex++);
695     }
696 chrisfen 417
697 gezelter 403 for(rb = mol->beginRigidBody(ri); rb != NULL;
698     rb = mol->nextRigidBody(ri)) {
699     rb->setGlobalIndex(beginRigidBodyIndex++);
700     }
701 chrisfen 417
702 gezelter 403 //local index of cutoff group is trivial, it only depends on the order of travesing
703     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
704     cg = mol->nextCutoffGroup(ci)) {
705     cg->setGlobalIndex(beginCutoffGroupIndex++);
706     }
707     }
708 chrisfen 417
709 gezelter 403 //fill globalGroupMembership
710     std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
711     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
712     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
713 chrisfen 417
714 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
715     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
716     }
717 chrisfen 417
718 gezelter 403 }
719     }
720 chrisfen 417
721 gezelter 403 #ifdef IS_MPI
722     // Since the globalGroupMembership has been zero filled and we've only
723     // poked values into the atoms we know, we can do an Allreduce
724     // to get the full globalGroupMembership array (We think).
725     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
726     // docs said we could.
727     std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
728     MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
729     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
730     info->setGlobalGroupMembership(tmpGroupMembership);
731     #else
732     info->setGlobalGroupMembership(globalGroupMembership);
733     #endif
734 chrisfen 417
735 gezelter 403 //fill molMembership
736     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
737    
738     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
739 chrisfen 417
740 gezelter 403 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     std::vector<int> tmpMolMembership(nGlobalAtoms, 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     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     int myGlobalIndex = mol->getGlobalIndex();
783     int globalIO = startingIOIndexForMol[myGlobalIndex];
784     for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
785     integrableObject = mol->nextIntegrableObject(ioi)) {
786     integrableObject->setGlobalIntegrableObjectIndex(globalIO);
787     IOIndexToIntegrableObject[globalIO] = integrableObject;
788     globalIO++;
789     }
790     }
791    
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     simParams = info->getSimParams();
799    
800 chrisfen 417
801 tim 1024 DumpReader reader(info, mdFileName);
802 gezelter 403 int nframes = reader.getNFrames();
803 chrisfen 417
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 chrisfen 417
815 gezelter 403 //copy the current snapshot to previous snapshot
816     info->getSnapshotManager()->advance();
817     }
818 chrisfen 417
819 gezelter 403 } //end namespace oopse
820    
821