ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-2.0/src/brains/SimCreator.cpp
Revision: 1886
Committed: Wed Dec 15 16:13:35 2004 UTC (19 years, 7 months ago) by tim
File size: 19584 byte(s)
Log Message:
MPI in NVE is working

File Contents

# Content
1 /*
2 * Copyright (C) 2000-2004 Object Oriented Parallel Simulation Engine (OOPSE) project
3 *
4 * Contact: oopse@oopse.org
5 *
6 * This program is free software; you can redistribute it and/or
7 * modify it under the terms of the GNU Lesser General Public License
8 * as published by the Free Software Foundation; either version 2.1
9 * of the License, or (at your option) any later version.
10 * All we ask is that proper credit is given for our work, which includes
11 * - but is not limited to - adding the above copyright notice to the beginning
12 * of your source code files, and to any copyright notice that you may distribute
13 * with programs based on this work.
14 *
15 * This program is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 * GNU Lesser General Public License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with this program; if not, write to the Free Software
22 * Foundation, Inc., 59 Temple Place - Suite 330, Boston, MA 02111-1307, USA.
23 *
24 */
25
26 /**
27 * @file SimCreator.cpp
28 * @author tlin
29 * @date 11/03/2004
30 * @time 13:51am
31 * @version 1.0
32 */
33
34 #include <sprng.h>
35
36 #include "brains/MoleculeCreator.hpp"
37 #include "brains/SimCreator.hpp"
38 #include "brains/SimSnapshotManager.hpp"
39 #include "io/DumpReader.hpp"
40 #include "io/parse_me.h"
41 #include "UseTheForce/ForceFieldFactory.hpp"
42 #include "utils/simError.h"
43 #include "utils/StringUtils.hpp"
44 #ifdef IS_MPI
45 #include "io/mpiBASS.h"
46 #include "math/randomSPRNG.hpp"
47 #endif
48
49 namespace oopse {
50
51 void SimCreator::parseFile(const std::string mdFileName, MakeStamps* stamps, Globals* simParams){
52
53 #ifdef IS_MPI
54
55 if (worldRank == 0) {
56 #endif // is_mpi
57
58 simParams->initalize();
59 set_interface_stamps(stamps, simParams);
60
61 #ifdef IS_MPI
62
63 mpiEventInit();
64
65 #endif
66
67 yacc_BASS(mdFileName.c_str());
68
69 #ifdef IS_MPI
70
71 throwMPIEvent(NULL);
72 } else {
73 set_interface_stamps(stamps, simParams);
74 mpiEventInit();
75 MPIcheckPoint();
76 mpiEventLoop();
77 }
78
79 #endif
80
81 }
82
83 SimInfo* SimCreator::createSim(const std::string & mdFileName) {
84
85 MakeStamps * stamps = new MakeStamps();
86
87 Globals * simParams = new Globals();
88
89 //parse meta-data file
90 parseFile(mdFileName, stamps, simParams);
91
92 //create the force field
93 ForceField * ff = ForceFieldFactory::getInstance()->createForceField(
94 simParams->getForceField());
95
96 if (ff == NULL) {
97 sprintf(painCave.errMsg, "ForceField Factory can not create %s force field\n",
98 simParams->getForceField());
99 painCave.isFatal = 1;
100 simError();
101 }
102
103 std::string forcefieldFileName;
104 forcefieldFileName = ff->getForceFieldFileName();
105
106 if (simParams->haveForceFieldVariant()) {
107 //If the force field has variant, the variant force field name will be
108 //Base.variant.frc. For exampel EAM.u6.frc
109
110 std::string variant = simParams->getForceFieldVariant();
111
112 std::string::size_type pos = forcefieldFileName.rfind(".frc");
113 variant = "." + variant;
114 if (pos != std::string::npos) {
115 forcefieldFileName.insert(pos, variant);
116 } else {
117 //If the default force field file name does not containt .frc suffix, just append the .variant
118 forcefieldFileName.append(variant);
119 }
120 }
121
122 ff->parse(forcefieldFileName);
123
124 //extract the molecule stamps
125 std::vector < std::pair<MoleculeStamp *, int> > moleculeStampPairs;
126 compList(stamps, simParams, moleculeStampPairs);
127
128 //create SimInfo
129 SimInfo * info = new SimInfo(moleculeStampPairs, ff, simParams);
130
131 //gather parameters (SimCreator only retrieves part of the parameters)
132 gatherParameters(info, mdFileName);
133
134 //divide the molecules and determine the global index of molecules
135 #ifdef IS_MPI
136 divideMolecules(info);
137 #endif
138
139 //create the molecules
140 createMolecules(info);
141
142
143 //allocate memory for DataStorage(circular reference, need to break it)
144 info->setSnapshotManager(new SimSnapshotManager(info));
145
146 //set the global index of atoms, rigidbodies and cutoffgroups (only need to be set once, the
147 //global index will never change again). Local indices of atoms and rigidbodies are already set by
148 //MoleculeCreator class which actually delegates the responsibility to LocalIndexManager.
149 setGlobalIndex(info);
150
151 //Alought addExculdePairs is called inside SimInfo's addMolecule method, at that point
152 //atoms don't have the global index yet (their global index are all initialized to -1).
153 //Therefore we have to call addExcludePairs explicitly here. A way to work around is that
154 //we can determine the beginning global indices of atoms before they get created.
155 SimInfo::MoleculeIterator mi;
156 Molecule* mol;
157 for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
158 info->addExcludePairs(mol);
159 }
160
161
162 //load initial coordinates, some extra information are pushed into SimInfo's property map ( such as
163 //eta, chi for NPT integrator)
164 loadCoordinates(info);
165 return info;
166 }
167
168 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
169
170 //setup seed for random number generator
171 int seedValue;
172 Globals * simParams = info->getSimParams();
173
174 if (simParams->haveSeed()) {
175 seedValue = simParams->getSeed();
176
177 if (seedValue / 1000000000 == 0) {
178 sprintf(painCave.errMsg,
179 "Seed for sprng library should contain at least 9 digits\n"
180 "OOPSE will generate a seed for user\n");
181
182 painCave.isFatal = 0;
183 simError();
184
185 //using seed generated by system instead of invalid seed set by user
186
187 #ifndef IS_MPI
188
189 seedValue = make_sprng_seed();
190
191 #else
192
193 if (worldRank == 0) {
194 seedValue = make_sprng_seed();
195 }
196
197 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
198
199 #endif
200
201 } //end if (seedValue /1000000000 == 0)
202 } else {
203
204 #ifndef IS_MPI
205
206 seedValue = make_sprng_seed();
207
208 #else
209
210 if (worldRank == 0) {
211 seedValue = make_sprng_seed();
212 }
213
214 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
215
216 #endif
217
218 } //end of simParams->haveSeed()
219
220 info->setSeed(seedValue);
221
222
223 //figure out the ouput file names
224 std::string prefix;
225
226 #ifdef IS_MPI
227
228 if (worldRank == 0) {
229 #endif // is_mpi
230
231 if (simParams->haveFinalConfig()) {
232 prefix = getPrefix(simParams->getFinalConfig());
233 } else {
234 prefix = getPrefix(mdfile);
235 }
236
237 info->setFinalConfigFileName(prefix + ".eor");
238 info->setDumpFileName(prefix + ".dump");
239 info->setStatFileName(prefix + ".stat");
240
241 #ifdef IS_MPI
242
243 }
244
245 #endif
246
247 }
248
249 #ifdef IS_MPI
250 void SimCreator::divideMolecules(SimInfo *info) {
251 double numerator;
252 double denominator;
253 double precast;
254 double x;
255 double y;
256 double a;
257 int old_atoms;
258 int add_atoms;
259 int new_atoms;
260 int nTarget;
261 int done;
262 int i;
263 int j;
264 int loops;
265 int which_proc;
266 int nProcessors;
267 std::vector<int> atomsPerProc;
268 randomSPRNG myRandom(info->getSeed());
269 int nGlobalMols = info->getNGlobalMolecules();
270 std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
271
272 MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
273
274 if (nProcessors > nGlobalMols) {
275 sprintf(painCave.errMsg,
276 "nProcessors (%d) > nMol (%d)\n"
277 "\tThe number of processors is larger than\n"
278 "\tthe number of molecules. This will not result in a \n"
279 "\tusable division of atoms for force decomposition.\n"
280 "\tEither try a smaller number of processors, or run the\n"
281 "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
282
283 painCave.isFatal = 1;
284 simError();
285 }
286
287 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
288
289 //initialize atomsPerProc
290 atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
291
292 if (worldRank == 0) {
293 numerator = info->getNGlobalAtoms();
294 denominator = nProcessors;
295 precast = numerator / denominator;
296 nTarget = (int)(precast + 0.5);
297
298 for(i = 0; i < nGlobalMols; i++) {
299 done = 0;
300 loops = 0;
301
302 while (!done) {
303 loops++;
304
305 // Pick a processor at random
306
307 which_proc = (int) (myRandom.getRandom() * nProcessors);
308
309 //get the molecule stamp first
310 int stampId = info->getMoleculeStampId(i);
311 MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
312
313 // How many atoms does this processor have so far?
314 old_atoms = atomsPerProc[which_proc];
315 add_atoms = moleculeStamp->getNAtoms();
316 new_atoms = old_atoms + add_atoms;
317
318 // If we've been through this loop too many times, we need
319 // to just give up and assign the molecule to this processor
320 // and be done with it.
321
322 if (loops > 100) {
323 sprintf(painCave.errMsg,
324 "I've tried 100 times to assign molecule %d to a "
325 " processor, but can't find a good spot.\n"
326 "I'm assigning it at random to processor %d.\n",
327 i, which_proc);
328
329 painCave.isFatal = 0;
330 simError();
331
332 molToProcMap[i] = which_proc;
333 atomsPerProc[which_proc] += add_atoms;
334
335 done = 1;
336 continue;
337 }
338
339 // If we can add this molecule to this processor without sending
340 // it above nTarget, then go ahead and do it:
341
342 if (new_atoms <= nTarget) {
343 molToProcMap[i] = which_proc;
344 atomsPerProc[which_proc] += add_atoms;
345
346 done = 1;
347 continue;
348 }
349
350 // The only situation left is when new_atoms > nTarget. We
351 // want to accept this with some probability that dies off the
352 // farther we are from nTarget
353
354 // roughly: x = new_atoms - nTarget
355 // Pacc(x) = exp(- a * x)
356 // where a = penalty / (average atoms per molecule)
357
358 x = (double)(new_atoms - nTarget);
359 y = myRandom.getRandom();
360
361 if (y < exp(- a * x)) {
362 molToProcMap[i] = which_proc;
363 atomsPerProc[which_proc] += add_atoms;
364
365 done = 1;
366 continue;
367 } else {
368 continue;
369 }
370 }
371 }
372
373 // Spray out this nonsense to all other processors:
374
375 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
376 } else {
377
378 // Listen to your marching orders from processor 0:
379
380 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
381 }
382
383 info->setMolToProcMap(molToProcMap);
384 sprintf(checkPointMsg,
385 "Successfully divided the molecules among the processors.\n");
386 MPIcheckPoint();
387 }
388
389 #endif
390
391 void SimCreator::createMolecules(SimInfo *info) {
392 MoleculeCreator molCreator;
393 int stampId;
394
395 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
396
397 #ifdef IS_MPI
398
399 if (info->getMolToProc(i) == worldRank) {
400 #endif
401
402 stampId = info->getMoleculeStampId(i);
403 Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
404 stampId, i, info->getLocalIndexManager());
405
406 info->addMolecule(mol);
407
408 #ifdef IS_MPI
409
410 }
411
412 #endif
413
414 } //end for(int i=0)
415 }
416
417 void SimCreator::compList(MakeStamps *stamps, Globals* simParams,
418 std::vector < std::pair<MoleculeStamp *, int> > &moleculeStampPairs) {
419 int i;
420 char * id;
421 MoleculeStamp * currentStamp;
422 Component** the_components = simParams->getComponents();
423 int n_components = simParams->getNComponents();
424
425 if (!simParams->haveNMol()) {
426 // we don't have the total number of molecules, so we assume it is
427 // given in each component
428
429 for(i = 0; i < n_components; i++) {
430 if (!the_components[i]->haveNMol()) {
431 // we have a problem
432 sprintf(painCave.errMsg,
433 "SimCreator Error. No global NMol or component NMol given.\n"
434 "\tCannot calculate the number of atoms.\n");
435
436 painCave.isFatal = 1;
437 simError();
438 }
439
440 id = the_components[i]->getType();
441 currentStamp = (stamps->extractMolStamp(id))->getStamp();
442
443 if (currentStamp == NULL) {
444 sprintf(painCave.errMsg,
445 "SimCreator error: Component \"%s\" was not found in the "
446 "list of declared molecules\n", id);
447
448 painCave.isFatal = 1;
449 simError();
450 }
451
452 moleculeStampPairs.push_back(
453 std::make_pair(currentStamp, the_components[i]->getNMol()));
454 } //end for (i = 0; i < n_components; i++)
455 } else {
456 sprintf(painCave.errMsg, "SimSetup error.\n"
457 "\tSorry, the ability to specify total"
458 " nMols and then give molfractions in the components\n"
459 "\tis not currently supported."
460 " Please give nMol in the components.\n");
461
462 painCave.isFatal = 1;
463 simError();
464 }
465
466 #ifdef IS_MPI
467
468 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
469 MPIcheckPoint();
470
471 #endif // is_mpi
472
473 }
474
475 void SimCreator::setGlobalIndex(SimInfo *info) {
476 SimInfo::MoleculeIterator mi;
477 Molecule::AtomIterator ai;
478 Molecule::RigidBodyIterator ri;
479 Molecule::CutoffGroupIterator ci;
480 Molecule * mol;
481 Atom * atom;
482 RigidBody * rb;
483 CutoffGroup * cg;
484 int beginAtomIndex;
485 int beginRigidBodyIndex;
486 int beginCutoffGroupIndex;
487 int nGlobalAtoms = info->getNGlobalAtoms();
488
489 #ifndef IS_MPI
490
491 beginAtomIndex = 0;
492 beginRigidBodyIndex = 0;
493 beginCutoffGroupIndex = 0;
494
495 #else
496
497 int nproc;
498 int myNode;
499
500 myNode = worldRank;
501 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
502
503 std::vector < int > tmpAtomsInProc(nproc, 0);
504 std::vector < int > tmpRigidBodiesInProc(nproc, 0);
505 std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
506 std::vector < int > NumAtomsInProc(nproc, 0);
507 std::vector < int > NumRigidBodiesInProc(nproc, 0);
508 std::vector < int > NumCutoffGroupsInProc(nproc, 0);
509
510 tmpAtomsInProc[myNode] = info->getNAtoms();
511 tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
512 tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
513
514 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
515 MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
516 MPI_SUM, MPI_COMM_WORLD);
517 MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
518 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
519 MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
520 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
521
522 beginAtomIndex = 0;
523 beginRigidBodyIndex = 0;
524 beginCutoffGroupIndex = 0;
525
526 for(int i = 0; i < myNode; i++) {
527 beginAtomIndex += NumAtomsInProc[i];
528 beginRigidBodyIndex += NumRigidBodiesInProc[i];
529 beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
530 }
531
532 #endif
533
534 for(mol = info->beginMolecule(mi); mol != NULL;
535 mol = info->nextMolecule(mi)) {
536
537 //local index(index in DataStorge) of atom is important
538 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
539 atom->setGlobalIndex(beginAtomIndex++);
540 }
541
542 for(rb = mol->beginRigidBody(ri); rb != NULL;
543 rb = mol->nextRigidBody(ri)) {
544 rb->setGlobalIndex(beginRigidBodyIndex++);
545 }
546
547 //local index of cutoff group is trivial, it only depends on the order of travesing
548 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
549 cg = mol->nextCutoffGroup(ci)) {
550 cg->setGlobalIndex(beginCutoffGroupIndex++);
551 }
552 }
553
554 //fill globalGroupMembership
555 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
556 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
557 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
558
559 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
560 globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
561 }
562
563 }
564 }
565
566 #ifdef IS_MPI
567 // Since the globalGroupMembership has been zero filled and we've only
568 // poked values into the atoms we know, we can do an Allreduce
569 // to get the full globalGroupMembership array (We think).
570 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
571 // docs said we could.
572 std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
573 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
574 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
575 info->setGlobalGroupMembership(tmpGroupMembership);
576 #else
577 info->setGlobalGroupMembership(globalGroupMembership);
578 #endif
579
580 //fill molMembership
581 std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
582
583 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
584
585 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
586 globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
587 }
588 }
589
590 #ifdef IS_MPI
591 std::vector<int> tmpMolMembership(nGlobalAtoms, 0);
592
593 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
594 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
595
596 info->setGlobalMolMembership(tmpMolMembership);
597 #else
598 info->setGlobalMolMembership(globalMolMembership);
599 #endif
600
601 }
602
603 void SimCreator::loadCoordinates(SimInfo* info) {
604 Globals* simParams;
605 simParams = info->getSimParams();
606
607 if (!simParams->haveInitialConfig()) {
608 sprintf(painCave.errMsg,
609 "Cannot intialize a simulation without an initial configuration file.\n");
610 painCave.isFatal = 1;;
611 simError();
612 }
613
614 DumpReader reader(info, simParams->getInitialConfig());
615 int nframes = reader.getNFrames();
616
617 if (nframes > 0) {
618 reader.readFrame(nframes - 1);
619 } else {
620 //invalid initial coordinate file
621 sprintf(painCave.errMsg, "Initial configuration file %s should at least contain one frame\n",
622 simParams->getInitialConfig());
623 painCave.isFatal = 1;
624 simError();
625 }
626
627 }
628
629 } //end namespace oopse
630
631