ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 2469
Committed: Fri Dec 2 15:38:03 2005 UTC (18 years, 7 months ago) by tim
File size: 19946 byte(s)
Log Message:
End of the Link --> List
Return of the Oject-Oriented
replace yacc/lex parser with antlr parser

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