ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1802
Committed: Wed Oct 3 14:07:28 2012 UTC (12 years, 6 months ago) by gezelter
File size: 31571 byte(s)
Log Message:
Parallel fixes for selection syntax (e.g. DistanceFinder).

File Contents

# User Rev Content
1 gezelter 403 /*
2 gezelter 1793 * copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3 gezelter 403 *
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 gezelter 1390 * 1. Redistributions of source code must retain the above copyright
10 gezelter 403 * notice, this list of conditions and the following disclaimer.
11     *
12 gezelter 1390 * 2. Redistributions in binary form must reproduce the above copyright
13 gezelter 403 * notice, this list of conditions and the following disclaimer in the
14     * documentation and/or other materials provided with the
15     * distribution.
16     *
17     * This software is provided "AS IS," without a warranty of any
18     * kind. All express or implied conditions, representations and
19     * warranties, including any implied warranty of merchantability,
20     * fitness for a particular purpose or non-infringement, are hereby
21     * excluded. The University of Notre Dame and its licensors shall not
22     * be liable for any damages suffered by licensee as a result of
23     * using, modifying or distributing the software or its
24     * derivatives. In no event will the University of Notre Dame or its
25     * licensors be liable for any lost revenue, profit or data, or for
26     * direct, indirect, special, consequential, incidental or punitive
27     * damages, however caused and regardless of the theory of liability,
28     * arising out of the use of or inability to use software, even if the
29     * University of Notre Dame has been advised of the possibility of
30     * such damages.
31 gezelter 1390 *
32     * SUPPORT OPEN SCIENCE! If you use OpenMD or its source code in your
33     * research, please cite the appropriate papers when you publish your
34     * work. Good starting points are:
35     *
36     * [1] Meineke, et al., J. Comp. Chem. 26, 252-271 (2005).
37     * [2] Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006).
38     * [3] Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008).
39 gezelter 1782 * [4] Kuang & Gezelter, J. Chem. Phys. 133, 164101 (2010).
40     * [5] Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011).
41 gezelter 403 */
42    
43     /**
44     * @file SimCreator.cpp
45     * @author tlin
46     * @date 11/03/2004
47     * @time 13:51am
48     * @version 1.0
49     */
50 tim 845 #include <exception>
51 tim 770 #include <iostream>
52     #include <sstream>
53     #include <string>
54    
55 gezelter 403 #include "brains/MoleculeCreator.hpp"
56     #include "brains/SimCreator.hpp"
57     #include "brains/SimSnapshotManager.hpp"
58     #include "io/DumpReader.hpp"
59 gezelter 1782 #include "brains/ForceField.hpp"
60 gezelter 403 #include "utils/simError.h"
61     #include "utils/StringUtils.hpp"
62     #include "math/SeqRandNumGen.hpp"
63 tim 770 #include "mdParser/MDLexer.hpp"
64     #include "mdParser/MDParser.hpp"
65     #include "mdParser/MDTreeParser.hpp"
66     #include "mdParser/SimplePreprocessor.hpp"
67 tim 816 #include "antlr/ANTLRException.hpp"
68     #include "antlr/TokenStreamRecognitionException.hpp"
69     #include "antlr/TokenStreamIOException.hpp"
70     #include "antlr/TokenStreamException.hpp"
71     #include "antlr/RecognitionException.hpp"
72     #include "antlr/CharStreamException.hpp"
73 tim 770
74 tim 816 #include "antlr/MismatchedCharException.hpp"
75     #include "antlr/MismatchedTokenException.hpp"
76     #include "antlr/NoViableAltForCharException.hpp"
77     #include "antlr/NoViableAltException.hpp"
78 tim 770
79 gezelter 1782 #include "types/DirectionalAdapter.hpp"
80     #include "types/MultipoleAdapter.hpp"
81     #include "types/EAMAdapter.hpp"
82     #include "types/SuttonChenAdapter.hpp"
83     #include "types/PolarizableAdapter.hpp"
84     #include "types/FixedChargeAdapter.hpp"
85     #include "types/FluctuatingChargeAdapter.hpp"
86    
87 gezelter 403 #ifdef IS_MPI
88 gezelter 1782 #include "mpi.h"
89 gezelter 403 #include "math/ParallelRandNumGen.hpp"
90     #endif
91    
92 gezelter 1390 namespace OpenMD {
93 chrisfen 417
94 gezelter 1782 Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int mdFileVersion, int startOfMetaDataBlock ){
95 tim 1024 Globals* simParams = NULL;
96     try {
97 tim 770
98 tim 1024 // Create a preprocessor that preprocesses md file into an ostringstream
99     std::stringstream ppStream;
100 tim 770 #ifdef IS_MPI
101 tim 1024 int streamSize;
102     const int masterNode = 0;
103 gezelter 1793
104 tim 1024 if (worldRank == masterNode) {
105 gezelter 1793 MPI::COMM_WORLD.Bcast(&mdFileVersion, 1, MPI::INT, masterNode);
106 gezelter 1782 #endif
107 tim 1024 SimplePreprocessor preprocessor;
108 gezelter 1793 preprocessor.preprocess(rawMetaDataStream, filename, startOfMetaDataBlock,
109     ppStream);
110 tim 770
111     #ifdef IS_MPI
112 tim 1024 //brocasting the stream size
113     streamSize = ppStream.str().size() +1;
114 gezelter 1793 MPI::COMM_WORLD.Bcast(&streamSize, 1, MPI::LONG, masterNode);
115     MPI::COMM_WORLD.Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())),
116     streamSize, MPI::CHAR, masterNode);
117 tim 770
118 tim 1024 } else {
119 gezelter 1782
120 gezelter 1793 MPI::COMM_WORLD.Bcast(&mdFileVersion, 1, MPI::INT, masterNode);
121 gezelter 1782
122 tim 1024 //get stream size
123 gezelter 1793 MPI::COMM_WORLD.Bcast(&streamSize, 1, MPI::LONG, masterNode);
124 gezelter 1313
125 tim 1024 char* buf = new char[streamSize];
126     assert(buf);
127 tim 770
128 tim 1024 //receive file content
129 gezelter 1793 MPI::COMM_WORLD.Bcast(buf, streamSize, MPI::CHAR, masterNode);
130 tim 770
131 tim 1024 ppStream.str(buf);
132 gezelter 1313 delete [] buf;
133 tim 770
134 tim 1024 }
135 tim 770 #endif
136 tim 1024 // Create a scanner that reads from the input stream
137     MDLexer lexer(ppStream);
138     lexer.setFilename(filename);
139     lexer.initDeferredLineCount();
140 chrisfen 417
141 tim 1024 // Create a parser that reads from the scanner
142     MDParser parser(lexer);
143     parser.setFilename(filename);
144 tim 770
145 tim 1024 // Create an observer that synchorizes file name change
146     FilenameObserver observer;
147     observer.setLexer(&lexer);
148     observer.setParser(&parser);
149     lexer.setObserver(&observer);
150 chrisfen 417
151 tim 1024 antlr::ASTFactory factory;
152     parser.initializeASTFactory(factory);
153     parser.setASTFactory(&factory);
154     parser.mdfile();
155 tim 770
156 tim 1024 // Create a tree parser that reads information into Globals
157     MDTreeParser treeParser;
158     treeParser.initializeASTFactory(factory);
159     treeParser.setASTFactory(&factory);
160     simParams = treeParser.walkTree(parser.getAST());
161     }
162 tim 845
163 tim 816
164 tim 1024 catch(antlr::MismatchedCharException& e) {
165     sprintf(painCave.errMsg,
166     "parser exception: %s %s:%d:%d\n",
167     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
168     painCave.isFatal = 1;
169     simError();
170     }
171     catch(antlr::MismatchedTokenException &e) {
172     sprintf(painCave.errMsg,
173     "parser exception: %s %s:%d:%d\n",
174     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
175     painCave.isFatal = 1;
176     simError();
177     }
178     catch(antlr::NoViableAltForCharException &e) {
179     sprintf(painCave.errMsg,
180     "parser exception: %s %s:%d:%d\n",
181     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
182     painCave.isFatal = 1;
183     simError();
184     }
185     catch(antlr::NoViableAltException &e) {
186     sprintf(painCave.errMsg,
187     "parser exception: %s %s:%d:%d\n",
188     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
189     painCave.isFatal = 1;
190     simError();
191     }
192 tim 845
193 tim 1024 catch(antlr::TokenStreamRecognitionException& e) {
194     sprintf(painCave.errMsg,
195     "parser exception: %s %s:%d:%d\n",
196     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
197     painCave.isFatal = 1;
198     simError();
199     }
200 tim 845
201 tim 1024 catch(antlr::TokenStreamIOException& e) {
202     sprintf(painCave.errMsg,
203     "parser exception: %s\n",
204     e.getMessage().c_str());
205     painCave.isFatal = 1;
206     simError();
207     }
208 tim 845
209 tim 1024 catch(antlr::TokenStreamException& e) {
210     sprintf(painCave.errMsg,
211     "parser exception: %s\n",
212     e.getMessage().c_str());
213     painCave.isFatal = 1;
214     simError();
215     }
216     catch (antlr::RecognitionException& e) {
217     sprintf(painCave.errMsg,
218     "parser exception: %s %s:%d:%d\n",
219     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
220     painCave.isFatal = 1;
221     simError();
222     }
223     catch (antlr::CharStreamException& e) {
224     sprintf(painCave.errMsg,
225     "parser exception: %s\n",
226     e.getMessage().c_str());
227     painCave.isFatal = 1;
228     simError();
229     }
230 gezelter 1390 catch (OpenMDException& e) {
231 tim 1024 sprintf(painCave.errMsg,
232     "%s\n",
233     e.getMessage().c_str());
234     painCave.isFatal = 1;
235     simError();
236     }
237     catch (std::exception& e) {
238     sprintf(painCave.errMsg,
239     "parser exception: %s\n",
240     e.what());
241     painCave.isFatal = 1;
242     simError();
243     }
244 tim 770
245 gezelter 1782 simParams->setMDfileVersion(mdFileVersion);
246 tim 1024 return simParams;
247 gezelter 403 }
248 chrisfen 417
249 chrisfen 514 SimInfo* SimCreator::createSim(const std::string & mdFileName,
250     bool loadInitCoords) {
251 gezelter 1782
252 tim 1024 const int bufferSize = 65535;
253     char buffer[bufferSize];
254     int lineNo = 0;
255     std::string mdRawData;
256     int metaDataBlockStart = -1;
257     int metaDataBlockEnd = -1;
258     int i;
259 gezelter 1790 streamoff mdOffset(0);
260 gezelter 1782 int mdFileVersion;
261 tim 1024
262 gezelter 1790
263 tim 1024 #ifdef IS_MPI
264     const int masterNode = 0;
265     if (worldRank == masterNode) {
266     #endif
267    
268 gezelter 1790 std::ifstream mdFile_;
269     mdFile_.open(mdFileName.c_str(), ifstream::in | ifstream::binary);
270 tim 1024
271     if (mdFile_.fail()) {
272     sprintf(painCave.errMsg,
273     "SimCreator: Cannot open file: %s\n",
274     mdFileName.c_str());
275     painCave.isFatal = 1;
276     simError();
277     }
278    
279     mdFile_.getline(buffer, bufferSize);
280     ++lineNo;
281     std::string line = trimLeftCopy(buffer);
282 gezelter 1390 i = CaseInsensitiveFind(line, "<OpenMD");
283 gezelter 1287 if (static_cast<size_t>(i) == string::npos) {
284 gezelter 1390 // try the older file strings to see if that works:
285     i = CaseInsensitiveFind(line, "<OOPSE");
286     }
287    
288     if (static_cast<size_t>(i) == string::npos) {
289     // still no luck!
290 tim 1024 sprintf(painCave.errMsg,
291 gezelter 1390 "SimCreator: File: %s is not a valid OpenMD file!\n",
292 tim 1024 mdFileName.c_str());
293     painCave.isFatal = 1;
294     simError();
295     }
296 gezelter 1782
297     // found the correct opening string, now try to get the file
298     // format version number.
299 tim 1024
300 gezelter 1782 StringTokenizer tokenizer(line, "=<> \t\n\r");
301     std::string fileType = tokenizer.nextToken();
302     toUpper(fileType);
303    
304     mdFileVersion = 0;
305    
306     if (fileType == "OPENMD") {
307     while (tokenizer.hasMoreTokens()) {
308     std::string token(tokenizer.nextToken());
309     toUpper(token);
310     if (token == "VERSION") {
311     mdFileVersion = tokenizer.nextTokenAsInt();
312     break;
313     }
314     }
315     }
316    
317 tim 1024 //scan through the input stream and find MetaData tag
318     while(mdFile_.getline(buffer, bufferSize)) {
319     ++lineNo;
320    
321     std::string line = trimLeftCopy(buffer);
322     if (metaDataBlockStart == -1) {
323     i = CaseInsensitiveFind(line, "<MetaData>");
324     if (i != string::npos) {
325     metaDataBlockStart = lineNo;
326     mdOffset = mdFile_.tellg();
327     }
328     } else {
329     i = CaseInsensitiveFind(line, "</MetaData>");
330     if (i != string::npos) {
331     metaDataBlockEnd = lineNo;
332     }
333     }
334     }
335    
336     if (metaDataBlockStart == -1) {
337     sprintf(painCave.errMsg,
338     "SimCreator: File: %s did not contain a <MetaData> tag!\n",
339     mdFileName.c_str());
340     painCave.isFatal = 1;
341     simError();
342     }
343     if (metaDataBlockEnd == -1) {
344     sprintf(painCave.errMsg,
345     "SimCreator: File: %s did not contain a closed MetaData block!\n",
346     mdFileName.c_str());
347     painCave.isFatal = 1;
348     simError();
349     }
350    
351     mdFile_.clear();
352     mdFile_.seekg(0);
353     mdFile_.seekg(mdOffset);
354    
355     mdRawData.clear();
356    
357     for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
358     mdFile_.getline(buffer, bufferSize);
359     mdRawData += buffer;
360     mdRawData += "\n";
361     }
362    
363     mdFile_.close();
364    
365     #ifdef IS_MPI
366     }
367     #endif
368    
369     std::stringstream rawMetaDataStream(mdRawData);
370    
371 gezelter 403 //parse meta-data file
372 gezelter 1782 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
373     metaDataBlockStart + 1);
374 chrisfen 417
375 gezelter 403 //create the force field
376 gezelter 1782 ForceField * ff = new ForceField(simParams->getForceField());
377 gezelter 1277
378 gezelter 403 if (ff == NULL) {
379 chrisfen 514 sprintf(painCave.errMsg,
380     "ForceField Factory can not create %s force field\n",
381 tim 665 simParams->getForceField().c_str());
382 gezelter 403 painCave.isFatal = 1;
383     simError();
384     }
385 chrisfen 417
386 gezelter 403 if (simParams->haveForceFieldFileName()) {
387     ff->setForceFieldFileName(simParams->getForceFieldFileName());
388     }
389    
390     std::string forcefieldFileName;
391     forcefieldFileName = ff->getForceFieldFileName();
392 chrisfen 417
393 gezelter 403 if (simParams->haveForceFieldVariant()) {
394     //If the force field has variant, the variant force field name will be
395     //Base.variant.frc. For exampel EAM.u6.frc
396 chrisfen 417
397 gezelter 403 std::string variant = simParams->getForceFieldVariant();
398 chrisfen 417
399 gezelter 403 std::string::size_type pos = forcefieldFileName.rfind(".frc");
400     variant = "." + variant;
401     if (pos != std::string::npos) {
402     forcefieldFileName.insert(pos, variant);
403     } else {
404     //If the default force field file name does not containt .frc suffix, just append the .variant
405     forcefieldFileName.append(variant);
406     }
407     }
408    
409     ff->parse(forcefieldFileName);
410     //create SimInfo
411 tim 770 SimInfo * info = new SimInfo(ff, simParams);
412 tim 1024
413     info->setRawMetaData(mdRawData);
414 tim 490
415 gezelter 945 //gather parameters (SimCreator only retrieves part of the
416     //parameters)
417 gezelter 403 gatherParameters(info, mdFileName);
418 chrisfen 417
419 gezelter 403 //divide the molecules and determine the global index of molecules
420     #ifdef IS_MPI
421     divideMolecules(info);
422     #endif
423 chrisfen 417
424 gezelter 403 //create the molecules
425     createMolecules(info);
426 chrisfen 417
427 gezelter 1782 //find the storage layout
428    
429     int storageLayout = computeStorageLayout(info);
430    
431 gezelter 945 //allocate memory for DataStorage(circular reference, need to
432     //break it)
433 gezelter 1782 info->setSnapshotManager(new SimSnapshotManager(info, storageLayout));
434 gezelter 403
435 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
436     //(only need to be set once, the global index will never change
437     //again). Local indices of atoms and rigidbodies are already set
438     //by MoleculeCreator class which actually delegates the
439     //responsibility to LocalIndexManager.
440 gezelter 403 setGlobalIndex(info);
441 chrisfen 417
442 gezelter 1287 //Although addInteractionPairs is called inside SimInfo's addMolecule
443 gezelter 945 //method, at that point atoms don't have the global index yet
444     //(their global index are all initialized to -1). Therefore we
445 gezelter 1287 //have to call addInteractionPairs explicitly here. A way to work
446 gezelter 945 //around is that we can determine the beginning global indices of
447     //atoms before they get created.
448 gezelter 403 SimInfo::MoleculeIterator mi;
449     Molecule* mol;
450     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
451 gezelter 1287 info->addInteractionPairs(mol);
452 gezelter 403 }
453    
454     if (loadInitCoords)
455 tim 1024 loadCoordinates(info, mdFileName);
456 gezelter 403 return info;
457     }
458 chrisfen 417
459 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
460 chrisfen 417
461 tim 749 //figure out the output file names
462 gezelter 403 std::string prefix;
463 chrisfen 417
464 gezelter 403 #ifdef IS_MPI
465 chrisfen 417
466 gezelter 403 if (worldRank == 0) {
467     #endif // is_mpi
468     Globals * simParams = info->getSimParams();
469     if (simParams->haveFinalConfig()) {
470     prefix = getPrefix(simParams->getFinalConfig());
471     } else {
472     prefix = getPrefix(mdfile);
473     }
474 chrisfen 417
475 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
476     info->setDumpFileName(prefix + ".dump");
477     info->setStatFileName(prefix + ".stat");
478 chrisfen 417 info->setRestFileName(prefix + ".zang");
479    
480 gezelter 403 #ifdef IS_MPI
481 chrisfen 417
482 gezelter 403 }
483 chrisfen 417
484 gezelter 403 #endif
485 chrisfen 417
486 gezelter 403 }
487 chrisfen 417
488 gezelter 403 #ifdef IS_MPI
489     void SimCreator::divideMolecules(SimInfo *info) {
490 tim 963 RealType numerator;
491     RealType denominator;
492     RealType precast;
493     RealType x;
494     RealType y;
495     RealType a;
496 gezelter 403 int old_atoms;
497     int add_atoms;
498     int new_atoms;
499     int nTarget;
500     int done;
501     int i;
502     int loops;
503     int which_proc;
504     int nProcessors;
505     std::vector<int> atomsPerProc;
506     int nGlobalMols = info->getNGlobalMolecules();
507     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
508    
509 gezelter 1796 nProcessors = MPI::COMM_WORLD.Get_size();
510 chrisfen 417
511 gezelter 403 if (nProcessors > nGlobalMols) {
512     sprintf(painCave.errMsg,
513     "nProcessors (%d) > nMol (%d)\n"
514     "\tThe number of processors is larger than\n"
515     "\tthe number of molecules. This will not result in a \n"
516     "\tusable division of atoms for force decomposition.\n"
517     "\tEither try a smaller number of processors, or run the\n"
518 gezelter 1390 "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols);
519 chrisfen 417
520 gezelter 403 painCave.isFatal = 1;
521     simError();
522     }
523 chrisfen 417
524 gezelter 403 int seedValue;
525     Globals * simParams = info->getSimParams();
526     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
527     if (simParams->haveSeed()) {
528     seedValue = simParams->getSeed();
529     myRandom = new SeqRandNumGen(seedValue);
530     }else {
531     myRandom = new SeqRandNumGen();
532     }
533 chrisfen 417
534    
535 gezelter 403 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
536 chrisfen 417
537 gezelter 403 //initialize atomsPerProc
538     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
539 chrisfen 417
540 gezelter 403 if (worldRank == 0) {
541     numerator = info->getNGlobalAtoms();
542     denominator = nProcessors;
543     precast = numerator / denominator;
544     nTarget = (int)(precast + 0.5);
545 chrisfen 417
546 gezelter 403 for(i = 0; i < nGlobalMols; i++) {
547 gezelter 1801
548 gezelter 403 done = 0;
549     loops = 0;
550 chrisfen 417
551 gezelter 403 while (!done) {
552     loops++;
553 chrisfen 417
554 gezelter 403 // Pick a processor at random
555 chrisfen 417
556 gezelter 403 which_proc = (int) (myRandom->rand() * nProcessors);
557 chrisfen 417
558 gezelter 403 //get the molecule stamp first
559     int stampId = info->getMoleculeStampId(i);
560     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
561 chrisfen 417
562 gezelter 403 // How many atoms does this processor have so far?
563     old_atoms = atomsPerProc[which_proc];
564     add_atoms = moleculeStamp->getNAtoms();
565     new_atoms = old_atoms + add_atoms;
566 chrisfen 417
567 gezelter 403 // If we've been through this loop too many times, we need
568     // to just give up and assign the molecule to this processor
569     // and be done with it.
570 chrisfen 417
571 gezelter 403 if (loops > 100) {
572 gezelter 1801
573 gezelter 403 sprintf(painCave.errMsg,
574 gezelter 1801 "There have been 100 attempts to assign molecule %d to an\n"
575     "\tunderworked processor, but there's no good place to\n"
576     "\tleave it. OpenMD is assigning it at random to processor %d.\n",
577 gezelter 403 i, which_proc);
578 gezelter 1801
579 gezelter 403 painCave.isFatal = 0;
580 gezelter 1801 painCave.severity = OPENMD_INFO;
581 gezelter 403 simError();
582 chrisfen 417
583 gezelter 403 molToProcMap[i] = which_proc;
584     atomsPerProc[which_proc] += add_atoms;
585 chrisfen 417
586 gezelter 403 done = 1;
587     continue;
588     }
589 chrisfen 417
590 gezelter 403 // If we can add this molecule to this processor without sending
591     // it above nTarget, then go ahead and do it:
592 chrisfen 417
593 gezelter 403 if (new_atoms <= nTarget) {
594     molToProcMap[i] = which_proc;
595     atomsPerProc[which_proc] += add_atoms;
596 chrisfen 417
597 gezelter 403 done = 1;
598     continue;
599     }
600 chrisfen 417
601 gezelter 403 // The only situation left is when new_atoms > nTarget. We
602     // want to accept this with some probability that dies off the
603     // farther we are from nTarget
604 chrisfen 417
605 gezelter 403 // roughly: x = new_atoms - nTarget
606     // Pacc(x) = exp(- a * x)
607     // where a = penalty / (average atoms per molecule)
608 chrisfen 417
609 tim 963 x = (RealType)(new_atoms - nTarget);
610 gezelter 403 y = myRandom->rand();
611 chrisfen 417
612 gezelter 403 if (y < exp(- a * x)) {
613     molToProcMap[i] = which_proc;
614     atomsPerProc[which_proc] += add_atoms;
615 chrisfen 417
616 gezelter 403 done = 1;
617     continue;
618     } else {
619     continue;
620     }
621     }
622     }
623 chrisfen 417
624 gezelter 403 delete myRandom;
625 gezelter 1801
626 gezelter 403 // Spray out this nonsense to all other processors:
627 gezelter 1796 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
628 gezelter 403 } else {
629 chrisfen 417
630 gezelter 403 // Listen to your marching orders from processor 0:
631 gezelter 1796 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
632 gezelter 1801
633 gezelter 403 }
634 chrisfen 417
635 gezelter 403 info->setMolToProcMap(molToProcMap);
636     sprintf(checkPointMsg,
637     "Successfully divided the molecules among the processors.\n");
638 gezelter 1241 errorCheckPoint();
639 gezelter 403 }
640 chrisfen 417
641 gezelter 403 #endif
642 chrisfen 417
643 gezelter 403 void SimCreator::createMolecules(SimInfo *info) {
644     MoleculeCreator molCreator;
645     int stampId;
646 chrisfen 417
647 gezelter 403 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
648 chrisfen 417
649 gezelter 403 #ifdef IS_MPI
650 chrisfen 417
651 gezelter 403 if (info->getMolToProc(i) == worldRank) {
652     #endif
653 chrisfen 417
654 gezelter 403 stampId = info->getMoleculeStampId(i);
655 gezelter 1782 Molecule * mol = molCreator.createMolecule(info->getForceField(),
656     info->getMoleculeStamp(stampId),
657     stampId, i,
658     info->getLocalIndexManager());
659 chrisfen 417
660 gezelter 403 info->addMolecule(mol);
661 chrisfen 417
662 gezelter 403 #ifdef IS_MPI
663 chrisfen 417
664 gezelter 403 }
665 chrisfen 417
666 gezelter 403 #endif
667 chrisfen 417
668 gezelter 403 } //end for(int i=0)
669     }
670 chrisfen 417
671 gezelter 1782 int SimCreator::computeStorageLayout(SimInfo* info) {
672    
673     Globals* simParams = info->getSimParams();
674     int nRigidBodies = info->getNGlobalRigidBodies();
675     set<AtomType*> atomTypes = info->getSimulatedAtomTypes();
676     set<AtomType*>::iterator i;
677     bool hasDirectionalAtoms = false;
678     bool hasFixedCharge = false;
679     bool hasMultipoles = false;
680     bool hasPolarizable = false;
681     bool hasFluctuatingCharge = false;
682     bool hasMetallic = false;
683     int storageLayout = 0;
684     storageLayout |= DataStorage::dslPosition;
685     storageLayout |= DataStorage::dslVelocity;
686     storageLayout |= DataStorage::dslForce;
687    
688     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
689    
690     DirectionalAdapter da = DirectionalAdapter( (*i) );
691     MultipoleAdapter ma = MultipoleAdapter( (*i) );
692     EAMAdapter ea = EAMAdapter( (*i) );
693     SuttonChenAdapter sca = SuttonChenAdapter( (*i) );
694     PolarizableAdapter pa = PolarizableAdapter( (*i) );
695     FixedChargeAdapter fca = FixedChargeAdapter( (*i) );
696     FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) );
697    
698     if (da.isDirectional()){
699     hasDirectionalAtoms = true;
700     }
701     if (ma.isMultipole()){
702     hasMultipoles = true;
703     }
704     if (ea.isEAM() || sca.isSuttonChen()){
705     hasMetallic = true;
706     }
707     if ( fca.isFixedCharge() ){
708     hasFixedCharge = true;
709     }
710     if ( fqa.isFluctuatingCharge() ){
711     hasFluctuatingCharge = true;
712     }
713     if ( pa.isPolarizable() ){
714     hasPolarizable = true;
715     }
716     }
717    
718     if (nRigidBodies > 0 || hasDirectionalAtoms) {
719     storageLayout |= DataStorage::dslAmat;
720     if(storageLayout & DataStorage::dslVelocity) {
721     storageLayout |= DataStorage::dslAngularMomentum;
722     }
723     if (storageLayout & DataStorage::dslForce) {
724     storageLayout |= DataStorage::dslTorque;
725     }
726     }
727     if (hasMultipoles) {
728     storageLayout |= DataStorage::dslElectroFrame;
729     }
730     if (hasFixedCharge || hasFluctuatingCharge) {
731     storageLayout |= DataStorage::dslSkippedCharge;
732     }
733     if (hasMetallic) {
734     storageLayout |= DataStorage::dslDensity;
735     storageLayout |= DataStorage::dslFunctional;
736     storageLayout |= DataStorage::dslFunctionalDerivative;
737     }
738     if (hasPolarizable) {
739     storageLayout |= DataStorage::dslElectricField;
740     }
741     if (hasFluctuatingCharge){
742     storageLayout |= DataStorage::dslFlucQPosition;
743     if(storageLayout & DataStorage::dslVelocity) {
744     storageLayout |= DataStorage::dslFlucQVelocity;
745     }
746     if (storageLayout & DataStorage::dslForce) {
747     storageLayout |= DataStorage::dslFlucQForce;
748     }
749     }
750    
751     // if the user has asked for them, make sure we've got the memory for the
752     // objects defined.
753    
754     if (simParams->getOutputParticlePotential()) {
755     storageLayout |= DataStorage::dslParticlePot;
756     }
757    
758     if (simParams->havePrintHeatFlux()) {
759     if (simParams->getPrintHeatFlux()) {
760     storageLayout |= DataStorage::dslParticlePot;
761     }
762     }
763    
764     if (simParams->getOutputElectricField()) {
765     storageLayout |= DataStorage::dslElectricField;
766     }
767     if (simParams->getOutputFluctuatingCharges()) {
768     storageLayout |= DataStorage::dslFlucQPosition;
769     storageLayout |= DataStorage::dslFlucQVelocity;
770     storageLayout |= DataStorage::dslFlucQForce;
771     }
772    
773     return storageLayout;
774     }
775    
776 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
777     SimInfo::MoleculeIterator mi;
778     Molecule::AtomIterator ai;
779     Molecule::RigidBodyIterator ri;
780     Molecule::CutoffGroupIterator ci;
781 tim 1024 Molecule::IntegrableObjectIterator ioi;
782 gezelter 403 Molecule * mol;
783     Atom * atom;
784     RigidBody * rb;
785     CutoffGroup * cg;
786     int beginAtomIndex;
787     int beginRigidBodyIndex;
788     int beginCutoffGroupIndex;
789     int nGlobalAtoms = info->getNGlobalAtoms();
790 gezelter 1802 int nGlobalRigidBodies = info->getNGlobalRigidBodies();
791 chrisfen 417
792 gezelter 403 beginAtomIndex = 0;
793 gezelter 1782 //rigidbody's index begins right after atom's
794     beginRigidBodyIndex = info->getNGlobalAtoms();
795 gezelter 403 beginCutoffGroupIndex = 0;
796 gezelter 1782
797     for(int i = 0; i < info->getNGlobalMolecules(); i++) {
798 chrisfen 417
799 gezelter 1782 #ifdef IS_MPI
800     if (info->getMolToProc(i) == worldRank) {
801     #endif
802     // stuff to do if I own this molecule
803     mol = info->getMoleculeByGlobalIndex(i);
804    
805     //local index(index in DataStorge) of atom is important
806     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
807     atom->setGlobalIndex(beginAtomIndex++);
808     }
809    
810     for(rb = mol->beginRigidBody(ri); rb != NULL;
811     rb = mol->nextRigidBody(ri)) {
812     rb->setGlobalIndex(beginRigidBodyIndex++);
813     }
814    
815     //local index of cutoff group is trivial, it only depends on
816     //the order of travesing
817     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
818     cg = mol->nextCutoffGroup(ci)) {
819     cg->setGlobalIndex(beginCutoffGroupIndex++);
820     }
821    
822     #ifdef IS_MPI
823     } else {
824    
825     // stuff to do if I don't own this molecule
826    
827     int stampId = info->getMoleculeStampId(i);
828     MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
829    
830     beginAtomIndex += stamp->getNAtoms();
831     beginRigidBodyIndex += stamp->getNRigidBodies();
832     beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
833 gezelter 403 }
834 gezelter 1782 #endif
835    
836     } //end for(int i=0)
837    
838 gezelter 403 //fill globalGroupMembership
839     std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
840     for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
841     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
842 chrisfen 417
843 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
844     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
845     }
846 chrisfen 417
847 gezelter 403 }
848     }
849 gezelter 1782
850 gezelter 403 #ifdef IS_MPI
851     // Since the globalGroupMembership has been zero filled and we've only
852     // poked values into the atoms we know, we can do an Allreduce
853     // to get the full globalGroupMembership array (We think).
854     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
855     // docs said we could.
856 gezelter 1313 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
857 gezelter 1796 MPI::COMM_WORLD.Allreduce(&globalGroupMembership[0],
858     &tmpGroupMembership[0], nGlobalAtoms,
859     MPI::INT, MPI::SUM);
860 gezelter 403 info->setGlobalGroupMembership(tmpGroupMembership);
861     #else
862     info->setGlobalGroupMembership(globalGroupMembership);
863     #endif
864 chrisfen 417
865 gezelter 403 //fill molMembership
866 gezelter 1802 std::vector<int> globalMolMembership(info->getNGlobalAtoms() +
867     info->getNGlobalRigidBodies(), 0);
868 gezelter 403
869 gezelter 1802 for(mol = info->beginMolecule(mi); mol != NULL;
870     mol = info->nextMolecule(mi)) {
871 gezelter 403 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
872     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
873     }
874 gezelter 1802 for (rb = mol->beginRigidBody(ri); rb != NULL;
875     rb = mol->nextRigidBody(ri)) {
876     globalMolMembership[rb->getGlobalIndex()] = mol->getGlobalIndex();
877     }
878 gezelter 403 }
879 chrisfen 417
880 gezelter 403 #ifdef IS_MPI
881 gezelter 1802 std::vector<int> tmpMolMembership(info->getNGlobalAtoms() +
882     info->getNGlobalRigidBodies(), 0);
883 gezelter 1796 MPI::COMM_WORLD.Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
884 gezelter 1802 nGlobalAtoms + nGlobalRigidBodies,
885 gezelter 1796 MPI::INT, MPI::SUM);
886 chrisfen 417
887 gezelter 403 info->setGlobalMolMembership(tmpMolMembership);
888     #else
889     info->setGlobalMolMembership(globalMolMembership);
890     #endif
891 tim 1024
892     // nIOPerMol holds the number of integrable objects per molecule
893     // here the molecules are listed by their global indices.
894    
895     std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
896 gezelter 1802 for (mol = info->beginMolecule(mi); mol != NULL;
897     mol = info->nextMolecule(mi)) {
898 tim 1024 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
899     }
900 chrisfen 417
901 tim 1024 #ifdef IS_MPI
902     std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
903 gezelter 1796 MPI::COMM_WORLD.Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
904     info->getNGlobalMolecules(), MPI::INT, MPI::SUM);
905 tim 1024 #else
906     std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
907     #endif
908    
909 gezelter 1313 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
910    
911     int startingIndex = 0;
912     for (int i = 0; i < info->getNGlobalMolecules(); i++) {
913     startingIOIndexForMol[i] = startingIndex;
914     startingIndex += numIntegrableObjectsPerMol[i];
915     }
916    
917     std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
918 gezelter 1802 for (mol = info->beginMolecule(mi); mol != NULL;
919     mol = info->nextMolecule(mi)) {
920 tim 1024 int myGlobalIndex = mol->getGlobalIndex();
921     int globalIO = startingIOIndexForMol[myGlobalIndex];
922 gezelter 1782 for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
923     sd = mol->nextIntegrableObject(ioi)) {
924     sd->setGlobalIntegrableObjectIndex(globalIO);
925     IOIndexToIntegrableObject[globalIO] = sd;
926 gezelter 1313 globalIO++;
927 tim 1024 }
928     }
929 gezelter 1782
930 gezelter 1313 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
931    
932 gezelter 403 }
933 chrisfen 417
934 tim 1024 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
935 gezelter 1782
936 tim 1024 DumpReader reader(info, mdFileName);
937 gezelter 403 int nframes = reader.getNFrames();
938 gezelter 1782
939 gezelter 403 if (nframes > 0) {
940     reader.readFrame(nframes - 1);
941     } else {
942     //invalid initial coordinate file
943 chrisfen 417 sprintf(painCave.errMsg,
944     "Initial configuration file %s should at least contain one frame\n",
945 tim 1024 mdFileName.c_str());
946 gezelter 403 painCave.isFatal = 1;
947     simError();
948     }
949     //copy the current snapshot to previous snapshot
950     info->getSnapshotManager()->advance();
951     }
952 chrisfen 417
953 gezelter 1390 } //end namespace OpenMD
954 gezelter 403
955    

Properties

Name Value
svn:keywords Author Id Revision Date