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: 1856
Committed: Mon Dec 6 04:49:53 2004 UTC (19 years, 7 months ago) by tim
File size: 19494 byte(s)
Log Message:
fix a bug in Exclude List

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