ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 1927
Committed: Wed Jan 12 17:14:35 2005 UTC (21 years, 4 months ago) by tim
File size: 20659 byte(s)
Log Message:
change license

File Contents

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