ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1808
Committed: Mon Oct 22 20:42:10 2012 UTC (12 years, 6 months ago) by gezelter
File size: 31893 byte(s)
Log Message:
A bug fix in the electric field for the new electrostatic code.  Also comment fixes for Doxygen 

File Contents

# User Rev Content
1 gezelter 403 /*
2     * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved.
3     *
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 1665 * [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     * @version 1.0
48     */
49 tim 845 #include <exception>
50 tim 770 #include <iostream>
51     #include <sstream>
52     #include <string>
53    
54 gezelter 403 #include "brains/MoleculeCreator.hpp"
55     #include "brains/SimCreator.hpp"
56     #include "brains/SimSnapshotManager.hpp"
57     #include "io/DumpReader.hpp"
58 gezelter 1725 #include "brains/ForceField.hpp"
59 gezelter 403 #include "utils/simError.h"
60     #include "utils/StringUtils.hpp"
61     #include "math/SeqRandNumGen.hpp"
62 tim 770 #include "mdParser/MDLexer.hpp"
63     #include "mdParser/MDParser.hpp"
64     #include "mdParser/MDTreeParser.hpp"
65     #include "mdParser/SimplePreprocessor.hpp"
66 tim 816 #include "antlr/ANTLRException.hpp"
67     #include "antlr/TokenStreamRecognitionException.hpp"
68     #include "antlr/TokenStreamIOException.hpp"
69     #include "antlr/TokenStreamException.hpp"
70     #include "antlr/RecognitionException.hpp"
71     #include "antlr/CharStreamException.hpp"
72 tim 770
73 tim 816 #include "antlr/MismatchedCharException.hpp"
74     #include "antlr/MismatchedTokenException.hpp"
75     #include "antlr/NoViableAltForCharException.hpp"
76     #include "antlr/NoViableAltException.hpp"
77 tim 770
78 gezelter 1710 #include "types/DirectionalAdapter.hpp"
79     #include "types/MultipoleAdapter.hpp"
80     #include "types/EAMAdapter.hpp"
81     #include "types/SuttonChenAdapter.hpp"
82     #include "types/PolarizableAdapter.hpp"
83     #include "types/FixedChargeAdapter.hpp"
84     #include "types/FluctuatingChargeAdapter.hpp"
85    
86 gezelter 403 #ifdef IS_MPI
87 gezelter 1627 #include "mpi.h"
88 gezelter 403 #include "math/ParallelRandNumGen.hpp"
89     #endif
90    
91 gezelter 1390 namespace OpenMD {
92 chrisfen 417
93 gezelter 1613 Globals* SimCreator::parseFile(std::istream& rawMetaDataStream, const std::string& filename, int mdFileVersion, int startOfMetaDataBlock ){
94 tim 1024 Globals* simParams = NULL;
95     try {
96 tim 770
97 tim 1024 // Create a preprocessor that preprocesses md file into an ostringstream
98     std::stringstream ppStream;
99 tim 770 #ifdef IS_MPI
100 tim 1024 int streamSize;
101     const int masterNode = 0;
102     int commStatus;
103     if (worldRank == masterNode) {
104 gezelter 1613 commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
105     #endif
106 tim 1024 SimplePreprocessor preprocessor;
107     preprocessor.preprocess(rawMetaDataStream, filename, startOfMetaDataBlock, ppStream);
108 tim 770
109     #ifdef IS_MPI
110 tim 1024 //brocasting the stream size
111     streamSize = ppStream.str().size() +1;
112     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
113 tim 770
114 tim 1024 commStatus = MPI_Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())), streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
115 tim 770
116    
117 tim 1024 } else {
118 gezelter 1613
119     commStatus = MPI_Bcast(&mdFileVersion, 1, MPI_INT, masterNode, MPI_COMM_WORLD);
120    
121 tim 1024 //get stream size
122     commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD);
123 gezelter 1313
124 tim 1024 char* buf = new char[streamSize];
125     assert(buf);
126 tim 770
127 tim 1024 //receive file content
128     commStatus = MPI_Bcast(buf, streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD);
129 tim 770
130 tim 1024 ppStream.str(buf);
131 gezelter 1313 delete [] buf;
132 tim 770
133 tim 1024 }
134 tim 770 #endif
135 tim 1024 // Create a scanner that reads from the input stream
136     MDLexer lexer(ppStream);
137     lexer.setFilename(filename);
138     lexer.initDeferredLineCount();
139 chrisfen 417
140 tim 1024 // Create a parser that reads from the scanner
141     MDParser parser(lexer);
142     parser.setFilename(filename);
143 tim 770
144 tim 1024 // Create an observer that synchorizes file name change
145     FilenameObserver observer;
146     observer.setLexer(&lexer);
147     observer.setParser(&parser);
148     lexer.setObserver(&observer);
149 chrisfen 417
150 tim 1024 antlr::ASTFactory factory;
151     parser.initializeASTFactory(factory);
152     parser.setASTFactory(&factory);
153     parser.mdfile();
154 tim 770
155 tim 1024 // Create a tree parser that reads information into Globals
156     MDTreeParser treeParser;
157     treeParser.initializeASTFactory(factory);
158     treeParser.setASTFactory(&factory);
159     simParams = treeParser.walkTree(parser.getAST());
160     }
161 tim 845
162 tim 816
163 tim 1024 catch(antlr::MismatchedCharException& e) {
164     sprintf(painCave.errMsg,
165     "parser exception: %s %s:%d:%d\n",
166     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
167     painCave.isFatal = 1;
168     simError();
169     }
170     catch(antlr::MismatchedTokenException &e) {
171     sprintf(painCave.errMsg,
172     "parser exception: %s %s:%d:%d\n",
173     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
174     painCave.isFatal = 1;
175     simError();
176     }
177     catch(antlr::NoViableAltForCharException &e) {
178     sprintf(painCave.errMsg,
179     "parser exception: %s %s:%d:%d\n",
180     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
181     painCave.isFatal = 1;
182     simError();
183     }
184     catch(antlr::NoViableAltException &e) {
185     sprintf(painCave.errMsg,
186     "parser exception: %s %s:%d:%d\n",
187     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
188     painCave.isFatal = 1;
189     simError();
190     }
191 tim 845
192 tim 1024 catch(antlr::TokenStreamRecognitionException& e) {
193     sprintf(painCave.errMsg,
194     "parser exception: %s %s:%d:%d\n",
195     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
196     painCave.isFatal = 1;
197     simError();
198     }
199 tim 845
200 tim 1024 catch(antlr::TokenStreamIOException& e) {
201     sprintf(painCave.errMsg,
202     "parser exception: %s\n",
203     e.getMessage().c_str());
204     painCave.isFatal = 1;
205     simError();
206     }
207 tim 845
208 tim 1024 catch(antlr::TokenStreamException& e) {
209     sprintf(painCave.errMsg,
210     "parser exception: %s\n",
211     e.getMessage().c_str());
212     painCave.isFatal = 1;
213     simError();
214     }
215     catch (antlr::RecognitionException& e) {
216     sprintf(painCave.errMsg,
217     "parser exception: %s %s:%d:%d\n",
218     e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn());
219     painCave.isFatal = 1;
220     simError();
221     }
222     catch (antlr::CharStreamException& e) {
223     sprintf(painCave.errMsg,
224     "parser exception: %s\n",
225     e.getMessage().c_str());
226     painCave.isFatal = 1;
227     simError();
228     }
229 gezelter 1390 catch (OpenMDException& e) {
230 tim 1024 sprintf(painCave.errMsg,
231     "%s\n",
232     e.getMessage().c_str());
233     painCave.isFatal = 1;
234     simError();
235     }
236     catch (std::exception& e) {
237     sprintf(painCave.errMsg,
238     "parser exception: %s\n",
239     e.what());
240     painCave.isFatal = 1;
241     simError();
242     }
243 tim 770
244 gezelter 1613 simParams->setMDfileVersion(mdFileVersion);
245 tim 1024 return simParams;
246 gezelter 403 }
247 chrisfen 417
248 chrisfen 514 SimInfo* SimCreator::createSim(const std::string & mdFileName,
249     bool loadInitCoords) {
250 gezelter 1469
251 tim 1024 const int bufferSize = 65535;
252     char buffer[bufferSize];
253     int lineNo = 0;
254     std::string mdRawData;
255     int metaDataBlockStart = -1;
256     int metaDataBlockEnd = -1;
257     int i;
258 gezelter 1767 streamoff mdOffset;
259 gezelter 1613 int mdFileVersion;
260 tim 1024
261 gezelter 1794
262 tim 1024 #ifdef IS_MPI
263     const int masterNode = 0;
264     if (worldRank == masterNode) {
265     #endif
266    
267 gezelter 1794 std::ifstream mdFile_;
268     mdFile_.open(mdFileName.c_str(), ifstream::in | ifstream::binary);
269 tim 1024
270     if (mdFile_.fail()) {
271     sprintf(painCave.errMsg,
272     "SimCreator: Cannot open file: %s\n",
273     mdFileName.c_str());
274     painCave.isFatal = 1;
275     simError();
276     }
277    
278     mdFile_.getline(buffer, bufferSize);
279     ++lineNo;
280     std::string line = trimLeftCopy(buffer);
281 gezelter 1390 i = CaseInsensitiveFind(line, "<OpenMD");
282 gezelter 1287 if (static_cast<size_t>(i) == string::npos) {
283 gezelter 1390 // try the older file strings to see if that works:
284     i = CaseInsensitiveFind(line, "<OOPSE");
285     }
286    
287     if (static_cast<size_t>(i) == string::npos) {
288     // still no luck!
289 tim 1024 sprintf(painCave.errMsg,
290 gezelter 1390 "SimCreator: File: %s is not a valid OpenMD file!\n",
291 tim 1024 mdFileName.c_str());
292     painCave.isFatal = 1;
293     simError();
294     }
295 gezelter 1613
296     // found the correct opening string, now try to get the file
297     // format version number.
298 tim 1024
299 gezelter 1613 StringTokenizer tokenizer(line, "=<> \t\n\r");
300     std::string fileType = tokenizer.nextToken();
301     toUpper(fileType);
302    
303     mdFileVersion = 0;
304    
305     if (fileType == "OPENMD") {
306     while (tokenizer.hasMoreTokens()) {
307     std::string token(tokenizer.nextToken());
308     toUpper(token);
309     if (token == "VERSION") {
310     mdFileVersion = tokenizer.nextTokenAsInt();
311     break;
312     }
313     }
314     }
315    
316 tim 1024 //scan through the input stream and find MetaData tag
317     while(mdFile_.getline(buffer, bufferSize)) {
318     ++lineNo;
319    
320     std::string line = trimLeftCopy(buffer);
321     if (metaDataBlockStart == -1) {
322     i = CaseInsensitiveFind(line, "<MetaData>");
323     if (i != string::npos) {
324     metaDataBlockStart = lineNo;
325     mdOffset = mdFile_.tellg();
326     }
327     } else {
328     i = CaseInsensitiveFind(line, "</MetaData>");
329     if (i != string::npos) {
330     metaDataBlockEnd = lineNo;
331     }
332     }
333     }
334    
335     if (metaDataBlockStart == -1) {
336     sprintf(painCave.errMsg,
337     "SimCreator: File: %s did not contain a <MetaData> tag!\n",
338     mdFileName.c_str());
339     painCave.isFatal = 1;
340     simError();
341     }
342     if (metaDataBlockEnd == -1) {
343     sprintf(painCave.errMsg,
344     "SimCreator: File: %s did not contain a closed MetaData block!\n",
345     mdFileName.c_str());
346     painCave.isFatal = 1;
347     simError();
348     }
349    
350     mdFile_.clear();
351     mdFile_.seekg(0);
352     mdFile_.seekg(mdOffset);
353    
354     mdRawData.clear();
355    
356     for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
357     mdFile_.getline(buffer, bufferSize);
358     mdRawData += buffer;
359     mdRawData += "\n";
360     }
361    
362     mdFile_.close();
363    
364     #ifdef IS_MPI
365     }
366     #endif
367    
368     std::stringstream rawMetaDataStream(mdRawData);
369    
370 gezelter 403 //parse meta-data file
371 gezelter 1613 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
372     metaDataBlockStart + 1);
373 chrisfen 417
374 gezelter 403 //create the force field
375 gezelter 1725 ForceField * ff = new ForceField(simParams->getForceField());
376 gezelter 1277
377 gezelter 403 if (ff == NULL) {
378 chrisfen 514 sprintf(painCave.errMsg,
379     "ForceField Factory can not create %s force field\n",
380 tim 665 simParams->getForceField().c_str());
381 gezelter 403 painCave.isFatal = 1;
382     simError();
383     }
384 chrisfen 417
385 gezelter 403 if (simParams->haveForceFieldFileName()) {
386     ff->setForceFieldFileName(simParams->getForceFieldFileName());
387     }
388    
389     std::string forcefieldFileName;
390     forcefieldFileName = ff->getForceFieldFileName();
391 chrisfen 417
392 gezelter 403 if (simParams->haveForceFieldVariant()) {
393     //If the force field has variant, the variant force field name will be
394     //Base.variant.frc. For exampel EAM.u6.frc
395 chrisfen 417
396 gezelter 403 std::string variant = simParams->getForceFieldVariant();
397 chrisfen 417
398 gezelter 403 std::string::size_type pos = forcefieldFileName.rfind(".frc");
399     variant = "." + variant;
400     if (pos != std::string::npos) {
401     forcefieldFileName.insert(pos, variant);
402     } else {
403     //If the default force field file name does not containt .frc suffix, just append the .variant
404     forcefieldFileName.append(variant);
405     }
406     }
407    
408     ff->parse(forcefieldFileName);
409     //create SimInfo
410 tim 770 SimInfo * info = new SimInfo(ff, simParams);
411 tim 1024
412     info->setRawMetaData(mdRawData);
413 tim 490
414 gezelter 945 //gather parameters (SimCreator only retrieves part of the
415     //parameters)
416 gezelter 403 gatherParameters(info, mdFileName);
417 chrisfen 417
418 gezelter 403 //divide the molecules and determine the global index of molecules
419     #ifdef IS_MPI
420     divideMolecules(info);
421     #endif
422 chrisfen 417
423 gezelter 403 //create the molecules
424     createMolecules(info);
425 chrisfen 417
426 gezelter 1710 //find the storage layout
427    
428     int storageLayout = computeStorageLayout(info);
429    
430 gezelter 945 //allocate memory for DataStorage(circular reference, need to
431     //break it)
432 gezelter 1710 info->setSnapshotManager(new SimSnapshotManager(info, storageLayout));
433 gezelter 403
434 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
435     //(only need to be set once, the global index will never change
436     //again). Local indices of atoms and rigidbodies are already set
437     //by MoleculeCreator class which actually delegates the
438     //responsibility to LocalIndexManager.
439 gezelter 403 setGlobalIndex(info);
440 chrisfen 417
441 gezelter 1287 //Although addInteractionPairs is called inside SimInfo's addMolecule
442 gezelter 945 //method, at that point atoms don't have the global index yet
443     //(their global index are all initialized to -1). Therefore we
444 gezelter 1287 //have to call addInteractionPairs explicitly here. A way to work
445 gezelter 945 //around is that we can determine the beginning global indices of
446     //atoms before they get created.
447 gezelter 403 SimInfo::MoleculeIterator mi;
448     Molecule* mol;
449     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
450 gezelter 1287 info->addInteractionPairs(mol);
451 gezelter 403 }
452    
453     if (loadInitCoords)
454 tim 1024 loadCoordinates(info, mdFileName);
455 gezelter 403 return info;
456     }
457 chrisfen 417
458 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
459 chrisfen 417
460 tim 749 //figure out the output file names
461 gezelter 403 std::string prefix;
462 chrisfen 417
463 gezelter 403 #ifdef IS_MPI
464 chrisfen 417
465 gezelter 403 if (worldRank == 0) {
466     #endif // is_mpi
467     Globals * simParams = info->getSimParams();
468     if (simParams->haveFinalConfig()) {
469     prefix = getPrefix(simParams->getFinalConfig());
470     } else {
471     prefix = getPrefix(mdfile);
472     }
473 chrisfen 417
474 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
475     info->setDumpFileName(prefix + ".dump");
476     info->setStatFileName(prefix + ".stat");
477 chrisfen 417 info->setRestFileName(prefix + ".zang");
478    
479 gezelter 403 #ifdef IS_MPI
480 chrisfen 417
481 gezelter 403 }
482 chrisfen 417
483 gezelter 403 #endif
484 chrisfen 417
485 gezelter 403 }
486 chrisfen 417
487 gezelter 403 #ifdef IS_MPI
488     void SimCreator::divideMolecules(SimInfo *info) {
489 tim 963 RealType numerator;
490     RealType denominator;
491     RealType precast;
492     RealType x;
493     RealType y;
494     RealType a;
495 gezelter 403 int old_atoms;
496     int add_atoms;
497     int new_atoms;
498     int nTarget;
499     int done;
500     int i;
501     int j;
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 1798 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 1803
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 1803
573 gezelter 403 sprintf(painCave.errMsg,
574 gezelter 1803 "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 1803
579 gezelter 403 painCave.isFatal = 0;
580 gezelter 1803 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 1803
626 gezelter 403 // Spray out this nonsense to all other processors:
627 gezelter 1798 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 1798 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
632 gezelter 1803
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 1600 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 1710 int SimCreator::computeStorageLayout(SimInfo* info) {
672    
673 gezelter 1711 Globals* simParams = info->getSimParams();
674 gezelter 1710 int nRigidBodies = info->getNGlobalRigidBodies();
675     set<AtomType*> atomTypes = info->getSimulatedAtomTypes();
676     set<AtomType*>::iterator i;
677     bool hasDirectionalAtoms = false;
678     bool hasFixedCharge = false;
679 gezelter 1794 bool hasDipoles = false;
680     bool hasQuadrupoles = false;
681 gezelter 1710 bool hasPolarizable = false;
682     bool hasFluctuatingCharge = false;
683     bool hasMetallic = false;
684     int storageLayout = 0;
685     storageLayout |= DataStorage::dslPosition;
686     storageLayout |= DataStorage::dslVelocity;
687     storageLayout |= DataStorage::dslForce;
688    
689     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
690    
691     DirectionalAdapter da = DirectionalAdapter( (*i) );
692     MultipoleAdapter ma = MultipoleAdapter( (*i) );
693     EAMAdapter ea = EAMAdapter( (*i) );
694     SuttonChenAdapter sca = SuttonChenAdapter( (*i) );
695     PolarizableAdapter pa = PolarizableAdapter( (*i) );
696     FixedChargeAdapter fca = FixedChargeAdapter( (*i) );
697     FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) );
698    
699     if (da.isDirectional()){
700     hasDirectionalAtoms = true;
701     }
702 gezelter 1794 if (ma.isDipole()){
703     hasDipoles = true;
704 gezelter 1710 }
705 gezelter 1794 if (ma.isQuadrupole()){
706     hasQuadrupoles = true;
707     }
708 gezelter 1710 if (ea.isEAM() || sca.isSuttonChen()){
709     hasMetallic = true;
710     }
711     if ( fca.isFixedCharge() ){
712     hasFixedCharge = true;
713     }
714     if ( fqa.isFluctuatingCharge() ){
715     hasFluctuatingCharge = true;
716     }
717     if ( pa.isPolarizable() ){
718     hasPolarizable = true;
719     }
720     }
721    
722     if (nRigidBodies > 0 || hasDirectionalAtoms) {
723     storageLayout |= DataStorage::dslAmat;
724     if(storageLayout & DataStorage::dslVelocity) {
725     storageLayout |= DataStorage::dslAngularMomentum;
726     }
727     if (storageLayout & DataStorage::dslForce) {
728     storageLayout |= DataStorage::dslTorque;
729     }
730     }
731 gezelter 1794 if (hasDipoles) {
732     storageLayout |= DataStorage::dslDipole;
733     }
734     if (hasQuadrupoles) {
735     storageLayout |= DataStorage::dslQuadrupole;
736     }
737 gezelter 1710 if (hasFixedCharge || hasFluctuatingCharge) {
738     storageLayout |= DataStorage::dslSkippedCharge;
739     }
740     if (hasMetallic) {
741     storageLayout |= DataStorage::dslDensity;
742     storageLayout |= DataStorage::dslFunctional;
743     storageLayout |= DataStorage::dslFunctionalDerivative;
744     }
745     if (hasPolarizable) {
746     storageLayout |= DataStorage::dslElectricField;
747     }
748     if (hasFluctuatingCharge){
749     storageLayout |= DataStorage::dslFlucQPosition;
750     if(storageLayout & DataStorage::dslVelocity) {
751     storageLayout |= DataStorage::dslFlucQVelocity;
752     }
753     if (storageLayout & DataStorage::dslForce) {
754     storageLayout |= DataStorage::dslFlucQForce;
755     }
756     }
757 gezelter 1714
758     // if the user has asked for them, make sure we've got the memory for the
759     // objects defined.
760 gezelter 1711
761     if (simParams->getOutputParticlePotential()) {
762     storageLayout |= DataStorage::dslParticlePot;
763     }
764 gezelter 1723
765     if (simParams->havePrintHeatFlux()) {
766     if (simParams->getPrintHeatFlux()) {
767     storageLayout |= DataStorage::dslParticlePot;
768     }
769     }
770    
771 gezelter 1714 if (simParams->getOutputElectricField()) {
772     storageLayout |= DataStorage::dslElectricField;
773     }
774 gezelter 1787
775 gezelter 1714 if (simParams->getOutputFluctuatingCharges()) {
776     storageLayout |= DataStorage::dslFlucQPosition;
777     storageLayout |= DataStorage::dslFlucQVelocity;
778     storageLayout |= DataStorage::dslFlucQForce;
779     }
780 gezelter 1711
781 gezelter 1710 return storageLayout;
782     }
783    
784 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
785     SimInfo::MoleculeIterator mi;
786     Molecule::AtomIterator ai;
787     Molecule::RigidBodyIterator ri;
788     Molecule::CutoffGroupIterator ci;
789 tim 1024 Molecule::IntegrableObjectIterator ioi;
790 gezelter 403 Molecule * mol;
791     Atom * atom;
792     RigidBody * rb;
793     CutoffGroup * cg;
794     int beginAtomIndex;
795     int beginRigidBodyIndex;
796     int beginCutoffGroupIndex;
797     int nGlobalAtoms = info->getNGlobalAtoms();
798 gezelter 1803 int nGlobalRigidBodies = info->getNGlobalRigidBodies();
799 chrisfen 417
800 gezelter 403 beginAtomIndex = 0;
801 gezelter 1757 //rigidbody's index begins right after atom's
802     beginRigidBodyIndex = info->getNGlobalAtoms();
803 gezelter 403 beginCutoffGroupIndex = 0;
804 gezelter 1600
805     for(int i = 0; i < info->getNGlobalMolecules(); i++) {
806 chrisfen 417
807 gezelter 1600 #ifdef IS_MPI
808     if (info->getMolToProc(i) == worldRank) {
809     #endif
810     // stuff to do if I own this molecule
811     mol = info->getMoleculeByGlobalIndex(i);
812    
813     //local index(index in DataStorge) of atom is important
814     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
815     atom->setGlobalIndex(beginAtomIndex++);
816     }
817    
818     for(rb = mol->beginRigidBody(ri); rb != NULL;
819     rb = mol->nextRigidBody(ri)) {
820     rb->setGlobalIndex(beginRigidBodyIndex++);
821     }
822    
823     //local index of cutoff group is trivial, it only depends on
824     //the order of travesing
825     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
826     cg = mol->nextCutoffGroup(ci)) {
827     cg->setGlobalIndex(beginCutoffGroupIndex++);
828     }
829    
830     #ifdef IS_MPI
831     } else {
832    
833     // stuff to do if I don't own this molecule
834    
835     int stampId = info->getMoleculeStampId(i);
836     MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
837    
838     beginAtomIndex += stamp->getNAtoms();
839     beginRigidBodyIndex += stamp->getNRigidBodies();
840     beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
841 gezelter 403 }
842 gezelter 1600 #endif
843    
844     } //end for(int i=0)
845    
846 gezelter 403 //fill globalGroupMembership
847 gezelter 1600 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
848 gezelter 403 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
849     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
850 chrisfen 417
851 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
852     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
853     }
854 chrisfen 417
855 gezelter 403 }
856     }
857 gezelter 1577
858 gezelter 403 #ifdef IS_MPI
859     // Since the globalGroupMembership has been zero filled and we've only
860     // poked values into the atoms we know, we can do an Allreduce
861     // to get the full globalGroupMembership array (We think).
862     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
863     // docs said we could.
864 gezelter 1313 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
865 gezelter 1798 MPI::COMM_WORLD.Allreduce(&globalGroupMembership[0],
866     &tmpGroupMembership[0], nGlobalAtoms,
867     MPI::INT, MPI::SUM);
868 gezelter 403 info->setGlobalGroupMembership(tmpGroupMembership);
869     #else
870     info->setGlobalGroupMembership(globalGroupMembership);
871     #endif
872 chrisfen 417
873 gezelter 403 //fill molMembership
874 gezelter 1803 std::vector<int> globalMolMembership(info->getNGlobalAtoms() +
875     info->getNGlobalRigidBodies(), 0);
876 gezelter 403
877 gezelter 1803 for(mol = info->beginMolecule(mi); mol != NULL;
878     mol = info->nextMolecule(mi)) {
879 gezelter 403 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
880     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
881     }
882 gezelter 1803 for (rb = mol->beginRigidBody(ri); rb != NULL;
883     rb = mol->nextRigidBody(ri)) {
884     globalMolMembership[rb->getGlobalIndex()] = mol->getGlobalIndex();
885     }
886 gezelter 403 }
887 chrisfen 417
888 gezelter 403 #ifdef IS_MPI
889 gezelter 1803 std::vector<int> tmpMolMembership(info->getNGlobalAtoms() +
890     info->getNGlobalRigidBodies(), 0);
891 gezelter 1798 MPI::COMM_WORLD.Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
892 gezelter 1803 nGlobalAtoms + nGlobalRigidBodies,
893 gezelter 1798 MPI::INT, MPI::SUM);
894 chrisfen 417
895 gezelter 403 info->setGlobalMolMembership(tmpMolMembership);
896     #else
897     info->setGlobalMolMembership(globalMolMembership);
898     #endif
899 tim 1024
900     // nIOPerMol holds the number of integrable objects per molecule
901     // here the molecules are listed by their global indices.
902    
903     std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
904 gezelter 1803 for (mol = info->beginMolecule(mi); mol != NULL;
905     mol = info->nextMolecule(mi)) {
906 tim 1024 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
907     }
908 chrisfen 417
909 tim 1024 #ifdef IS_MPI
910     std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
911 gezelter 1798 MPI::COMM_WORLD.Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
912     info->getNGlobalMolecules(), MPI::INT, MPI::SUM);
913 tim 1024 #else
914     std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
915     #endif
916    
917 gezelter 1313 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
918    
919     int startingIndex = 0;
920     for (int i = 0; i < info->getNGlobalMolecules(); i++) {
921     startingIOIndexForMol[i] = startingIndex;
922     startingIndex += numIntegrableObjectsPerMol[i];
923     }
924    
925     std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
926 gezelter 1803 for (mol = info->beginMolecule(mi); mol != NULL;
927     mol = info->nextMolecule(mi)) {
928 tim 1024 int myGlobalIndex = mol->getGlobalIndex();
929     int globalIO = startingIOIndexForMol[myGlobalIndex];
930 gezelter 1769 for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
931     sd = mol->nextIntegrableObject(ioi)) {
932     sd->setGlobalIntegrableObjectIndex(globalIO);
933     IOIndexToIntegrableObject[globalIO] = sd;
934 gezelter 1313 globalIO++;
935 tim 1024 }
936     }
937 gezelter 1597
938 gezelter 1313 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
939    
940 gezelter 403 }
941 chrisfen 417
942 tim 1024 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
943 gezelter 403 Globals* simParams;
944 gezelter 1540
945 gezelter 403 simParams = info->getSimParams();
946    
947 tim 1024 DumpReader reader(info, mdFileName);
948 gezelter 403 int nframes = reader.getNFrames();
949 gezelter 1540
950 gezelter 403 if (nframes > 0) {
951     reader.readFrame(nframes - 1);
952     } else {
953     //invalid initial coordinate file
954 chrisfen 417 sprintf(painCave.errMsg,
955     "Initial configuration file %s should at least contain one frame\n",
956 tim 1024 mdFileName.c_str());
957 gezelter 403 painCave.isFatal = 1;
958     simError();
959     }
960     //copy the current snapshot to previous snapshot
961     info->getSnapshotManager()->advance();
962     }
963 chrisfen 417
964 gezelter 1390 } //end namespace OpenMD
965 gezelter 403
966    

Properties

Name Value
svn:keywords Author Id Revision Date