ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2448
Committed: Wed Nov 16 23:10:02 2005 UTC (18 years, 7 months ago) by tim
File size: 19974 byte(s)
Log Message:
OptionSectionParser get compiled

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    
50     #include "brains/MoleculeCreator.hpp"
51     #include "brains/SimCreator.hpp"
52     #include "brains/SimSnapshotManager.hpp"
53     #include "io/DumpReader.hpp"
54     #include "io/parse_me.h"
55     #include "UseTheForce/ForceFieldFactory.hpp"
56     #include "utils/simError.h"
57     #include "utils/StringUtils.hpp"
58     #include "math/SeqRandNumGen.hpp"
59     #ifdef IS_MPI
60     #include "io/mpiBASS.h"
61     #include "math/ParallelRandNumGen.hpp"
62     #endif
63    
64     namespace oopse {
65 chrisfen 2101
66     void SimCreator::parseFile(const std::string mdFileName, MakeStamps* stamps,
67     Globals* simParams){
68    
69 gezelter 2087 #ifdef IS_MPI
70 chrisfen 2101
71 gezelter 2087 if (worldRank == 0) {
72     #endif // is_mpi
73 chrisfen 2101
74 gezelter 2087 set_interface_stamps(stamps, simParams);
75 chrisfen 2101
76 gezelter 2087 #ifdef IS_MPI
77 chrisfen 2101
78 gezelter 2087 mpiEventInit();
79 chrisfen 2101
80 gezelter 2087 #endif
81 chrisfen 2101
82 gezelter 2087 yacc_BASS(mdFileName.c_str());
83 chrisfen 2101
84 gezelter 2087 #ifdef IS_MPI
85 chrisfen 2101
86 gezelter 2087 throwMPIEvent(NULL);
87     } else {
88     set_interface_stamps(stamps, simParams);
89     mpiEventInit();
90     MPIcheckPoint();
91     mpiEventLoop();
92     }
93 chrisfen 2101
94 gezelter 2087 #endif
95 chrisfen 2101
96 gezelter 2087 }
97 chrisfen 2101
98 chrisfen 2211 SimInfo* SimCreator::createSim(const std::string & mdFileName,
99     bool loadInitCoords) {
100 gezelter 2087
101     MakeStamps * stamps = new MakeStamps();
102 chrisfen 2101
103 gezelter 2087 Globals * simParams = new Globals();
104 chrisfen 2101
105 gezelter 2087 //parse meta-data file
106     parseFile(mdFileName, stamps, simParams);
107 chrisfen 2101
108 gezelter 2087 //create the force field
109 chrisfen 2211 ForceField * ff = ForceFieldFactory::getInstance()
110     ->createForceField(simParams->getForceField());
111 gezelter 2087
112     if (ff == NULL) {
113 chrisfen 2211 sprintf(painCave.errMsg,
114     "ForceField Factory can not create %s force field\n",
115 tim 2364 simParams->getForceField().c_str());
116 gezelter 2087 painCave.isFatal = 1;
117     simError();
118     }
119 chrisfen 2101
120 gezelter 2087 if (simParams->haveForceFieldFileName()) {
121     ff->setForceFieldFileName(simParams->getForceFieldFileName());
122     }
123    
124     std::string forcefieldFileName;
125     forcefieldFileName = ff->getForceFieldFileName();
126 chrisfen 2101
127 gezelter 2087 if (simParams->haveForceFieldVariant()) {
128     //If the force field has variant, the variant force field name will be
129     //Base.variant.frc. For exampel EAM.u6.frc
130 chrisfen 2101
131 gezelter 2087 std::string variant = simParams->getForceFieldVariant();
132 chrisfen 2101
133 gezelter 2087 std::string::size_type pos = forcefieldFileName.rfind(".frc");
134     variant = "." + variant;
135     if (pos != std::string::npos) {
136     forcefieldFileName.insert(pos, variant);
137     } else {
138     //If the default force field file name does not containt .frc suffix, just append the .variant
139     forcefieldFileName.append(variant);
140     }
141     }
142    
143     ff->parse(forcefieldFileName);
144    
145     //extract the molecule stamps
146     std::vector < std::pair<MoleculeStamp *, int> > moleculeStampPairs;
147     compList(stamps, simParams, moleculeStampPairs);
148 chrisfen 2101
149 gezelter 2087 //create SimInfo
150 tim 2187 SimInfo * info = new SimInfo(stamps, moleculeStampPairs, ff, simParams);
151    
152 gezelter 2087 //gather parameters (SimCreator only retrieves part of the parameters)
153     gatherParameters(info, mdFileName);
154 chrisfen 2101
155 gezelter 2087 //divide the molecules and determine the global index of molecules
156     #ifdef IS_MPI
157     divideMolecules(info);
158     #endif
159 chrisfen 2101
160 gezelter 2087 //create the molecules
161     createMolecules(info);
162 chrisfen 2101
163    
164 gezelter 2087 //allocate memory for DataStorage(circular reference, need to break it)
165     info->setSnapshotManager(new SimSnapshotManager(info));
166    
167     //set the global index of atoms, rigidbodies and cutoffgroups (only need to be set once, the
168     //global index will never change again). Local indices of atoms and rigidbodies are already set by
169     //MoleculeCreator class which actually delegates the responsibility to LocalIndexManager.
170     setGlobalIndex(info);
171 chrisfen 2101
172 gezelter 2087 //Alought addExculdePairs is called inside SimInfo's addMolecule method, at that point
173     //atoms don't have the global index yet (their global index are all initialized to -1).
174     //Therefore we have to call addExcludePairs explicitly here. A way to work around is that
175     //we can determine the beginning global indices of atoms before they get created.
176     SimInfo::MoleculeIterator mi;
177     Molecule* mol;
178     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
179     info->addExcludePairs(mol);
180     }
181    
182     if (loadInitCoords)
183     loadCoordinates(info);
184    
185     return info;
186     }
187 chrisfen 2101
188 gezelter 2087 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
189 chrisfen 2101
190 tim 2448 //figure out the output file names
191 gezelter 2087 std::string prefix;
192 chrisfen 2101
193 gezelter 2087 #ifdef IS_MPI
194 chrisfen 2101
195 gezelter 2087 if (worldRank == 0) {
196     #endif // is_mpi
197     Globals * simParams = info->getSimParams();
198     if (simParams->haveFinalConfig()) {
199     prefix = getPrefix(simParams->getFinalConfig());
200     } else {
201     prefix = getPrefix(mdfile);
202     }
203 chrisfen 2101
204 gezelter 2087 info->setFinalConfigFileName(prefix + ".eor");
205     info->setDumpFileName(prefix + ".dump");
206     info->setStatFileName(prefix + ".stat");
207 chrisfen 2101 info->setRestFileName(prefix + ".zang");
208    
209 gezelter 2087 #ifdef IS_MPI
210 chrisfen 2101
211 gezelter 2087 }
212 chrisfen 2101
213 gezelter 2087 #endif
214 chrisfen 2101
215 gezelter 2087 }
216 chrisfen 2101
217 gezelter 2087 #ifdef IS_MPI
218     void SimCreator::divideMolecules(SimInfo *info) {
219     double numerator;
220     double denominator;
221     double precast;
222     double x;
223     double y;
224     double a;
225     int old_atoms;
226     int add_atoms;
227     int new_atoms;
228     int nTarget;
229     int done;
230     int i;
231     int j;
232     int loops;
233     int which_proc;
234     int nProcessors;
235     std::vector<int> atomsPerProc;
236     int nGlobalMols = info->getNGlobalMolecules();
237     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
238    
239     MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
240 chrisfen 2101
241 gezelter 2087 if (nProcessors > nGlobalMols) {
242     sprintf(painCave.errMsg,
243     "nProcessors (%d) > nMol (%d)\n"
244     "\tThe number of processors is larger than\n"
245     "\tthe number of molecules. This will not result in a \n"
246     "\tusable division of atoms for force decomposition.\n"
247     "\tEither try a smaller number of processors, or run the\n"
248     "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
249 chrisfen 2101
250 gezelter 2087 painCave.isFatal = 1;
251     simError();
252     }
253 chrisfen 2101
254 gezelter 2087 int seedValue;
255     Globals * simParams = info->getSimParams();
256     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
257     if (simParams->haveSeed()) {
258     seedValue = simParams->getSeed();
259     myRandom = new SeqRandNumGen(seedValue);
260     }else {
261     myRandom = new SeqRandNumGen();
262     }
263 chrisfen 2101
264    
265 gezelter 2087 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
266 chrisfen 2101
267 gezelter 2087 //initialize atomsPerProc
268     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
269 chrisfen 2101
270 gezelter 2087 if (worldRank == 0) {
271     numerator = info->getNGlobalAtoms();
272     denominator = nProcessors;
273     precast = numerator / denominator;
274     nTarget = (int)(precast + 0.5);
275 chrisfen 2101
276 gezelter 2087 for(i = 0; i < nGlobalMols; i++) {
277     done = 0;
278     loops = 0;
279 chrisfen 2101
280 gezelter 2087 while (!done) {
281     loops++;
282 chrisfen 2101
283 gezelter 2087 // Pick a processor at random
284 chrisfen 2101
285 gezelter 2087 which_proc = (int) (myRandom->rand() * nProcessors);
286 chrisfen 2101
287 gezelter 2087 //get the molecule stamp first
288     int stampId = info->getMoleculeStampId(i);
289     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
290 chrisfen 2101
291 gezelter 2087 // How many atoms does this processor have so far?
292     old_atoms = atomsPerProc[which_proc];
293     add_atoms = moleculeStamp->getNAtoms();
294     new_atoms = old_atoms + add_atoms;
295 chrisfen 2101
296 gezelter 2087 // If we've been through this loop too many times, we need
297     // to just give up and assign the molecule to this processor
298     // and be done with it.
299 chrisfen 2101
300 gezelter 2087 if (loops > 100) {
301     sprintf(painCave.errMsg,
302     "I've tried 100 times to assign molecule %d to a "
303     " processor, but can't find a good spot.\n"
304     "I'm assigning it at random to processor %d.\n",
305     i, which_proc);
306 chrisfen 2101
307 gezelter 2087 painCave.isFatal = 0;
308     simError();
309 chrisfen 2101
310 gezelter 2087 molToProcMap[i] = which_proc;
311     atomsPerProc[which_proc] += add_atoms;
312 chrisfen 2101
313 gezelter 2087 done = 1;
314     continue;
315     }
316 chrisfen 2101
317 gezelter 2087 // If we can add this molecule to this processor without sending
318     // it above nTarget, then go ahead and do it:
319 chrisfen 2101
320 gezelter 2087 if (new_atoms <= nTarget) {
321     molToProcMap[i] = which_proc;
322     atomsPerProc[which_proc] += add_atoms;
323 chrisfen 2101
324 gezelter 2087 done = 1;
325     continue;
326     }
327 chrisfen 2101
328 gezelter 2087 // The only situation left is when new_atoms > nTarget. We
329     // want to accept this with some probability that dies off the
330     // farther we are from nTarget
331 chrisfen 2101
332 gezelter 2087 // roughly: x = new_atoms - nTarget
333     // Pacc(x) = exp(- a * x)
334     // where a = penalty / (average atoms per molecule)
335 chrisfen 2101
336 gezelter 2087 x = (double)(new_atoms - nTarget);
337     y = myRandom->rand();
338 chrisfen 2101
339 gezelter 2087 if (y < exp(- a * x)) {
340     molToProcMap[i] = which_proc;
341     atomsPerProc[which_proc] += add_atoms;
342 chrisfen 2101
343 gezelter 2087 done = 1;
344     continue;
345     } else {
346     continue;
347     }
348     }
349     }
350 chrisfen 2101
351 gezelter 2087 delete myRandom;
352 chrisfen 2101
353 gezelter 2087 // Spray out this nonsense to all other processors:
354 chrisfen 2101
355 gezelter 2087 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
356     } else {
357 chrisfen 2101
358 gezelter 2087 // Listen to your marching orders from processor 0:
359 chrisfen 2101
360 gezelter 2087 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
361     }
362 chrisfen 2101
363 gezelter 2087 info->setMolToProcMap(molToProcMap);
364     sprintf(checkPointMsg,
365     "Successfully divided the molecules among the processors.\n");
366     MPIcheckPoint();
367     }
368 chrisfen 2101
369 gezelter 2087 #endif
370 chrisfen 2101
371 gezelter 2087 void SimCreator::createMolecules(SimInfo *info) {
372     MoleculeCreator molCreator;
373     int stampId;
374 chrisfen 2101
375 gezelter 2087 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
376 chrisfen 2101
377 gezelter 2087 #ifdef IS_MPI
378 chrisfen 2101
379 gezelter 2087 if (info->getMolToProc(i) == worldRank) {
380     #endif
381 chrisfen 2101
382 gezelter 2087 stampId = info->getMoleculeStampId(i);
383     Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
384     stampId, i, info->getLocalIndexManager());
385 chrisfen 2101
386 gezelter 2087 info->addMolecule(mol);
387 chrisfen 2101
388 gezelter 2087 #ifdef IS_MPI
389 chrisfen 2101
390 gezelter 2087 }
391 chrisfen 2101
392 gezelter 2087 #endif
393 chrisfen 2101
394 gezelter 2087 } //end for(int i=0)
395     }
396 chrisfen 2101
397 gezelter 2087 void SimCreator::compList(MakeStamps *stamps, Globals* simParams,
398     std::vector < std::pair<MoleculeStamp *, int> > &moleculeStampPairs) {
399     int i;
400     char * id;
401     MoleculeStamp * currentStamp;
402     Component** the_components = simParams->getComponents();
403     int n_components = simParams->getNComponents();
404 chrisfen 2101
405 gezelter 2087 if (!simParams->haveNMol()) {
406     // we don't have the total number of molecules, so we assume it is
407     // given in each component
408 chrisfen 2101
409 gezelter 2087 for(i = 0; i < n_components; i++) {
410     if (!the_components[i]->haveNMol()) {
411     // we have a problem
412     sprintf(painCave.errMsg,
413     "SimCreator Error. No global NMol or component NMol given.\n"
414     "\tCannot calculate the number of atoms.\n");
415 chrisfen 2101
416 gezelter 2087 painCave.isFatal = 1;
417     simError();
418     }
419 chrisfen 2101
420 gezelter 2087 id = the_components[i]->getType();
421 tim 2187
422     currentStamp = stamps->getMolStamp(id);
423     if (currentStamp == NULL) {
424 gezelter 2087 sprintf(painCave.errMsg,
425     "SimCreator error: Component \"%s\" was not found in the "
426     "list of declared molecules\n", id);
427 chrisfen 2101
428 gezelter 2087 painCave.isFatal = 1;
429     simError();
430     }
431 chrisfen 2101
432 gezelter 2087 moleculeStampPairs.push_back(
433     std::make_pair(currentStamp, the_components[i]->getNMol()));
434     } //end for (i = 0; i < n_components; i++)
435     } else {
436     sprintf(painCave.errMsg, "SimSetup error.\n"
437     "\tSorry, the ability to specify total"
438     " nMols and then give molfractions in the components\n"
439     "\tis not currently supported."
440     " Please give nMol in the components.\n");
441 chrisfen 2101
442 gezelter 2087 painCave.isFatal = 1;
443     simError();
444     }
445 chrisfen 2101
446 gezelter 2087 #ifdef IS_MPI
447 chrisfen 2101
448 gezelter 2087 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
449     MPIcheckPoint();
450 chrisfen 2101
451 gezelter 2087 #endif // is_mpi
452 chrisfen 2101
453 gezelter 2087 }
454 chrisfen 2101
455 gezelter 2087 void SimCreator::setGlobalIndex(SimInfo *info) {
456     SimInfo::MoleculeIterator mi;
457     Molecule::AtomIterator ai;
458     Molecule::RigidBodyIterator ri;
459     Molecule::CutoffGroupIterator ci;
460     Molecule * mol;
461     Atom * atom;
462     RigidBody * rb;
463     CutoffGroup * cg;
464     int beginAtomIndex;
465     int beginRigidBodyIndex;
466     int beginCutoffGroupIndex;
467     int nGlobalAtoms = info->getNGlobalAtoms();
468    
469     #ifndef IS_MPI
470 chrisfen 2101
471 gezelter 2087 beginAtomIndex = 0;
472     beginRigidBodyIndex = 0;
473     beginCutoffGroupIndex = 0;
474 chrisfen 2101
475 gezelter 2087 #else
476 chrisfen 2101
477 gezelter 2087 int nproc;
478     int myNode;
479 chrisfen 2101
480 gezelter 2087 myNode = worldRank;
481     MPI_Comm_size(MPI_COMM_WORLD, &nproc);
482 chrisfen 2101
483 gezelter 2087 std::vector < int > tmpAtomsInProc(nproc, 0);
484     std::vector < int > tmpRigidBodiesInProc(nproc, 0);
485     std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
486     std::vector < int > NumAtomsInProc(nproc, 0);
487     std::vector < int > NumRigidBodiesInProc(nproc, 0);
488     std::vector < int > NumCutoffGroupsInProc(nproc, 0);
489 chrisfen 2101
490 gezelter 2087 tmpAtomsInProc[myNode] = info->getNAtoms();
491     tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
492     tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
493 chrisfen 2101
494 gezelter 2087 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
495     MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
496     MPI_SUM, MPI_COMM_WORLD);
497     MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
498     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
499     MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
500     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
501 chrisfen 2101
502 gezelter 2087 beginAtomIndex = 0;
503     beginRigidBodyIndex = 0;
504     beginCutoffGroupIndex = 0;
505 chrisfen 2101
506 gezelter 2087 for(int i = 0; i < myNode; i++) {
507     beginAtomIndex += NumAtomsInProc[i];
508     beginRigidBodyIndex += NumRigidBodiesInProc[i];
509     beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
510     }
511 chrisfen 2101
512 gezelter 2087 #endif
513 chrisfen 2101
514 gezelter 2087 //rigidbody's index begins right after atom's
515     beginRigidBodyIndex += info->getNGlobalAtoms();
516 chrisfen 2101
517 gezelter 2087 for(mol = info->beginMolecule(mi); mol != NULL;
518     mol = info->nextMolecule(mi)) {
519 chrisfen 2101
520 gezelter 2087 //local index(index in DataStorge) of atom is important
521     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
522     atom->setGlobalIndex(beginAtomIndex++);
523     }
524 chrisfen 2101
525 gezelter 2087 for(rb = mol->beginRigidBody(ri); rb != NULL;
526     rb = mol->nextRigidBody(ri)) {
527     rb->setGlobalIndex(beginRigidBodyIndex++);
528     }
529 chrisfen 2101
530 gezelter 2087 //local index of cutoff group is trivial, it only depends on the order of travesing
531     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
532     cg = mol->nextCutoffGroup(ci)) {
533     cg->setGlobalIndex(beginCutoffGroupIndex++);
534     }
535     }
536 chrisfen 2101
537 gezelter 2087 //fill globalGroupMembership
538     std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
539     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
540     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
541 chrisfen 2101
542 gezelter 2087 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
543     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
544     }
545 chrisfen 2101
546 gezelter 2087 }
547     }
548 chrisfen 2101
549 gezelter 2087 #ifdef IS_MPI
550     // Since the globalGroupMembership has been zero filled and we've only
551     // poked values into the atoms we know, we can do an Allreduce
552     // to get the full globalGroupMembership array (We think).
553     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
554     // docs said we could.
555     std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
556     MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
557     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
558     info->setGlobalGroupMembership(tmpGroupMembership);
559     #else
560     info->setGlobalGroupMembership(globalGroupMembership);
561     #endif
562 chrisfen 2101
563 gezelter 2087 //fill molMembership
564     std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
565    
566     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
567 chrisfen 2101
568 gezelter 2087 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
569     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
570     }
571     }
572 chrisfen 2101
573 gezelter 2087 #ifdef IS_MPI
574     std::vector<int> tmpMolMembership(nGlobalAtoms, 0);
575 chrisfen 2101
576 gezelter 2087 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
577     MPI_INT, MPI_SUM, MPI_COMM_WORLD);
578    
579     info->setGlobalMolMembership(tmpMolMembership);
580     #else
581     info->setGlobalMolMembership(globalMolMembership);
582     #endif
583 chrisfen 2101
584 gezelter 2087 }
585 chrisfen 2101
586 gezelter 2087 void SimCreator::loadCoordinates(SimInfo* info) {
587     Globals* simParams;
588     simParams = info->getSimParams();
589    
590     if (!simParams->haveInitialConfig()) {
591     sprintf(painCave.errMsg,
592     "Cannot intialize a simulation without an initial configuration file.\n");
593     painCave.isFatal = 1;;
594     simError();
595     }
596 chrisfen 2101
597 gezelter 2087 DumpReader reader(info, simParams->getInitialConfig());
598     int nframes = reader.getNFrames();
599 chrisfen 2101
600 gezelter 2087 if (nframes > 0) {
601     reader.readFrame(nframes - 1);
602     } else {
603     //invalid initial coordinate file
604 chrisfen 2101 sprintf(painCave.errMsg,
605     "Initial configuration file %s should at least contain one frame\n",
606 tim 2364 simParams->getInitialConfig().c_str());
607 gezelter 2087 painCave.isFatal = 1;
608     simError();
609     }
610 chrisfen 2101
611 gezelter 2087 //copy the current snapshot to previous snapshot
612     info->getSnapshotManager()->advance();
613     }
614 chrisfen 2101
615 gezelter 2087 } //end namespace oopse
616    
617