ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2544
Committed: Wed Jan 11 19:01:20 2006 UTC (18 years, 6 months ago) by tim
File size: 23292 byte(s)
Log Message:
instead of printing to std::cout, throwing an exception when error is found.

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