ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/branches/new_design/OOPSE-4/src/brains/SimCreator.cpp
Revision: 1903
Committed: Thu Jan 6 00:16:07 2005 UTC (19 years, 7 months ago) by tim
File size: 19633 byte(s)
Log Message:
simpleBuilder in progress

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, bool loadInitCoords) {
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 if (loadInitCoords)
165 loadCoordinates(info);
166
167 return info;
168 }
169
170 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
171
172 //setup seed for random number generator
173 int seedValue;
174 Globals * simParams = info->getSimParams();
175
176 if (simParams->haveSeed()) {
177 seedValue = simParams->getSeed();
178
179 if (seedValue < 100000000 ) {
180 sprintf(painCave.errMsg,
181 "Seed for sprng library should contain at least 9 digits\n"
182 "OOPSE will generate a seed for user\n");
183
184 painCave.isFatal = 0;
185 simError();
186
187 //using seed generated by system instead of invalid seed set by user
188
189 #ifndef IS_MPI
190
191 seedValue = make_sprng_seed();
192
193 #else
194
195 if (worldRank == 0) {
196 seedValue = make_sprng_seed();
197 }
198
199 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
200
201 #endif
202
203 } //end if (seedValue /1000000000 == 0)
204 } else {
205
206 #ifndef IS_MPI
207
208 seedValue = make_sprng_seed();
209
210 #else
211
212 if (worldRank == 0) {
213 seedValue = make_sprng_seed();
214 }
215
216 MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD);
217
218 #endif
219
220 } //end of simParams->haveSeed()
221
222 info->setSeed(seedValue);
223
224
225 //figure out the ouput file names
226 std::string prefix;
227
228 #ifdef IS_MPI
229
230 if (worldRank == 0) {
231 #endif // is_mpi
232
233 if (simParams->haveFinalConfig()) {
234 prefix = getPrefix(simParams->getFinalConfig());
235 } else {
236 prefix = getPrefix(mdfile);
237 }
238
239 info->setFinalConfigFileName(prefix + ".eor");
240 info->setDumpFileName(prefix + ".dump");
241 info->setStatFileName(prefix + ".stat");
242
243 #ifdef IS_MPI
244
245 }
246
247 #endif
248
249 }
250
251 #ifdef IS_MPI
252 void SimCreator::divideMolecules(SimInfo *info) {
253 double numerator;
254 double denominator;
255 double precast;
256 double x;
257 double y;
258 double a;
259 int old_atoms;
260 int add_atoms;
261 int new_atoms;
262 int nTarget;
263 int done;
264 int i;
265 int j;
266 int loops;
267 int which_proc;
268 int nProcessors;
269 std::vector<int> atomsPerProc;
270 randomSPRNG myRandom(info->getSeed());
271 int nGlobalMols = info->getNGlobalMolecules();
272 std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
273
274 MPI_Comm_size(MPI_COMM_WORLD, &nProcessors);
275
276 if (nProcessors > nGlobalMols) {
277 sprintf(painCave.errMsg,
278 "nProcessors (%d) > nMol (%d)\n"
279 "\tThe number of processors is larger than\n"
280 "\tthe number of molecules. This will not result in a \n"
281 "\tusable division of atoms for force decomposition.\n"
282 "\tEither try a smaller number of processors, or run the\n"
283 "\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols);
284
285 painCave.isFatal = 1;
286 simError();
287 }
288
289 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
290
291 //initialize atomsPerProc
292 atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
293
294 if (worldRank == 0) {
295 numerator = info->getNGlobalAtoms();
296 denominator = nProcessors;
297 precast = numerator / denominator;
298 nTarget = (int)(precast + 0.5);
299
300 for(i = 0; i < nGlobalMols; i++) {
301 done = 0;
302 loops = 0;
303
304 while (!done) {
305 loops++;
306
307 // Pick a processor at random
308
309 which_proc = (int) (myRandom.getRandom() * nProcessors);
310
311 //get the molecule stamp first
312 int stampId = info->getMoleculeStampId(i);
313 MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
314
315 // How many atoms does this processor have so far?
316 old_atoms = atomsPerProc[which_proc];
317 add_atoms = moleculeStamp->getNAtoms();
318 new_atoms = old_atoms + add_atoms;
319
320 // If we've been through this loop too many times, we need
321 // to just give up and assign the molecule to this processor
322 // and be done with it.
323
324 if (loops > 100) {
325 sprintf(painCave.errMsg,
326 "I've tried 100 times to assign molecule %d to a "
327 " processor, but can't find a good spot.\n"
328 "I'm assigning it at random to processor %d.\n",
329 i, which_proc);
330
331 painCave.isFatal = 0;
332 simError();
333
334 molToProcMap[i] = which_proc;
335 atomsPerProc[which_proc] += add_atoms;
336
337 done = 1;
338 continue;
339 }
340
341 // If we can add this molecule to this processor without sending
342 // it above nTarget, then go ahead and do it:
343
344 if (new_atoms <= nTarget) {
345 molToProcMap[i] = which_proc;
346 atomsPerProc[which_proc] += add_atoms;
347
348 done = 1;
349 continue;
350 }
351
352 // The only situation left is when new_atoms > nTarget. We
353 // want to accept this with some probability that dies off the
354 // farther we are from nTarget
355
356 // roughly: x = new_atoms - nTarget
357 // Pacc(x) = exp(- a * x)
358 // where a = penalty / (average atoms per molecule)
359
360 x = (double)(new_atoms - nTarget);
361 y = myRandom.getRandom();
362
363 if (y < exp(- a * x)) {
364 molToProcMap[i] = which_proc;
365 atomsPerProc[which_proc] += add_atoms;
366
367 done = 1;
368 continue;
369 } else {
370 continue;
371 }
372 }
373 }
374
375 // Spray out this nonsense to all other processors:
376
377 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
378 } else {
379
380 // Listen to your marching orders from processor 0:
381
382 MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD);
383 }
384
385 info->setMolToProcMap(molToProcMap);
386 sprintf(checkPointMsg,
387 "Successfully divided the molecules among the processors.\n");
388 MPIcheckPoint();
389 }
390
391 #endif
392
393 void SimCreator::createMolecules(SimInfo *info) {
394 MoleculeCreator molCreator;
395 int stampId;
396
397 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
398
399 #ifdef IS_MPI
400
401 if (info->getMolToProc(i) == worldRank) {
402 #endif
403
404 stampId = info->getMoleculeStampId(i);
405 Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId),
406 stampId, i, info->getLocalIndexManager());
407
408 info->addMolecule(mol);
409
410 #ifdef IS_MPI
411
412 }
413
414 #endif
415
416 } //end for(int i=0)
417 }
418
419 void SimCreator::compList(MakeStamps *stamps, Globals* simParams,
420 std::vector < std::pair<MoleculeStamp *, int> > &moleculeStampPairs) {
421 int i;
422 char * id;
423 MoleculeStamp * currentStamp;
424 Component** the_components = simParams->getComponents();
425 int n_components = simParams->getNComponents();
426
427 if (!simParams->haveNMol()) {
428 // we don't have the total number of molecules, so we assume it is
429 // given in each component
430
431 for(i = 0; i < n_components; i++) {
432 if (!the_components[i]->haveNMol()) {
433 // we have a problem
434 sprintf(painCave.errMsg,
435 "SimCreator Error. No global NMol or component NMol given.\n"
436 "\tCannot calculate the number of atoms.\n");
437
438 painCave.isFatal = 1;
439 simError();
440 }
441
442 id = the_components[i]->getType();
443 currentStamp = (stamps->extractMolStamp(id))->getStamp();
444
445 if (currentStamp == NULL) {
446 sprintf(painCave.errMsg,
447 "SimCreator error: Component \"%s\" was not found in the "
448 "list of declared molecules\n", id);
449
450 painCave.isFatal = 1;
451 simError();
452 }
453
454 moleculeStampPairs.push_back(
455 std::make_pair(currentStamp, the_components[i]->getNMol()));
456 } //end for (i = 0; i < n_components; i++)
457 } else {
458 sprintf(painCave.errMsg, "SimSetup error.\n"
459 "\tSorry, the ability to specify total"
460 " nMols and then give molfractions in the components\n"
461 "\tis not currently supported."
462 " Please give nMol in the components.\n");
463
464 painCave.isFatal = 1;
465 simError();
466 }
467
468 #ifdef IS_MPI
469
470 strcpy(checkPointMsg, "Component stamps successfully extracted\n");
471 MPIcheckPoint();
472
473 #endif // is_mpi
474
475 }
476
477 void SimCreator::setGlobalIndex(SimInfo *info) {
478 SimInfo::MoleculeIterator mi;
479 Molecule::AtomIterator ai;
480 Molecule::RigidBodyIterator ri;
481 Molecule::CutoffGroupIterator ci;
482 Molecule * mol;
483 Atom * atom;
484 RigidBody * rb;
485 CutoffGroup * cg;
486 int beginAtomIndex;
487 int beginRigidBodyIndex;
488 int beginCutoffGroupIndex;
489 int nGlobalAtoms = info->getNGlobalAtoms();
490
491 #ifndef IS_MPI
492
493 beginAtomIndex = 0;
494 beginRigidBodyIndex = 0;
495 beginCutoffGroupIndex = 0;
496
497 #else
498
499 int nproc;
500 int myNode;
501
502 myNode = worldRank;
503 MPI_Comm_size(MPI_COMM_WORLD, &nproc);
504
505 std::vector < int > tmpAtomsInProc(nproc, 0);
506 std::vector < int > tmpRigidBodiesInProc(nproc, 0);
507 std::vector < int > tmpCutoffGroupsInProc(nproc, 0);
508 std::vector < int > NumAtomsInProc(nproc, 0);
509 std::vector < int > NumRigidBodiesInProc(nproc, 0);
510 std::vector < int > NumCutoffGroupsInProc(nproc, 0);
511
512 tmpAtomsInProc[myNode] = info->getNAtoms();
513 tmpRigidBodiesInProc[myNode] = info->getNRigidBodies();
514 tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups();
515
516 //do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups
517 MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT,
518 MPI_SUM, MPI_COMM_WORLD);
519 MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc,
520 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
521 MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc,
522 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
523
524 beginAtomIndex = 0;
525 beginRigidBodyIndex = 0;
526 beginCutoffGroupIndex = 0;
527
528 for(int i = 0; i < myNode; i++) {
529 beginAtomIndex += NumAtomsInProc[i];
530 beginRigidBodyIndex += NumRigidBodiesInProc[i];
531 beginCutoffGroupIndex += NumCutoffGroupsInProc[i];
532 }
533
534 #endif
535
536 for(mol = info->beginMolecule(mi); mol != NULL;
537 mol = info->nextMolecule(mi)) {
538
539 //local index(index in DataStorge) of atom is important
540 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
541 atom->setGlobalIndex(beginAtomIndex++);
542 }
543
544 for(rb = mol->beginRigidBody(ri); rb != NULL;
545 rb = mol->nextRigidBody(ri)) {
546 rb->setGlobalIndex(beginRigidBodyIndex++);
547 }
548
549 //local index of cutoff group is trivial, it only depends on the order of travesing
550 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
551 cg = mol->nextCutoffGroup(ci)) {
552 cg->setGlobalIndex(beginCutoffGroupIndex++);
553 }
554 }
555
556 //fill globalGroupMembership
557 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
558 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
559 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
560
561 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
562 globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
563 }
564
565 }
566 }
567
568 #ifdef IS_MPI
569 // Since the globalGroupMembership has been zero filled and we've only
570 // poked values into the atoms we know, we can do an Allreduce
571 // to get the full globalGroupMembership array (We think).
572 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
573 // docs said we could.
574 std::vector<int> tmpGroupMembership(nGlobalAtoms, 0);
575 MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms,
576 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
577 info->setGlobalGroupMembership(tmpGroupMembership);
578 #else
579 info->setGlobalGroupMembership(globalGroupMembership);
580 #endif
581
582 //fill molMembership
583 std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
584
585 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
586
587 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
588 globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
589 }
590 }
591
592 #ifdef IS_MPI
593 std::vector<int> tmpMolMembership(nGlobalAtoms, 0);
594
595 MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms,
596 MPI_INT, MPI_SUM, MPI_COMM_WORLD);
597
598 info->setGlobalMolMembership(tmpMolMembership);
599 #else
600 info->setGlobalMolMembership(globalMolMembership);
601 #endif
602
603 }
604
605 void SimCreator::loadCoordinates(SimInfo* info) {
606 Globals* simParams;
607 simParams = info->getSimParams();
608
609 if (!simParams->haveInitialConfig()) {
610 sprintf(painCave.errMsg,
611 "Cannot intialize a simulation without an initial configuration file.\n");
612 painCave.isFatal = 1;;
613 simError();
614 }
615
616 DumpReader reader(info, simParams->getInitialConfig());
617 int nframes = reader.getNFrames();
618
619 if (nframes > 0) {
620 reader.readFrame(nframes - 1);
621 } else {
622 //invalid initial coordinate file
623 sprintf(painCave.errMsg, "Initial configuration file %s should at least contain one frame\n",
624 simParams->getInitialConfig());
625 painCave.isFatal = 1;
626 simError();
627 }
628
629 }
630
631 } //end namespace oopse
632
633