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

File Contents

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

Properties

Name Value
svn:keywords Author Id Revision Date