ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimCreator.cpp
Revision: 2759
Committed: Wed May 17 21:51:42 2006 UTC (18 years, 3 months ago) by tim
File size: 23341 byte(s)
Log Message:
Adding single precision capabilities to c++ side

File Contents

# User Rev Content
1 gezelter 2087 /*
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 2544 #include <exception>
50 tim 2469 #include <iostream>
51     #include <sstream>
52     #include <string>
53    
54 gezelter 2087 #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 2469 #include "mdParser/MDLexer.hpp"
63     #include "mdParser/MDParser.hpp"
64     #include "mdParser/MDTreeParser.hpp"
65     #include "mdParser/SimplePreprocessor.hpp"
66 tim 2515 #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 2469
73 tim 2515 #include "antlr/MismatchedCharException.hpp"
74     #include "antlr/MismatchedTokenException.hpp"
75     #include "antlr/NoViableAltForCharException.hpp"
76     #include "antlr/NoViableAltException.hpp"
77 tim 2469
78 gezelter 2087 #ifdef IS_MPI
79     #include "math/ParallelRandNumGen.hpp"
80     #endif
81    
82     namespace oopse {
83 chrisfen 2101
84 tim 2469 Globals* SimCreator::parseFile(const std::string mdFileName){
85     Globals* simParams = NULL;
86     try {
87    
88     // Create a preprocessor that preprocesses md file into an ostringstream
89     std::stringstream ppStream;
90     #ifdef IS_MPI
91     int streamSize;
92     const int masterNode = 0;
93     int commStatus;
94     if (worldRank == masterNode) {
95     #endif
96    
97     SimplePreprocessor preprocessor;
98     preprocessor.preprocess(mdFileName, ppStream);
99    
100     #ifdef IS_MPI
101     //brocasting the stream size
102     streamSize = ppStream.str().size() +1;
103     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
104    
105 tim 2531 commStatus = MPI_Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())), streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
106 tim 2469
107    
108     } else {
109     //get stream size
110     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
111    
112     char* buf = new char[streamSize];
113     assert(buf);
114    
115     //receive file content
116     commStatus = MPI_Bcast(buf, streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
117    
118     ppStream.str(buf);
119     delete buf;
120    
121     }
122     #endif
123     // Create a scanner that reads from the input stream
124     MDLexer lexer(ppStream);
125     lexer.setFilename(mdFileName);
126     lexer.initDeferredLineCount();
127 chrisfen 2101
128 tim 2469 // Create a parser that reads from the scanner
129     MDParser parser(lexer);
130     parser.setFilename(mdFileName);
131    
132     // 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 2101
138 tim 2469 antlr::ASTFactory factory;
139     parser.initializeASTFactory(factory);
140     parser.setASTFactory(&factory);
141     parser.mdfile();
142    
143     // 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     }
150 tim 2544
151 tim 2515
152     catch(antlr::MismatchedCharException& e) {
153 tim 2544 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 tim 2515 }
159     catch(antlr::MismatchedTokenException &e) {
160 tim 2544 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 tim 2515 }
166     catch(antlr::NoViableAltForCharException &e) {
167 tim 2544 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 tim 2515 }
173     catch(antlr::NoViableAltException &e) {
174 tim 2544 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 tim 2515 }
180 tim 2544
181 tim 2515 catch(antlr::TokenStreamRecognitionException& e) {
182 tim 2544 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 tim 2515 }
188 tim 2544
189 tim 2515 catch(antlr::TokenStreamIOException& e) {
190 tim 2544 sprintf(painCave.errMsg,
191     "parser exception: %s\n",
192     e.getMessage().c_str());
193     painCave.isFatal = 1;
194     simError();
195 tim 2515 }
196 tim 2544
197 tim 2515 catch(antlr::TokenStreamException& e) {
198 tim 2544 sprintf(painCave.errMsg,
199     "parser exception: %s\n",
200     e.getMessage().c_str());
201     painCave.isFatal = 1;
202     simError();
203 tim 2515 }
204     catch (antlr::RecognitionException& e) {
205 tim 2544 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 tim 2515 }
211     catch (antlr::CharStreamException& e) {
212 tim 2544 sprintf(painCave.errMsg,
213     "parser exception: %s\n",
214     e.getMessage().c_str());
215     painCave.isFatal = 1;
216     simError();
217 tim 2515 }
218 tim 2544 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 2469
233     return simParams;
234 gezelter 2087 }
235 chrisfen 2101
236 chrisfen 2211 SimInfo* SimCreator::createSim(const std::string & mdFileName,
237     bool loadInitCoords) {
238 tim 2469
239 gezelter 2087 //parse meta-data file
240 tim 2469 Globals* simParams = parseFile(mdFileName);
241 chrisfen 2101
242 gezelter 2087 //create the force field
243 chrisfen 2211 ForceField * ff = ForceFieldFactory::getInstance()
244     ->createForceField(simParams->getForceField());
245 gezelter 2087
246     if (ff == NULL) {
247 chrisfen 2211 sprintf(painCave.errMsg,
248     "ForceField Factory can not create %s force field\n",
249 tim 2364 simParams->getForceField().c_str());
250 gezelter 2087 painCave.isFatal = 1;
251     simError();
252     }
253 chrisfen 2101
254 gezelter 2087 if (simParams->haveForceFieldFileName()) {
255     ff->setForceFieldFileName(simParams->getForceFieldFileName());
256     }
257    
258     std::string forcefieldFileName;
259     forcefieldFileName = ff->getForceFieldFileName();
260 chrisfen 2101
261 gezelter 2087 if (simParams->haveForceFieldVariant()) {
262     //If the force field has variant, the variant force field name will be
263     //Base.variant.frc. For exampel EAM.u6.frc
264 chrisfen 2101
265 gezelter 2087 std::string variant = simParams->getForceFieldVariant();
266 chrisfen 2101
267 gezelter 2087 std::string::size_type pos = forcefieldFileName.rfind(".frc");
268     variant = "." + variant;
269     if (pos != std::string::npos) {
270     forcefieldFileName.insert(pos, variant);
271     } else {
272     //If the default force field file name does not containt .frc suffix, just append the .variant
273     forcefieldFileName.append(variant);
274     }
275     }
276    
277     ff->parse(forcefieldFileName);
278 chuckv 2522 ff->setFortranForceOptions();
279 gezelter 2087 //create SimInfo
280 tim 2469 SimInfo * info = new SimInfo(ff, simParams);
281 tim 2187
282 gezelter 2733 //gather parameters (SimCreator only retrieves part of the
283     //parameters)
284 gezelter 2087 gatherParameters(info, mdFileName);
285 chrisfen 2101
286 gezelter 2087 //divide the molecules and determine the global index of molecules
287     #ifdef IS_MPI
288     divideMolecules(info);
289     #endif
290 chrisfen 2101
291 gezelter 2087 //create the molecules
292     createMolecules(info);
293 chrisfen 2101
294    
295 gezelter 2733 //allocate memory for DataStorage(circular reference, need to
296     //break it)
297 gezelter 2087 info->setSnapshotManager(new SimSnapshotManager(info));
298    
299 gezelter 2733 //set the global index of atoms, rigidbodies and cutoffgroups
300     //(only need to be set once, the global index will never change
301     //again). Local indices of atoms and rigidbodies are already set
302     //by MoleculeCreator class which actually delegates the
303     //responsibility to LocalIndexManager.
304 gezelter 2087 setGlobalIndex(info);
305 chrisfen 2101
306 gezelter 2733 //Although addExcludePairs is called inside SimInfo's addMolecule
307     //method, at that point atoms don't have the global index yet
308     //(their global index are all initialized to -1). Therefore we
309     //have to call addExcludePairs explicitly here. A way to work
310     //around is that we can determine the beginning global indices of
311     //atoms before they get created.
312 gezelter 2087 SimInfo::MoleculeIterator mi;
313     Molecule* mol;
314     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
315     info->addExcludePairs(mol);
316     }
317    
318     if (loadInitCoords)
319     loadCoordinates(info);
320    
321     return info;
322     }
323 chrisfen 2101
324 gezelter 2087 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
325 chrisfen 2101
326 tim 2448 //figure out the output file names
327 gezelter 2087 std::string prefix;
328 chrisfen 2101
329 gezelter 2087 #ifdef IS_MPI
330 chrisfen 2101
331 gezelter 2087 if (worldRank == 0) {
332     #endif // is_mpi
333     Globals * simParams = info->getSimParams();
334     if (simParams->haveFinalConfig()) {
335     prefix = getPrefix(simParams->getFinalConfig());
336     } else {
337     prefix = getPrefix(mdfile);
338     }
339 chrisfen 2101
340 gezelter 2087 info->setFinalConfigFileName(prefix + ".eor");
341     info->setDumpFileName(prefix + ".dump");
342     info->setStatFileName(prefix + ".stat");
343 chrisfen 2101 info->setRestFileName(prefix + ".zang");
344    
345 gezelter 2087 #ifdef IS_MPI
346 chrisfen 2101
347 gezelter 2087 }
348 chrisfen 2101
349 gezelter 2087 #endif
350 chrisfen 2101
351 gezelter 2087 }
352 chrisfen 2101
353 gezelter 2087 #ifdef IS_MPI
354     void SimCreator::divideMolecules(SimInfo *info) {
355 tim 2759 RealType numerator;
356     RealType denominator;
357     RealType precast;
358     RealType x;
359     RealType y;
360     RealType a;
361 gezelter 2087 int old_atoms;
362     int add_atoms;
363     int new_atoms;
364     int nTarget;
365     int done;
366     int i;
367     int j;
368     int loops;
369     int which_proc;
370     int nProcessors;
371     std::vector<int> atomsPerProc;
372     int nGlobalMols = info->getNGlobalMolecules();
373     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
374    
375     MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
376 chrisfen 2101
377 gezelter 2087 if (nProcessors > nGlobalMols) {
378     sprintf(painCave.errMsg,
379     "nProcessors (%d) > nMol (%d)\n"
380     "\tThe number of processors is larger than\n"
381     "\tthe number of molecules. This will not result in a \n"
382     "\tusable division of atoms for force decomposition.\n"
383     "\tEither try a smaller number of processors, or run the\n"
384     "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
385 chrisfen 2101
386 gezelter 2087 painCave.isFatal = 1;
387     simError();
388     }
389 chrisfen 2101
390 gezelter 2087 int seedValue;
391     Globals * simParams = info->getSimParams();
392     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
393     if (simParams->haveSeed()) {
394     seedValue = simParams->getSeed();
395     myRandom = new SeqRandNumGen(seedValue);
396     }else {
397     myRandom = new SeqRandNumGen();
398     }
399 chrisfen 2101
400    
401 gezelter 2087 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
402 chrisfen 2101
403 gezelter 2087 //initialize atomsPerProc
404     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
405 chrisfen 2101
406 gezelter 2087 if (worldRank == 0) {
407     numerator = info->getNGlobalAtoms();
408     denominator = nProcessors;
409     precast = numerator / denominator;
410     nTarget = (int)(precast + 0.5);
411 chrisfen 2101
412 gezelter 2087 for(i = 0; i < nGlobalMols; i++) {
413     done = 0;
414     loops = 0;
415 chrisfen 2101
416 gezelter 2087 while (!done) {
417     loops++;
418 chrisfen 2101
419 gezelter 2087 // Pick a processor at random
420 chrisfen 2101
421 gezelter 2087 which_proc = (int) (myRandom->rand() * nProcessors);
422 chrisfen 2101
423 gezelter 2087 //get the molecule stamp first
424     int stampId = info->getMoleculeStampId(i);
425     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
426 chrisfen 2101
427 gezelter 2087 // How many atoms does this processor have so far?
428     old_atoms = atomsPerProc[which_proc];
429     add_atoms = moleculeStamp->getNAtoms();
430     new_atoms = old_atoms + add_atoms;
431 chrisfen 2101
432 gezelter 2087 // If we've been through this loop too many times, we need
433     // to just give up and assign the molecule to this processor
434     // and be done with it.
435 chrisfen 2101
436 gezelter 2087 if (loops > 100) {
437     sprintf(painCave.errMsg,
438     "I've tried 100 times to assign molecule %d to a "
439     " processor, but can't find a good spot.\n"
440     "I'm assigning it at random to processor %d.\n",
441     i, which_proc);
442 chrisfen 2101
443 gezelter 2087 painCave.isFatal = 0;
444     simError();
445 chrisfen 2101
446 gezelter 2087 molToProcMap[i] = which_proc;
447     atomsPerProc[which_proc] += add_atoms;
448 chrisfen 2101
449 gezelter 2087 done = 1;
450     continue;
451     }
452 chrisfen 2101
453 gezelter 2087 // If we can add this molecule to this processor without sending
454     // it above nTarget, then go ahead and do it:
455 chrisfen 2101
456 gezelter 2087 if (new_atoms <= nTarget) {
457     molToProcMap[i] = which_proc;
458     atomsPerProc[which_proc] += add_atoms;
459 chrisfen 2101
460 gezelter 2087 done = 1;
461     continue;
462     }
463 chrisfen 2101
464 gezelter 2087 // The only situation left is when new_atoms > nTarget. We
465     // want to accept this with some probability that dies off the
466     // farther we are from nTarget
467 chrisfen 2101
468 gezelter 2087 // roughly: x = new_atoms - nTarget
469     // Pacc(x) = exp(- a * x)
470     // where a = penalty / (average atoms per molecule)
471 chrisfen 2101
472 tim 2759 x = (RealType)(new_atoms - nTarget);
473 gezelter 2087 y = myRandom->rand();
474 chrisfen 2101
475 gezelter 2087 if (y < exp(- a * x)) {
476     molToProcMap[i] = which_proc;
477     atomsPerProc[which_proc] += add_atoms;
478 chrisfen 2101
479 gezelter 2087 done = 1;
480     continue;
481     } else {
482     continue;
483     }
484     }
485     }
486 chrisfen 2101
487 gezelter 2087 delete myRandom;
488 chrisfen 2101
489 gezelter 2087 // Spray out this nonsense to all other processors:
490 chrisfen 2101
491 gezelter 2087 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
492     } else {
493 chrisfen 2101
494 gezelter 2087 // Listen to your marching orders from processor 0:
495 chrisfen 2101
496 gezelter 2087 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
497     }
498 chrisfen 2101
499 gezelter 2087 info->setMolToProcMap(molToProcMap);
500     sprintf(checkPointMsg,
501     "Successfully divided the molecules among the processors.\n");
502     MPIcheckPoint();
503     }
504 chrisfen 2101
505 gezelter 2087 #endif
506 chrisfen 2101
507 gezelter 2087 void SimCreator::createMolecules(SimInfo *info) {
508     MoleculeCreator molCreator;
509     int stampId;
510 chrisfen 2101
511 gezelter 2087 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
512 chrisfen 2101
513 gezelter 2087 #ifdef IS_MPI
514 chrisfen 2101
515 gezelter 2087 if (info->getMolToProc(i) == worldRank) {
516     #endif
517 chrisfen 2101
518 gezelter 2087 stampId = info->getMoleculeStampId(i);
519     Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
520     stampId, i, info->getLocalIndexManager());
521 chrisfen 2101
522 gezelter 2087 info->addMolecule(mol);
523 chrisfen 2101
524 gezelter 2087 #ifdef IS_MPI
525 chrisfen 2101
526 gezelter 2087 }
527 chrisfen 2101
528 gezelter 2087 #endif
529 chrisfen 2101
530 gezelter 2087 } //end for(int i=0)
531     }
532 chrisfen 2101
533 gezelter 2087 void SimCreator::setGlobalIndex(SimInfo *info) {
534     SimInfo::MoleculeIterator mi;
535     Molecule::AtomIterator ai;
536     Molecule::RigidBodyIterator ri;
537     Molecule::CutoffGroupIterator ci;
538     Molecule * mol;
539     Atom * atom;
540     RigidBody * rb;
541     CutoffGroup * cg;
542     int beginAtomIndex;
543     int beginRigidBodyIndex;
544     int beginCutoffGroupIndex;
545     int nGlobalAtoms = info->getNGlobalAtoms();
546    
547     #ifndef IS_MPI
548 chrisfen 2101
549 gezelter 2087 beginAtomIndex = 0;
550     beginRigidBodyIndex = 0;
551     beginCutoffGroupIndex = 0;
552 chrisfen 2101
553 gezelter 2087 #else
554 chrisfen 2101
555 gezelter 2087 int nproc;
556     int myNode;
557 chrisfen 2101
558 gezelter 2087 myNode = worldRank;
559     MPI_Comm_size(MPI_COMM_WORLD, &nproc);
560 chrisfen 2101
561 gezelter 2087 std::vector < int > tmpAtomsInProc(nproc, 0);
562     std::vector < int > tmpRigidBodiesInProc(nproc, 0);
563     std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
564     std::vector < int > NumAtomsInProc(nproc, 0);
565     std::vector < int > NumRigidBodiesInProc(nproc, 0);
566     std::vector < int > NumCutoffGroupsInProc(nproc, 0);
567 chrisfen 2101
568 gezelter 2087 tmpAtomsInProc[myNode] = info->getNAtoms();
569     tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
570     tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
571 chrisfen 2101
572 gezelter 2087 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
573     MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
574     MPI_SUM, MPI_COMM_WORLD);
575     MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
576     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
577     MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
578     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
579 chrisfen 2101
580 gezelter 2087 beginAtomIndex = 0;
581     beginRigidBodyIndex = 0;
582     beginCutoffGroupIndex = 0;
583 chrisfen 2101
584 gezelter 2087 for(int i = 0; i < myNode; i++) {
585     beginAtomIndex += NumAtomsInProc[i];
586     beginRigidBodyIndex += NumRigidBodiesInProc[i];
587     beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
588     }
589 chrisfen 2101
590 gezelter 2087 #endif
591 chrisfen 2101
592 gezelter 2087 //rigidbody's index begins right after atom's
593     beginRigidBodyIndex += info->getNGlobalAtoms();
594 chrisfen 2101
595 gezelter 2087 for(mol = info->beginMolecule(mi); mol != NULL;
596     mol = info->nextMolecule(mi)) {
597 chrisfen 2101
598 gezelter 2087 //local index(index in DataStorge) of atom is important
599     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
600     atom->setGlobalIndex(beginAtomIndex++);
601     }
602 chrisfen 2101
603 gezelter 2087 for(rb = mol->beginRigidBody(ri); rb != NULL;
604     rb = mol->nextRigidBody(ri)) {
605     rb->setGlobalIndex(beginRigidBodyIndex++);
606     }
607 chrisfen 2101
608 gezelter 2087 //local index of cutoff group is trivial, it only depends on the order of travesing
609     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
610     cg = mol->nextCutoffGroup(ci)) {
611     cg->setGlobalIndex(beginCutoffGroupIndex++);
612     }
613     }
614 chrisfen 2101
615 gezelter 2087 //fill globalGroupMembership
616     std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
617     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
618     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
619 chrisfen 2101
620 gezelter 2087 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
621     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
622     }
623 chrisfen 2101
624 gezelter 2087 }
625     }
626 chrisfen 2101
627 gezelter 2087 #ifdef IS_MPI
628     // Since the globalGroupMembership has been zero filled and we've only
629     // poked values into the atoms we know, we can do an Allreduce
630     // to get the full globalGroupMembership array (We think).
631     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
632     // docs said we could.
633     std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
634     MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
635     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
636     info->setGlobalGroupMembership(tmpGroupMembership);
637     #else
638     info->setGlobalGroupMembership(globalGroupMembership);
639     #endif
640 chrisfen 2101
641 gezelter 2087 //fill molMembership
642     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
643    
644     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
645 chrisfen 2101
646 gezelter 2087 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
647     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
648     }
649     }
650 chrisfen 2101
651 gezelter 2087 #ifdef IS_MPI
652     std::vector<int> tmpMolMembership(nGlobalAtoms, 0);
653 chrisfen 2101
654 gezelter 2087 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
655     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
656    
657     info->setGlobalMolMembership(tmpMolMembership);
658     #else
659     info->setGlobalMolMembership(globalMolMembership);
660     #endif
661 chrisfen 2101
662 gezelter 2087 }
663 chrisfen 2101
664 gezelter 2087 void SimCreator::loadCoordinates(SimInfo* info) {
665     Globals* simParams;
666     simParams = info->getSimParams();
667    
668     if (!simParams->haveInitialConfig()) {
669     sprintf(painCave.errMsg,
670     "Cannot intialize a simulation without an initial configuration file.\n");
671     painCave.isFatal = 1;;
672     simError();
673     }
674 chrisfen 2101
675 gezelter 2087 DumpReader reader(info, simParams->getInitialConfig());
676     int nframes = reader.getNFrames();
677 chrisfen 2101
678 gezelter 2087 if (nframes > 0) {
679     reader.readFrame(nframes - 1);
680     } else {
681     //invalid initial coordinate file
682 chrisfen 2101 sprintf(painCave.errMsg,
683     "Initial configuration file %s should at least contain one frame\n",
684 tim 2364 simParams->getInitialConfig().c_str());
685 gezelter 2087 painCave.isFatal = 1;
686     simError();
687     }
688 chrisfen 2101
689 gezelter 2087 //copy the current snapshot to previous snapshot
690     info->getSnapshotManager()->advance();
691     }
692 chrisfen 2101
693 gezelter 2087 } //end namespace oopse
694    
695