ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2101
Committed: Thu Mar 10 15:10:24 2005 UTC (19 years, 4 months ago) by chrisfen
File size: 20267 byte(s)
Log Message:
First commit of the new restraints code

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