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: 1740
Committed: Mon Nov 15 23:00:32 2004 UTC (19 years, 8 months ago) by tim
File size: 17909 byte(s)
Log Message:
adding ForceFieldFactory and LJFF classes

File Contents

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