ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2187
Committed: Wed Apr 13 18:41:17 2005 UTC (19 years, 2 months ago) by tim
File size: 20004 byte(s)
Log Message:
more memory leak are fixed

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