ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1825
Committed: Wed Jan 9 19:27:52 2013 UTC (12 years, 3 months ago) by gezelter
File size: 32810 byte(s)
Log Message:
Deleting unused variables

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 gezelter 1812 int i, j;
258 gezelter 1767 streamoff mdOffset;
259 gezelter 1613 int mdFileVersion;
260 tim 1024
261 gezelter 1812 // Create a string for embedding the version information in the MetaData
262     std::string version;
263     version.assign("## Last run using OpenMD Version: ");
264     version.append(OPENMD_VERSION_MAJOR);
265     version.append(".");
266     version.append(OPENMD_VERSION_MINOR);
267 gezelter 1794
268 gezelter 1812 std::string svnrev;
269     //convert a macro from compiler to a string in c++
270     STR_DEFINE(svnrev, SVN_REV );
271     version.append(" Revision: ");
272     // If there's no SVN revision, just call this the RELEASE revision.
273     if (!svnrev.empty()) {
274     version.append(svnrev);
275     } else {
276     version.append("RELEASE");
277     }
278    
279 tim 1024 #ifdef IS_MPI
280     const int masterNode = 0;
281     if (worldRank == masterNode) {
282     #endif
283    
284 gezelter 1794 std::ifstream mdFile_;
285     mdFile_.open(mdFileName.c_str(), ifstream::in | ifstream::binary);
286 tim 1024
287     if (mdFile_.fail()) {
288     sprintf(painCave.errMsg,
289     "SimCreator: Cannot open file: %s\n",
290     mdFileName.c_str());
291     painCave.isFatal = 1;
292     simError();
293     }
294    
295     mdFile_.getline(buffer, bufferSize);
296     ++lineNo;
297     std::string line = trimLeftCopy(buffer);
298 gezelter 1390 i = CaseInsensitiveFind(line, "<OpenMD");
299 gezelter 1287 if (static_cast<size_t>(i) == string::npos) {
300 gezelter 1390 // try the older file strings to see if that works:
301     i = CaseInsensitiveFind(line, "<OOPSE");
302     }
303    
304     if (static_cast<size_t>(i) == string::npos) {
305     // still no luck!
306 tim 1024 sprintf(painCave.errMsg,
307 gezelter 1390 "SimCreator: File: %s is not a valid OpenMD file!\n",
308 tim 1024 mdFileName.c_str());
309     painCave.isFatal = 1;
310     simError();
311     }
312 gezelter 1613
313     // found the correct opening string, now try to get the file
314     // format version number.
315 tim 1024
316 gezelter 1613 StringTokenizer tokenizer(line, "=<> \t\n\r");
317     std::string fileType = tokenizer.nextToken();
318     toUpper(fileType);
319    
320     mdFileVersion = 0;
321    
322     if (fileType == "OPENMD") {
323     while (tokenizer.hasMoreTokens()) {
324     std::string token(tokenizer.nextToken());
325     toUpper(token);
326     if (token == "VERSION") {
327     mdFileVersion = tokenizer.nextTokenAsInt();
328     break;
329     }
330     }
331     }
332    
333 tim 1024 //scan through the input stream and find MetaData tag
334     while(mdFile_.getline(buffer, bufferSize)) {
335     ++lineNo;
336    
337     std::string line = trimLeftCopy(buffer);
338     if (metaDataBlockStart == -1) {
339     i = CaseInsensitiveFind(line, "<MetaData>");
340     if (i != string::npos) {
341     metaDataBlockStart = lineNo;
342     mdOffset = mdFile_.tellg();
343     }
344     } else {
345     i = CaseInsensitiveFind(line, "</MetaData>");
346     if (i != string::npos) {
347     metaDataBlockEnd = lineNo;
348     }
349     }
350     }
351    
352     if (metaDataBlockStart == -1) {
353     sprintf(painCave.errMsg,
354     "SimCreator: File: %s did not contain a <MetaData> tag!\n",
355     mdFileName.c_str());
356     painCave.isFatal = 1;
357     simError();
358     }
359     if (metaDataBlockEnd == -1) {
360     sprintf(painCave.errMsg,
361     "SimCreator: File: %s did not contain a closed MetaData block!\n",
362     mdFileName.c_str());
363     painCave.isFatal = 1;
364     simError();
365     }
366    
367     mdFile_.clear();
368     mdFile_.seekg(0);
369     mdFile_.seekg(mdOffset);
370    
371     mdRawData.clear();
372    
373 gezelter 1812 bool foundVersion = false;
374    
375 tim 1024 for (int i = 0; i < metaDataBlockEnd - metaDataBlockStart - 1; ++i) {
376     mdFile_.getline(buffer, bufferSize);
377 gezelter 1812 std::string line = trimLeftCopy(buffer);
378     j = CaseInsensitiveFind(line, "## Last run using OpenMD Version");
379     if (static_cast<size_t>(j) != string::npos) {
380     foundVersion = true;
381     mdRawData += version;
382     } else {
383     mdRawData += buffer;
384     }
385 tim 1024 mdRawData += "\n";
386     }
387 gezelter 1812
388     if (!foundVersion) mdRawData += version + "\n";
389    
390 tim 1024 mdFile_.close();
391    
392     #ifdef IS_MPI
393     }
394     #endif
395    
396     std::stringstream rawMetaDataStream(mdRawData);
397    
398 gezelter 403 //parse meta-data file
399 gezelter 1613 Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion,
400     metaDataBlockStart + 1);
401 chrisfen 417
402 gezelter 403 //create the force field
403 gezelter 1725 ForceField * ff = new ForceField(simParams->getForceField());
404 gezelter 1277
405 gezelter 403 if (ff == NULL) {
406 chrisfen 514 sprintf(painCave.errMsg,
407     "ForceField Factory can not create %s force field\n",
408 tim 665 simParams->getForceField().c_str());
409 gezelter 403 painCave.isFatal = 1;
410     simError();
411     }
412 chrisfen 417
413 gezelter 403 if (simParams->haveForceFieldFileName()) {
414     ff->setForceFieldFileName(simParams->getForceFieldFileName());
415     }
416    
417     std::string forcefieldFileName;
418     forcefieldFileName = ff->getForceFieldFileName();
419 chrisfen 417
420 gezelter 403 if (simParams->haveForceFieldVariant()) {
421     //If the force field has variant, the variant force field name will be
422     //Base.variant.frc. For exampel EAM.u6.frc
423 chrisfen 417
424 gezelter 403 std::string variant = simParams->getForceFieldVariant();
425 chrisfen 417
426 gezelter 403 std::string::size_type pos = forcefieldFileName.rfind(".frc");
427     variant = "." + variant;
428     if (pos != std::string::npos) {
429     forcefieldFileName.insert(pos, variant);
430     } else {
431     //If the default force field file name does not containt .frc suffix, just append the .variant
432     forcefieldFileName.append(variant);
433     }
434     }
435    
436     ff->parse(forcefieldFileName);
437     //create SimInfo
438 tim 770 SimInfo * info = new SimInfo(ff, simParams);
439 tim 1024
440     info->setRawMetaData(mdRawData);
441 tim 490
442 gezelter 945 //gather parameters (SimCreator only retrieves part of the
443     //parameters)
444 gezelter 403 gatherParameters(info, mdFileName);
445 chrisfen 417
446 gezelter 403 //divide the molecules and determine the global index of molecules
447     #ifdef IS_MPI
448     divideMolecules(info);
449     #endif
450 chrisfen 417
451 gezelter 403 //create the molecules
452     createMolecules(info);
453 chrisfen 417
454 gezelter 1710 //find the storage layout
455    
456     int storageLayout = computeStorageLayout(info);
457    
458 gezelter 945 //allocate memory for DataStorage(circular reference, need to
459     //break it)
460 gezelter 1710 info->setSnapshotManager(new SimSnapshotManager(info, storageLayout));
461 gezelter 403
462 gezelter 945 //set the global index of atoms, rigidbodies and cutoffgroups
463     //(only need to be set once, the global index will never change
464     //again). Local indices of atoms and rigidbodies are already set
465     //by MoleculeCreator class which actually delegates the
466     //responsibility to LocalIndexManager.
467 gezelter 403 setGlobalIndex(info);
468 chrisfen 417
469 gezelter 1287 //Although addInteractionPairs is called inside SimInfo's addMolecule
470 gezelter 945 //method, at that point atoms don't have the global index yet
471     //(their global index are all initialized to -1). Therefore we
472 gezelter 1287 //have to call addInteractionPairs explicitly here. A way to work
473 gezelter 945 //around is that we can determine the beginning global indices of
474     //atoms before they get created.
475 gezelter 403 SimInfo::MoleculeIterator mi;
476     Molecule* mol;
477     for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
478 gezelter 1287 info->addInteractionPairs(mol);
479 gezelter 403 }
480    
481     if (loadInitCoords)
482 tim 1024 loadCoordinates(info, mdFileName);
483 gezelter 403 return info;
484     }
485 chrisfen 417
486 gezelter 403 void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) {
487 chrisfen 417
488 tim 749 //figure out the output file names
489 gezelter 403 std::string prefix;
490 chrisfen 417
491 gezelter 403 #ifdef IS_MPI
492 chrisfen 417
493 gezelter 403 if (worldRank == 0) {
494     #endif // is_mpi
495     Globals * simParams = info->getSimParams();
496     if (simParams->haveFinalConfig()) {
497     prefix = getPrefix(simParams->getFinalConfig());
498     } else {
499     prefix = getPrefix(mdfile);
500     }
501 chrisfen 417
502 gezelter 403 info->setFinalConfigFileName(prefix + ".eor");
503     info->setDumpFileName(prefix + ".dump");
504     info->setStatFileName(prefix + ".stat");
505 chrisfen 417 info->setRestFileName(prefix + ".zang");
506    
507 gezelter 403 #ifdef IS_MPI
508 chrisfen 417
509 gezelter 403 }
510 chrisfen 417
511 gezelter 403 #endif
512 chrisfen 417
513 gezelter 403 }
514 chrisfen 417
515 gezelter 403 #ifdef IS_MPI
516     void SimCreator::divideMolecules(SimInfo *info) {
517 tim 963 RealType numerator;
518     RealType denominator;
519     RealType precast;
520     RealType x;
521     RealType y;
522     RealType a;
523 gezelter 403 int old_atoms;
524     int add_atoms;
525     int new_atoms;
526     int nTarget;
527     int done;
528     int i;
529     int j;
530     int loops;
531     int which_proc;
532     int nProcessors;
533     std::vector<int> atomsPerProc;
534     int nGlobalMols = info->getNGlobalMolecules();
535     std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition:
536    
537 gezelter 1798 nProcessors = MPI::COMM_WORLD.Get_size();
538 chrisfen 417
539 gezelter 403 if (nProcessors > nGlobalMols) {
540     sprintf(painCave.errMsg,
541     "nProcessors (%d) > nMol (%d)\n"
542     "\tThe number of processors is larger than\n"
543     "\tthe number of molecules. This will not result in a \n"
544     "\tusable division of atoms for force decomposition.\n"
545     "\tEither try a smaller number of processors, or run the\n"
546 gezelter 1390 "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols);
547 chrisfen 417
548 gezelter 403 painCave.isFatal = 1;
549     simError();
550     }
551 chrisfen 417
552 gezelter 403 int seedValue;
553     Globals * simParams = info->getSimParams();
554     SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator
555     if (simParams->haveSeed()) {
556     seedValue = simParams->getSeed();
557     myRandom = new SeqRandNumGen(seedValue);
558     }else {
559     myRandom = new SeqRandNumGen();
560     }
561 chrisfen 417
562    
563 gezelter 403 a = 3.0 * nGlobalMols / info->getNGlobalAtoms();
564 chrisfen 417
565 gezelter 403 //initialize atomsPerProc
566     atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0);
567 chrisfen 417
568 gezelter 403 if (worldRank == 0) {
569     numerator = info->getNGlobalAtoms();
570     denominator = nProcessors;
571     precast = numerator / denominator;
572     nTarget = (int)(precast + 0.5);
573 chrisfen 417
574 gezelter 403 for(i = 0; i < nGlobalMols; i++) {
575 gezelter 1803
576 gezelter 403 done = 0;
577     loops = 0;
578 chrisfen 417
579 gezelter 403 while (!done) {
580     loops++;
581 chrisfen 417
582 gezelter 403 // Pick a processor at random
583 chrisfen 417
584 gezelter 403 which_proc = (int) (myRandom->rand() * nProcessors);
585 chrisfen 417
586 gezelter 403 //get the molecule stamp first
587     int stampId = info->getMoleculeStampId(i);
588     MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId);
589 chrisfen 417
590 gezelter 403 // How many atoms does this processor have so far?
591     old_atoms = atomsPerProc[which_proc];
592     add_atoms = moleculeStamp->getNAtoms();
593     new_atoms = old_atoms + add_atoms;
594 chrisfen 417
595 gezelter 403 // If we've been through this loop too many times, we need
596     // to just give up and assign the molecule to this processor
597     // and be done with it.
598 chrisfen 417
599 gezelter 403 if (loops > 100) {
600 gezelter 1803
601 gezelter 403 sprintf(painCave.errMsg,
602 gezelter 1803 "There have been 100 attempts to assign molecule %d to an\n"
603     "\tunderworked processor, but there's no good place to\n"
604     "\tleave it. OpenMD is assigning it at random to processor %d.\n",
605 gezelter 403 i, which_proc);
606 gezelter 1803
607 gezelter 403 painCave.isFatal = 0;
608 gezelter 1803 painCave.severity = OPENMD_INFO;
609 gezelter 403 simError();
610 chrisfen 417
611 gezelter 403 molToProcMap[i] = which_proc;
612     atomsPerProc[which_proc] += add_atoms;
613 chrisfen 417
614 gezelter 403 done = 1;
615     continue;
616     }
617 chrisfen 417
618 gezelter 403 // If we can add this molecule to this processor without sending
619     // it above nTarget, then go ahead and do it:
620 chrisfen 417
621 gezelter 403 if (new_atoms <= nTarget) {
622     molToProcMap[i] = which_proc;
623     atomsPerProc[which_proc] += add_atoms;
624 chrisfen 417
625 gezelter 403 done = 1;
626     continue;
627     }
628 chrisfen 417
629 gezelter 403 // The only situation left is when new_atoms > nTarget. We
630     // want to accept this with some probability that dies off the
631     // farther we are from nTarget
632 chrisfen 417
633 gezelter 403 // roughly: x = new_atoms - nTarget
634     // Pacc(x) = exp(- a * x)
635     // where a = penalty / (average atoms per molecule)
636 chrisfen 417
637 tim 963 x = (RealType)(new_atoms - nTarget);
638 gezelter 403 y = myRandom->rand();
639 chrisfen 417
640 gezelter 403 if (y < exp(- a * x)) {
641     molToProcMap[i] = which_proc;
642     atomsPerProc[which_proc] += add_atoms;
643 chrisfen 417
644 gezelter 403 done = 1;
645     continue;
646     } else {
647     continue;
648     }
649     }
650     }
651 chrisfen 417
652 gezelter 403 delete myRandom;
653 gezelter 1803
654 gezelter 403 // Spray out this nonsense to all other processors:
655 gezelter 1798 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
656 gezelter 403 } else {
657 chrisfen 417
658 gezelter 403 // Listen to your marching orders from processor 0:
659 gezelter 1798 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
660 gezelter 1803
661 gezelter 403 }
662 chrisfen 417
663 gezelter 403 info->setMolToProcMap(molToProcMap);
664     sprintf(checkPointMsg,
665     "Successfully divided the molecules among the processors.\n");
666 gezelter 1241 errorCheckPoint();
667 gezelter 403 }
668 chrisfen 417
669 gezelter 403 #endif
670 chrisfen 417
671 gezelter 403 void SimCreator::createMolecules(SimInfo *info) {
672     MoleculeCreator molCreator;
673     int stampId;
674 chrisfen 417
675 gezelter 403 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
676 chrisfen 417
677 gezelter 403 #ifdef IS_MPI
678 chrisfen 417
679 gezelter 403 if (info->getMolToProc(i) == worldRank) {
680     #endif
681 chrisfen 417
682 gezelter 403 stampId = info->getMoleculeStampId(i);
683 gezelter 1600 Molecule * mol = molCreator.createMolecule(info->getForceField(),
684     info->getMoleculeStamp(stampId),
685     stampId, i,
686     info->getLocalIndexManager());
687 chrisfen 417
688 gezelter 403 info->addMolecule(mol);
689 chrisfen 417
690 gezelter 403 #ifdef IS_MPI
691 chrisfen 417
692 gezelter 403 }
693 chrisfen 417
694 gezelter 403 #endif
695 chrisfen 417
696 gezelter 403 } //end for(int i=0)
697     }
698 chrisfen 417
699 gezelter 1710 int SimCreator::computeStorageLayout(SimInfo* info) {
700    
701 gezelter 1711 Globals* simParams = info->getSimParams();
702 gezelter 1710 int nRigidBodies = info->getNGlobalRigidBodies();
703     set<AtomType*> atomTypes = info->getSimulatedAtomTypes();
704     set<AtomType*>::iterator i;
705     bool hasDirectionalAtoms = false;
706     bool hasFixedCharge = false;
707 gezelter 1794 bool hasDipoles = false;
708     bool hasQuadrupoles = false;
709 gezelter 1710 bool hasPolarizable = false;
710     bool hasFluctuatingCharge = false;
711     bool hasMetallic = false;
712     int storageLayout = 0;
713     storageLayout |= DataStorage::dslPosition;
714     storageLayout |= DataStorage::dslVelocity;
715     storageLayout |= DataStorage::dslForce;
716    
717     for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
718    
719     DirectionalAdapter da = DirectionalAdapter( (*i) );
720     MultipoleAdapter ma = MultipoleAdapter( (*i) );
721     EAMAdapter ea = EAMAdapter( (*i) );
722     SuttonChenAdapter sca = SuttonChenAdapter( (*i) );
723     PolarizableAdapter pa = PolarizableAdapter( (*i) );
724     FixedChargeAdapter fca = FixedChargeAdapter( (*i) );
725     FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) );
726    
727     if (da.isDirectional()){
728     hasDirectionalAtoms = true;
729     }
730 gezelter 1794 if (ma.isDipole()){
731     hasDipoles = true;
732 gezelter 1710 }
733 gezelter 1794 if (ma.isQuadrupole()){
734     hasQuadrupoles = true;
735     }
736 gezelter 1710 if (ea.isEAM() || sca.isSuttonChen()){
737     hasMetallic = true;
738     }
739     if ( fca.isFixedCharge() ){
740     hasFixedCharge = true;
741     }
742     if ( fqa.isFluctuatingCharge() ){
743     hasFluctuatingCharge = true;
744     }
745     if ( pa.isPolarizable() ){
746     hasPolarizable = true;
747     }
748     }
749    
750     if (nRigidBodies > 0 || hasDirectionalAtoms) {
751     storageLayout |= DataStorage::dslAmat;
752     if(storageLayout & DataStorage::dslVelocity) {
753     storageLayout |= DataStorage::dslAngularMomentum;
754     }
755     if (storageLayout & DataStorage::dslForce) {
756     storageLayout |= DataStorage::dslTorque;
757     }
758     }
759 gezelter 1794 if (hasDipoles) {
760     storageLayout |= DataStorage::dslDipole;
761     }
762     if (hasQuadrupoles) {
763     storageLayout |= DataStorage::dslQuadrupole;
764     }
765 gezelter 1710 if (hasFixedCharge || hasFluctuatingCharge) {
766     storageLayout |= DataStorage::dslSkippedCharge;
767     }
768     if (hasMetallic) {
769     storageLayout |= DataStorage::dslDensity;
770     storageLayout |= DataStorage::dslFunctional;
771     storageLayout |= DataStorage::dslFunctionalDerivative;
772     }
773     if (hasPolarizable) {
774     storageLayout |= DataStorage::dslElectricField;
775     }
776     if (hasFluctuatingCharge){
777     storageLayout |= DataStorage::dslFlucQPosition;
778     if(storageLayout & DataStorage::dslVelocity) {
779     storageLayout |= DataStorage::dslFlucQVelocity;
780     }
781     if (storageLayout & DataStorage::dslForce) {
782     storageLayout |= DataStorage::dslFlucQForce;
783     }
784     }
785 gezelter 1714
786     // if the user has asked for them, make sure we've got the memory for the
787     // objects defined.
788 gezelter 1711
789     if (simParams->getOutputParticlePotential()) {
790     storageLayout |= DataStorage::dslParticlePot;
791     }
792 gezelter 1723
793     if (simParams->havePrintHeatFlux()) {
794     if (simParams->getPrintHeatFlux()) {
795     storageLayout |= DataStorage::dslParticlePot;
796     }
797     }
798    
799 gezelter 1714 if (simParams->getOutputElectricField()) {
800     storageLayout |= DataStorage::dslElectricField;
801     }
802 gezelter 1787
803 gezelter 1714 if (simParams->getOutputFluctuatingCharges()) {
804     storageLayout |= DataStorage::dslFlucQPosition;
805     storageLayout |= DataStorage::dslFlucQVelocity;
806     storageLayout |= DataStorage::dslFlucQForce;
807     }
808 gezelter 1711
809 gezelter 1710 return storageLayout;
810     }
811    
812 gezelter 403 void SimCreator::setGlobalIndex(SimInfo *info) {
813     SimInfo::MoleculeIterator mi;
814     Molecule::AtomIterator ai;
815     Molecule::RigidBodyIterator ri;
816     Molecule::CutoffGroupIterator ci;
817 tim 1024 Molecule::IntegrableObjectIterator ioi;
818 gezelter 403 Molecule * mol;
819     Atom * atom;
820     RigidBody * rb;
821     CutoffGroup * cg;
822     int beginAtomIndex;
823     int beginRigidBodyIndex;
824     int beginCutoffGroupIndex;
825     int nGlobalAtoms = info->getNGlobalAtoms();
826 gezelter 1803 int nGlobalRigidBodies = info->getNGlobalRigidBodies();
827 chrisfen 417
828 gezelter 403 beginAtomIndex = 0;
829 gezelter 1757 //rigidbody's index begins right after atom's
830     beginRigidBodyIndex = info->getNGlobalAtoms();
831 gezelter 403 beginCutoffGroupIndex = 0;
832 gezelter 1600
833     for(int i = 0; i < info->getNGlobalMolecules(); i++) {
834 chrisfen 417
835 gezelter 1600 #ifdef IS_MPI
836     if (info->getMolToProc(i) == worldRank) {
837     #endif
838     // stuff to do if I own this molecule
839     mol = info->getMoleculeByGlobalIndex(i);
840    
841     //local index(index in DataStorge) of atom is important
842     for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
843     atom->setGlobalIndex(beginAtomIndex++);
844     }
845    
846     for(rb = mol->beginRigidBody(ri); rb != NULL;
847     rb = mol->nextRigidBody(ri)) {
848     rb->setGlobalIndex(beginRigidBodyIndex++);
849     }
850    
851     //local index of cutoff group is trivial, it only depends on
852     //the order of travesing
853     for(cg = mol->beginCutoffGroup(ci); cg != NULL;
854     cg = mol->nextCutoffGroup(ci)) {
855     cg->setGlobalIndex(beginCutoffGroupIndex++);
856     }
857    
858     #ifdef IS_MPI
859     } else {
860    
861     // stuff to do if I don't own this molecule
862    
863     int stampId = info->getMoleculeStampId(i);
864     MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
865    
866     beginAtomIndex += stamp->getNAtoms();
867     beginRigidBodyIndex += stamp->getNRigidBodies();
868     beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
869 gezelter 403 }
870 gezelter 1600 #endif
871    
872     } //end for(int i=0)
873    
874 gezelter 403 //fill globalGroupMembership
875 gezelter 1600 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
876 gezelter 403 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
877     for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
878 chrisfen 417
879 gezelter 403 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
880     globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
881     }
882 chrisfen 417
883 gezelter 403 }
884     }
885 gezelter 1577
886 gezelter 403 #ifdef IS_MPI
887     // Since the globalGroupMembership has been zero filled and we've only
888     // poked values into the atoms we know, we can do an Allreduce
889     // to get the full globalGroupMembership array (We think).
890     // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
891     // docs said we could.
892 gezelter 1313 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
893 gezelter 1798 MPI::COMM_WORLD.Allreduce(&globalGroupMembership[0],
894     &tmpGroupMembership[0], nGlobalAtoms,
895     MPI::INT, MPI::SUM);
896 gezelter 403 info->setGlobalGroupMembership(tmpGroupMembership);
897     #else
898     info->setGlobalGroupMembership(globalGroupMembership);
899     #endif
900 chrisfen 417
901 gezelter 403 //fill molMembership
902 gezelter 1803 std::vector<int> globalMolMembership(info->getNGlobalAtoms() +
903     info->getNGlobalRigidBodies(), 0);
904 gezelter 403
905 gezelter 1803 for(mol = info->beginMolecule(mi); mol != NULL;
906     mol = info->nextMolecule(mi)) {
907 gezelter 403 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
908     globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
909     }
910 gezelter 1803 for (rb = mol->beginRigidBody(ri); rb != NULL;
911     rb = mol->nextRigidBody(ri)) {
912     globalMolMembership[rb->getGlobalIndex()] = mol->getGlobalIndex();
913     }
914 gezelter 403 }
915 chrisfen 417
916 gezelter 403 #ifdef IS_MPI
917 gezelter 1803 std::vector<int> tmpMolMembership(info->getNGlobalAtoms() +
918     info->getNGlobalRigidBodies(), 0);
919 gezelter 1798 MPI::COMM_WORLD.Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
920 gezelter 1803 nGlobalAtoms + nGlobalRigidBodies,
921 gezelter 1798 MPI::INT, MPI::SUM);
922 chrisfen 417
923 gezelter 403 info->setGlobalMolMembership(tmpMolMembership);
924     #else
925     info->setGlobalMolMembership(globalMolMembership);
926     #endif
927 tim 1024
928     // nIOPerMol holds the number of integrable objects per molecule
929     // here the molecules are listed by their global indices.
930    
931     std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
932 gezelter 1803 for (mol = info->beginMolecule(mi); mol != NULL;
933     mol = info->nextMolecule(mi)) {
934 tim 1024 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
935     }
936 chrisfen 417
937 tim 1024 #ifdef IS_MPI
938     std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
939 gezelter 1798 MPI::COMM_WORLD.Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
940     info->getNGlobalMolecules(), MPI::INT, MPI::SUM);
941 tim 1024 #else
942     std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
943     #endif
944    
945 gezelter 1313 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
946    
947     int startingIndex = 0;
948     for (int i = 0; i < info->getNGlobalMolecules(); i++) {
949     startingIOIndexForMol[i] = startingIndex;
950     startingIndex += numIntegrableObjectsPerMol[i];
951     }
952    
953     std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
954 gezelter 1803 for (mol = info->beginMolecule(mi); mol != NULL;
955     mol = info->nextMolecule(mi)) {
956 tim 1024 int myGlobalIndex = mol->getGlobalIndex();
957     int globalIO = startingIOIndexForMol[myGlobalIndex];
958 gezelter 1769 for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
959     sd = mol->nextIntegrableObject(ioi)) {
960     sd->setGlobalIntegrableObjectIndex(globalIO);
961     IOIndexToIntegrableObject[globalIO] = sd;
962 gezelter 1313 globalIO++;
963 tim 1024 }
964     }
965 gezelter 1597
966 gezelter 1313 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
967    
968 gezelter 403 }
969 chrisfen 417
970 tim 1024 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
971 gezelter 403
972 tim 1024 DumpReader reader(info, mdFileName);
973 gezelter 403 int nframes = reader.getNFrames();
974 gezelter 1825
975 gezelter 403 if (nframes > 0) {
976     reader.readFrame(nframes - 1);
977     } else {
978     //invalid initial coordinate file
979 chrisfen 417 sprintf(painCave.errMsg,
980     "Initial configuration file %s should at least contain one frame\n",
981 tim 1024 mdFileName.c_str());
982 gezelter 403 painCave.isFatal = 1;
983     simError();
984     }
985     //copy the current snapshot to previous snapshot
986     info->getSnapshotManager()->advance();
987     }
988 chrisfen 417
989 gezelter 1390 } //end namespace OpenMD
990 gezelter 403
991    

Properties

Name Value
svn:keywords Author Id Revision Date