ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2522
Committed: Thu Dec 29 16:03:11 2005 UTC (18 years, 6 months ago) by chuckv
File size: 21940 byte(s)
Log Message:
Fixed issue with symbol naming and added call to set force options in SimCreator.

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