ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/trunk/src/brains/SimCreator.cpp
Revision: 1782
Committed: Wed Aug 22 02:28:28 2012 UTC (12 years, 8 months ago) by gezelter
File size: 31104 byte(s)
Log Message:
MERGE OpenMD development branch 1465:1781 into trunk

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

Properties

Name Value
svn:keywords Author Id Revision Date