ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1803
Committed: Wed Oct 3 14:20:07 2012 UTC (12 years, 6 months ago) by gezelter
File size: 31910 byte(s)
Log Message:
Merging trunk changes back to development branch

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date