ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/brains/SimCreator.cpp
Revision: 3470
Committed: Wed Oct 22 20:01:49 2008 UTC (16 years, 6 months ago) by gezelter
File size: 26986 byte(s)
Log Message:
General bug-fixes and other changes to make particle pots work with
the Helfand Energy correlation function

File Contents

# Content
1 /*
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 #include <exception>
50 #include <iostream>
51 #include <sstream>
52 #include <string>
53
54 #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 #include "mdParser/MDLexer.hpp"
63 #include "mdParser/MDParser.hpp"
64 #include "mdParser/MDTreeParser.hpp"
65 #include "mdParser/SimplePreprocessor.hpp"
66 #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
73 #include "antlr/MismatchedCharException.hpp"
74 #include "antlr/MismatchedTokenException.hpp"
75 #include "antlr/NoViableAltForCharException.hpp"
76 #include "antlr/NoViableAltException.hpp"
77
78 #ifdef IS_MPI
79 #include "math/ParallelRandNumGen.hpp"
80 #endif
81
82 namespace oopse {
83
84 Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int startOfMetaDataBlock ){
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(rawMetaDataStream, filename, startOfMetaDataBlock, 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(static_cast<void*>(const_cast<char*>(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(filename);
126 lexer.initDeferredLineCount();
127
128 // Create a parser that reads from the scanner
129 MDParser parser(lexer);
130 parser.setFilename(filename);
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
138 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
151 catch(antlr::MismatchedCharException& e) {
152 sprintf(painCave.errMsg,
153 "parser exception: %s %s:%d:%d\n",
154 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
155 painCave.isFatal = 1;
156 simError();
157 }
158 catch(antlr::MismatchedTokenException &e) {
159 sprintf(painCave.errMsg,
160 "parser exception: %s %s:%d:%d\n",
161 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
162 painCave.isFatal = 1;
163 simError();
164 }
165 catch(antlr::NoViableAltForCharException &e) {
166 sprintf(painCave.errMsg,
167 "parser exception: %s %s:%d:%d\n",
168 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
169 painCave.isFatal = 1;
170 simError();
171 }
172 catch(antlr::NoViableAltException &e) {
173 sprintf(painCave.errMsg,
174 "parser exception: %s %s:%d:%d\n",
175 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
176 painCave.isFatal = 1;
177 simError();
178 }
179
180 catch(antlr::TokenStreamRecognitionException& e) {
181 sprintf(painCave.errMsg,
182 "parser exception: %s %s:%d:%d\n",
183 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
184 painCave.isFatal = 1;
185 simError();
186 }
187
188 catch(antlr::TokenStreamIOException& e) {
189 sprintf(painCave.errMsg,
190 "parser exception: %s\n",
191 e.getMessage().c_str());
192 painCave.isFatal = 1;
193 simError();
194 }
195
196 catch(antlr::TokenStreamException& e) {
197 sprintf(painCave.errMsg,
198 "parser exception: %s\n",
199 e.getMessage().c_str());
200 painCave.isFatal = 1;
201 simError();
202 }
203 catch (antlr::RecognitionException& e) {
204 sprintf(painCave.errMsg,
205 "parser exception: %s %s:%d:%d\n",
206 e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
207 painCave.isFatal = 1;
208 simError();
209 }
210 catch (antlr::CharStreamException& e) {
211 sprintf(painCave.errMsg,
212 "parser exception: %s\n",
213 e.getMessage().c_str());
214 painCave.isFatal = 1;
215 simError();
216 }
217 catch (OOPSEException& e) {
218 sprintf(painCave.errMsg,
219 "%s\n",
220 e.getMessage().c_str());
221 painCave.isFatal = 1;
222 simError();
223 }
224 catch (std::exception& e) {
225 sprintf(painCave.errMsg,
226 "parser exception: %s\n",
227 e.what());
228 painCave.isFatal = 1;
229 simError();
230 }
231
232 return simParams;
233 }
234
235 SimInfo* SimCreator::createSim(const std::string & mdFileName,
236 bool loadInitCoords) {
237
238 const int bufferSize = 65535;
239 char buffer[bufferSize];
240 int lineNo = 0;
241 std::string mdRawData;
242 int metaDataBlockStart = -1;
243 int metaDataBlockEnd = -1;
244 int i;
245 int mdOffset;
246
247 #ifdef IS_MPI
248 const int masterNode = 0;
249 if (worldRank == masterNode) {
250 #endif
251
252 std::ifstream mdFile_(mdFileName.c_str());
253
254 if (mdFile_.fail()) {
255 sprintf(painCave.errMsg,
256 "SimCreator: Cannot open file: %s\n",
257 mdFileName.c_str());
258 painCave.isFatal = 1;
259 simError();
260 }
261
262 mdFile_.getline(buffer, bufferSize);
263 ++lineNo;
264 std::string line = trimLeftCopy(buffer);
265 i = CaseInsensitiveFind(line, "<OOPSE");
266 if (static_cast<size_t>(i) == string::npos) {
267 sprintf(painCave.errMsg,
268 "SimCreator: File: %s is not an OOPSE file!\n",
269 mdFileName.c_str());
270 painCave.isFatal = 1;
271 simError();
272 }
273
274 //scan through the input stream and find MetaData tag
275 while(mdFile_.getline(buffer, bufferSize)) {
276 ++lineNo;
277
278 std::string line = trimLeftCopy(buffer);
279 if (metaDataBlockStart == -1) {
280 i = CaseInsensitiveFind(line, "<MetaData>");
281 if (i != string::npos) {
282 metaDataBlockStart = lineNo;
283 mdOffset = mdFile_.tellg();
284 }
285 } else {
286 i = CaseInsensitiveFind(line, "</MetaData>");
287 if (i != string::npos) {
288 metaDataBlockEnd = lineNo;
289 }
290 }
291 }
292
293 if (metaDataBlockStart == -1) {
294 sprintf(painCave.errMsg,
295 "SimCreator: File: %s did not contain a <MetaData> tag!\n",
296 mdFileName.c_str());
297 painCave.isFatal = 1;
298 simError();
299 }
300 if (metaDataBlockEnd == -1) {
301 sprintf(painCave.errMsg,
302 "SimCreator: File: %s did not contain a closed MetaData block!\n",
303 mdFileName.c_str());
304 painCave.isFatal = 1;
305 simError();
306 }
307
308 mdFile_.clear();
309 mdFile_.seekg(0);
310 mdFile_.seekg(mdOffset);
311
312 mdRawData.clear();
313
314 for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
315 mdFile_.getline(buffer, bufferSize);
316 mdRawData += buffer;
317 mdRawData += "\n";
318 }
319
320 mdFile_.close();
321
322 #ifdef IS_MPI
323 }
324 #endif
325
326 std::stringstream rawMetaDataStream(mdRawData);
327
328 //parse meta-data file
329 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, metaDataBlockStart+1);
330
331 //create the force field
332 ForceField * ff = ForceFieldFactory::getInstance()->createForceField(simParams->getForceField());
333
334 if (ff == NULL) {
335 sprintf(painCave.errMsg,
336 "ForceField Factory can not create %s force field\n",
337 simParams->getForceField().c_str());
338 painCave.isFatal = 1;
339 simError();
340 }
341
342 if (simParams->haveForceFieldFileName()) {
343 ff->setForceFieldFileName(simParams->getForceFieldFileName());
344 }
345
346 std::string forcefieldFileName;
347 forcefieldFileName = ff->getForceFieldFileName();
348
349 if (simParams->haveForceFieldVariant()) {
350 //If the force field has variant, the variant force field name will be
351 //Base.variant.frc. For exampel EAM.u6.frc
352
353 std::string variant = simParams->getForceFieldVariant();
354
355 std::string::size_type pos = forcefieldFileName.rfind(".frc");
356 variant = "." + variant;
357 if (pos != std::string::npos) {
358 forcefieldFileName.insert(pos, variant);
359 } else {
360 //If the default force field file name does not containt .frc suffix, just append the .variant
361 forcefieldFileName.append(variant);
362 }
363 }
364
365 ff->parse(forcefieldFileName);
366 ff->setFortranForceOptions();
367 //create SimInfo
368 SimInfo * info = new SimInfo(ff, simParams);
369
370 info->setRawMetaData(mdRawData);
371
372 //gather parameters (SimCreator only retrieves part of the
373 //parameters)
374 gatherParameters(info, mdFileName);
375
376 //divide the molecules and determine the global index of molecules
377 #ifdef IS_MPI
378 divideMolecules(info);
379 #endif
380
381 //create the molecules
382 createMolecules(info);
383
384
385 //allocate memory for DataStorage(circular reference, need to
386 //break it)
387 info->setSnapshotManager(new SimSnapshotManager(info));
388
389 //set the global index of atoms, rigidbodies and cutoffgroups
390 //(only need to be set once, the global index will never change
391 //again). Local indices of atoms and rigidbodies are already set
392 //by MoleculeCreator class which actually delegates the
393 //responsibility to LocalIndexManager.
394 setGlobalIndex(info);
395
396 //Although addInteractionPairs is called inside SimInfo's addMolecule
397 //method, at that point atoms don't have the global index yet
398 //(their global index are all initialized to -1). Therefore we
399 //have to call addInteractionPairs explicitly here. A way to work
400 //around is that we can determine the beginning global indices of
401 //atoms before they get created.
402 SimInfo::MoleculeIterator mi;
403 Molecule* mol;
404 for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
405 info->addInteractionPairs(mol);
406 }
407
408 if (loadInitCoords)
409 loadCoordinates(info, mdFileName);
410
411 return info;
412 }
413
414 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
415
416 //figure out the output file names
417 std::string prefix;
418
419 #ifdef IS_MPI
420
421 if (worldRank == 0) {
422 #endif // is_mpi
423 Globals * simParams = info->getSimParams();
424 if (simParams->haveFinalConfig()) {
425 prefix = getPrefix(simParams->getFinalConfig());
426 } else {
427 prefix = getPrefix(mdfile);
428 }
429
430 info->setFinalConfigFileName(prefix + ".eor");
431 info->setDumpFileName(prefix + ".dump");
432 info->setStatFileName(prefix + ".stat");
433 info->setRestFileName(prefix + ".zang");
434
435 #ifdef IS_MPI
436
437 }
438
439 #endif
440
441 }
442
443 #ifdef IS_MPI
444 void SimCreator::divideMolecules(SimInfo *info) {
445 RealType numerator;
446 RealType denominator;
447 RealType precast;
448 RealType x;
449 RealType y;
450 RealType a;
451 int old_atoms;
452 int add_atoms;
453 int new_atoms;
454 int nTarget;
455 int done;
456 int i;
457 int j;
458 int loops;
459 int which_proc;
460 int nProcessors;
461 std::vector<int> atomsPerProc;
462 int nGlobalMols = info->getNGlobalMolecules();
463 std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
464
465 MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
466
467 if (nProcessors > nGlobalMols) {
468 sprintf(painCave.errMsg,
469 "nProcessors (%d) > nMol (%d)\n"
470 "\tThe number of processors is larger than\n"
471 "\tthe number of molecules. This will not result in a \n"
472 "\tusable division of atoms for force decomposition.\n"
473 "\tEither try a smaller number of processors, or run the\n"
474 "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
475
476 painCave.isFatal = 1;
477 simError();
478 }
479
480 int seedValue;
481 Globals * simParams = info->getSimParams();
482 SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
483 if (simParams->haveSeed()) {
484 seedValue = simParams->getSeed();
485 myRandom = new SeqRandNumGen(seedValue);
486 }else {
487 myRandom = new SeqRandNumGen();
488 }
489
490
491 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
492
493 //initialize atomsPerProc
494 atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
495
496 if (worldRank == 0) {
497 numerator = info->getNGlobalAtoms();
498 denominator = nProcessors;
499 precast = numerator / denominator;
500 nTarget = (int)(precast + 0.5);
501
502 for(i = 0; i < nGlobalMols; i++) {
503 done = 0;
504 loops = 0;
505
506 while (!done) {
507 loops++;
508
509 // Pick a processor at random
510
511 which_proc = (int) (myRandom->rand() * nProcessors);
512
513 //get the molecule stamp first
514 int stampId = info->getMoleculeStampId(i);
515 MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
516
517 // How many atoms does this processor have so far?
518 old_atoms = atomsPerProc[which_proc];
519 add_atoms = moleculeStamp->getNAtoms();
520 new_atoms = old_atoms + add_atoms;
521
522 // If we've been through this loop too many times, we need
523 // to just give up and assign the molecule to this processor
524 // and be done with it.
525
526 if (loops > 100) {
527 sprintf(painCave.errMsg,
528 "I've tried 100 times to assign molecule %d to a "
529 " processor, but can't find a good spot.\n"
530 "I'm assigning it at random to processor %d.\n",
531 i, which_proc);
532
533 painCave.isFatal = 0;
534 simError();
535
536 molToProcMap[i] = which_proc;
537 atomsPerProc[which_proc] += add_atoms;
538
539 done = 1;
540 continue;
541 }
542
543 // If we can add this molecule to this processor without sending
544 // it above nTarget, then go ahead and do it:
545
546 if (new_atoms <= nTarget) {
547 molToProcMap[i] = which_proc;
548 atomsPerProc[which_proc] += add_atoms;
549
550 done = 1;
551 continue;
552 }
553
554 // The only situation left is when new_atoms > nTarget. We
555 // want to accept this with some probability that dies off the
556 // farther we are from nTarget
557
558 // roughly: x = new_atoms - nTarget
559 // Pacc(x) = exp(- a * x)
560 // where a = penalty / (average atoms per molecule)
561
562 x = (RealType)(new_atoms - nTarget);
563 y = myRandom->rand();
564
565 if (y < exp(- a * x)) {
566 molToProcMap[i] = which_proc;
567 atomsPerProc[which_proc] += add_atoms;
568
569 done = 1;
570 continue;
571 } else {
572 continue;
573 }
574 }
575 }
576
577 delete myRandom;
578
579 // Spray out this nonsense to all other processors:
580
581 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
582 } else {
583
584 // Listen to your marching orders from processor 0:
585
586 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
587 }
588
589 info->setMolToProcMap(molToProcMap);
590 sprintf(checkPointMsg,
591 "Successfully divided the molecules among the processors.\n");
592 errorCheckPoint();
593 }
594
595 #endif
596
597 void SimCreator::createMolecules(SimInfo *info) {
598 MoleculeCreator molCreator;
599 int stampId;
600
601 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
602
603 #ifdef IS_MPI
604
605 if (info->getMolToProc(i) == worldRank) {
606 #endif
607
608 stampId = info->getMoleculeStampId(i);
609 Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
610 stampId, i, info->getLocalIndexManager());
611
612 info->addMolecule(mol);
613
614 #ifdef IS_MPI
615
616 }
617
618 #endif
619
620 } //end for(int i=0)
621 }
622
623 void SimCreator::setGlobalIndex(SimInfo *info) {
624 SimInfo::MoleculeIterator mi;
625 Molecule::AtomIterator ai;
626 Molecule::RigidBodyIterator ri;
627 Molecule::CutoffGroupIterator ci;
628 Molecule::IntegrableObjectIterator ioi;
629 Molecule * mol;
630 Atom * atom;
631 RigidBody * rb;
632 CutoffGroup * cg;
633 int beginAtomIndex;
634 int beginRigidBodyIndex;
635 int beginCutoffGroupIndex;
636 int nGlobalAtoms = info->getNGlobalAtoms();
637
638 /**@todo fixme */
639 #ifndef IS_MPI
640
641 beginAtomIndex = 0;
642 beginRigidBodyIndex = 0;
643 beginCutoffGroupIndex = 0;
644
645 #else
646
647 int nproc;
648 int myNode;
649
650 myNode = worldRank;
651 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
652
653 std::vector < int > tmpAtomsInProc(nproc, 0);
654 std::vector < int > tmpRigidBodiesInProc(nproc, 0);
655 std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
656 std::vector < int > NumAtomsInProc(nproc, 0);
657 std::vector < int > NumRigidBodiesInProc(nproc, 0);
658 std::vector < int > NumCutoffGroupsInProc(nproc, 0);
659
660 tmpAtomsInProc[myNode] = info->getNAtoms();
661 tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
662 tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
663
664 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
665 MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
666 MPI_SUM, MPI_COMM_WORLD);
667 MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
668 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
669 MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
670 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
671
672 beginAtomIndex = 0;
673 beginRigidBodyIndex = 0;
674 beginCutoffGroupIndex = 0;
675
676 for(int i = 0; i < myNode; i++) {
677 beginAtomIndex += NumAtomsInProc[i];
678 beginRigidBodyIndex += NumRigidBodiesInProc[i];
679 beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
680 }
681
682 #endif
683
684 //rigidbody's index begins right after atom's
685 beginRigidBodyIndex += info->getNGlobalAtoms();
686
687 for(mol = info->beginMolecule(mi); mol != NULL;
688 mol = info->nextMolecule(mi)) {
689
690 //local index(index in DataStorge) of atom is important
691 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
692 atom->setGlobalIndex(beginAtomIndex++);
693 }
694
695 for(rb = mol->beginRigidBody(ri); rb != NULL;
696 rb = mol->nextRigidBody(ri)) {
697 rb->setGlobalIndex(beginRigidBodyIndex++);
698 }
699
700 //local index of cutoff group is trivial, it only depends on the order of travesing
701 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
702 cg = mol->nextCutoffGroup(ci)) {
703 cg->setGlobalIndex(beginCutoffGroupIndex++);
704 }
705 }
706
707 //fill globalGroupMembership
708 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
709 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
710 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
711
712 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
713 globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
714 }
715
716 }
717 }
718
719 #ifdef IS_MPI
720 // Since the globalGroupMembership has been zero filled and we've only
721 // poked values into the atoms we know, we can do an Allreduce
722 // to get the full globalGroupMembership array (We think).
723 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
724 // docs said we could.
725 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
726 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
727 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
728 info->setGlobalGroupMembership(tmpGroupMembership);
729 #else
730 info->setGlobalGroupMembership(globalGroupMembership);
731 #endif
732
733 //fill molMembership
734 std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
735
736 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
737 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
738 globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
739 }
740 }
741
742 #ifdef IS_MPI
743 std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
744
745 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
746 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
747
748 info->setGlobalMolMembership(tmpMolMembership);
749 #else
750 info->setGlobalMolMembership(globalMolMembership);
751 #endif
752
753 // nIOPerMol holds the number of integrable objects per molecule
754 // here the molecules are listed by their global indices.
755
756 std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
757 for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
758 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
759 }
760
761 #ifdef IS_MPI
762 std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
763 MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
764 info->getNGlobalMolecules(), MPI_INT, MPI_SUM, MPI_COMM_WORLD);
765 #else
766 std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
767 #endif
768
769 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
770
771 int startingIndex = 0;
772 for (int i = 0; i < info->getNGlobalMolecules(); i++) {
773 startingIOIndexForMol[i] = startingIndex;
774 startingIndex += numIntegrableObjectsPerMol[i];
775 }
776
777 std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
778 for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
779 int myGlobalIndex = mol->getGlobalIndex();
780 int globalIO = startingIOIndexForMol[myGlobalIndex];
781 for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL;
782 integrableObject = mol->nextIntegrableObject(ioi)) {
783 integrableObject->setGlobalIntegrableObjectIndex(globalIO);
784 IOIndexToIntegrableObject[globalIO] = integrableObject;
785 globalIO++;
786 }
787 }
788
789 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
790
791 }
792
793 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
794 Globals* simParams;
795 simParams = info->getSimParams();
796
797
798 DumpReader reader(info, mdFileName);
799 int nframes = reader.getNFrames();
800
801 if (nframes > 0) {
802 reader.readFrame(nframes - 1);
803 } else {
804 //invalid initial coordinate file
805 sprintf(painCave.errMsg,
806 "Initial configuration file %s should at least contain one frame\n",
807 mdFileName.c_str());
808 painCave.isFatal = 1;
809 simError();
810 }
811
812 //copy the current snapshot to previous snapshot
813 info->getSnapshotManager()->advance();
814 }
815
816 } //end namespace oopse
817
818