ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/OpenMD/branches/development/src/brains/SimCreator.cpp
Revision: 1798
Committed: Thu Sep 13 14:10:11 2012 UTC (12 years, 7 months ago) by gezelter
File size: 31393 byte(s)
Log Message:
Merged trunk changes into the 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 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 sprintf(painCave.errMsg,
573 "I've tried 100 times to assign molecule %d to a "
574 " processor, but can't find a good spot.\n"
575 "I'm assigning it at random to processor %d.\n",
576 i, which_proc);
577
578 painCave.isFatal = 0;
579 simError();
580
581 molToProcMap[i] = which_proc;
582 atomsPerProc[which_proc] += add_atoms;
583
584 done = 1;
585 continue;
586 }
587
588 // If we can add this molecule to this processor without sending
589 // it above nTarget, then go ahead and do it:
590
591 if (new_atoms <= nTarget) {
592 molToProcMap[i] = which_proc;
593 atomsPerProc[which_proc] += add_atoms;
594
595 done = 1;
596 continue;
597 }
598
599 // The only situation left is when new_atoms > nTarget. We
600 // want to accept this with some probability that dies off the
601 // farther we are from nTarget
602
603 // roughly: x = new_atoms - nTarget
604 // Pacc(x) = exp(- a * x)
605 // where a = penalty / (average atoms per molecule)
606
607 x = (RealType)(new_atoms - nTarget);
608 y = myRandom->rand();
609
610 if (y < exp(- a * x)) {
611 molToProcMap[i] = which_proc;
612 atomsPerProc[which_proc] += add_atoms;
613
614 done = 1;
615 continue;
616 } else {
617 continue;
618 }
619 }
620 }
621
622 delete myRandom;
623
624 // Spray out this nonsense to all other processors:
625 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
626 } else {
627
628 // Listen to your marching orders from processor 0:
629 MPI::COMM_WORLD.Bcast(&molToProcMap[0], nGlobalMols, MPI::INT, 0);
630 }
631
632 info->setMolToProcMap(molToProcMap);
633 sprintf(checkPointMsg,
634 "Successfully divided the molecules among the processors.\n");
635 errorCheckPoint();
636 }
637
638 #endif
639
640 void SimCreator::createMolecules(SimInfo *info) {
641 MoleculeCreator molCreator;
642 int stampId;
643
644 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
645
646 #ifdef IS_MPI
647
648 if (info->getMolToProc(i) == worldRank) {
649 #endif
650
651 stampId = info->getMoleculeStampId(i);
652 Molecule * mol = molCreator.createMolecule(info->getForceField(),
653 info->getMoleculeStamp(stampId),
654 stampId, i,
655 info->getLocalIndexManager());
656
657 info->addMolecule(mol);
658
659 #ifdef IS_MPI
660
661 }
662
663 #endif
664
665 } //end for(int i=0)
666 }
667
668 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 hasDipoles = false;
677 bool hasQuadrupoles = false;
678 bool hasPolarizable = false;
679 bool hasFluctuatingCharge = false;
680 bool hasMetallic = false;
681 int storageLayout = 0;
682 storageLayout |= DataStorage::dslPosition;
683 storageLayout |= DataStorage::dslVelocity;
684 storageLayout |= DataStorage::dslForce;
685
686 for (i = atomTypes.begin(); i != atomTypes.end(); ++i) {
687
688 DirectionalAdapter da = DirectionalAdapter( (*i) );
689 MultipoleAdapter ma = MultipoleAdapter( (*i) );
690 EAMAdapter ea = EAMAdapter( (*i) );
691 SuttonChenAdapter sca = SuttonChenAdapter( (*i) );
692 PolarizableAdapter pa = PolarizableAdapter( (*i) );
693 FixedChargeAdapter fca = FixedChargeAdapter( (*i) );
694 FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) );
695
696 if (da.isDirectional()){
697 hasDirectionalAtoms = true;
698 }
699 if (ma.isDipole()){
700 hasDipoles = true;
701 }
702 if (ma.isQuadrupole()){
703 hasQuadrupoles = true;
704 }
705 if (ea.isEAM() || sca.isSuttonChen()){
706 hasMetallic = true;
707 }
708 if ( fca.isFixedCharge() ){
709 hasFixedCharge = true;
710 }
711 if ( fqa.isFluctuatingCharge() ){
712 hasFluctuatingCharge = true;
713 }
714 if ( pa.isPolarizable() ){
715 hasPolarizable = true;
716 }
717 }
718
719 if (nRigidBodies > 0 || hasDirectionalAtoms) {
720 storageLayout |= DataStorage::dslAmat;
721 if(storageLayout & DataStorage::dslVelocity) {
722 storageLayout |= DataStorage::dslAngularMomentum;
723 }
724 if (storageLayout & DataStorage::dslForce) {
725 storageLayout |= DataStorage::dslTorque;
726 }
727 }
728 if (hasDipoles) {
729 storageLayout |= DataStorage::dslDipole;
730 }
731 if (hasQuadrupoles) {
732 storageLayout |= DataStorage::dslQuadrupole;
733 }
734 if (hasFixedCharge || hasFluctuatingCharge) {
735 storageLayout |= DataStorage::dslSkippedCharge;
736 }
737 if (hasMetallic) {
738 storageLayout |= DataStorage::dslDensity;
739 storageLayout |= DataStorage::dslFunctional;
740 storageLayout |= DataStorage::dslFunctionalDerivative;
741 }
742 if (hasPolarizable) {
743 storageLayout |= DataStorage::dslElectricField;
744 }
745 if (hasFluctuatingCharge){
746 storageLayout |= DataStorage::dslFlucQPosition;
747 if(storageLayout & DataStorage::dslVelocity) {
748 storageLayout |= DataStorage::dslFlucQVelocity;
749 }
750 if (storageLayout & DataStorage::dslForce) {
751 storageLayout |= DataStorage::dslFlucQForce;
752 }
753 }
754
755 // if the user has asked for them, make sure we've got the memory for the
756 // objects defined.
757
758 if (simParams->getOutputParticlePotential()) {
759 storageLayout |= DataStorage::dslParticlePot;
760 }
761
762 if (simParams->havePrintHeatFlux()) {
763 if (simParams->getPrintHeatFlux()) {
764 storageLayout |= DataStorage::dslParticlePot;
765 }
766 }
767
768 if (simParams->getOutputElectricField()) {
769 storageLayout |= DataStorage::dslElectricField;
770 }
771
772 if (simParams->getOutputFluctuatingCharges()) {
773 storageLayout |= DataStorage::dslFlucQPosition;
774 storageLayout |= DataStorage::dslFlucQVelocity;
775 storageLayout |= DataStorage::dslFlucQForce;
776 }
777
778 return storageLayout;
779 }
780
781 void SimCreator::setGlobalIndex(SimInfo *info) {
782 SimInfo::MoleculeIterator mi;
783 Molecule::AtomIterator ai;
784 Molecule::RigidBodyIterator ri;
785 Molecule::CutoffGroupIterator ci;
786 Molecule::IntegrableObjectIterator ioi;
787 Molecule * mol;
788 Atom * atom;
789 RigidBody * rb;
790 CutoffGroup * cg;
791 int beginAtomIndex;
792 int beginRigidBodyIndex;
793 int beginCutoffGroupIndex;
794 int nGlobalAtoms = info->getNGlobalAtoms();
795
796 beginAtomIndex = 0;
797 //rigidbody's index begins right after atom's
798 beginRigidBodyIndex = info->getNGlobalAtoms();
799 beginCutoffGroupIndex = 0;
800
801 for(int i = 0; i < info->getNGlobalMolecules(); i++) {
802
803 #ifdef IS_MPI
804 if (info->getMolToProc(i) == worldRank) {
805 #endif
806 // stuff to do if I own this molecule
807 mol = info->getMoleculeByGlobalIndex(i);
808
809 //local index(index in DataStorge) of atom is important
810 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
811 atom->setGlobalIndex(beginAtomIndex++);
812 }
813
814 for(rb = mol->beginRigidBody(ri); rb != NULL;
815 rb = mol->nextRigidBody(ri)) {
816 rb->setGlobalIndex(beginRigidBodyIndex++);
817 }
818
819 //local index of cutoff group is trivial, it only depends on
820 //the order of travesing
821 for(cg = mol->beginCutoffGroup(ci); cg != NULL;
822 cg = mol->nextCutoffGroup(ci)) {
823 cg->setGlobalIndex(beginCutoffGroupIndex++);
824 }
825
826 #ifdef IS_MPI
827 } else {
828
829 // stuff to do if I don't own this molecule
830
831 int stampId = info->getMoleculeStampId(i);
832 MoleculeStamp* stamp = info->getMoleculeStamp(stampId);
833
834 beginAtomIndex += stamp->getNAtoms();
835 beginRigidBodyIndex += stamp->getNRigidBodies();
836 beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms();
837 }
838 #endif
839
840 } //end for(int i=0)
841
842 //fill globalGroupMembership
843 std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0);
844 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
845 for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) {
846
847 for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) {
848 globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex();
849 }
850
851 }
852 }
853
854 #ifdef IS_MPI
855 // Since the globalGroupMembership has been zero filled and we've only
856 // poked values into the atoms we know, we can do an Allreduce
857 // to get the full globalGroupMembership array (We think).
858 // This would be prettier if we could use MPI_IN_PLACE like the MPI-2
859 // docs said we could.
860 std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0);
861 MPI::COMM_WORLD.Allreduce(&globalGroupMembership[0],
862 &tmpGroupMembership[0], nGlobalAtoms,
863 MPI::INT, MPI::SUM);
864 info->setGlobalGroupMembership(tmpGroupMembership);
865 #else
866 info->setGlobalGroupMembership(globalGroupMembership);
867 #endif
868
869 //fill molMembership
870 std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0);
871
872 for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
873 for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) {
874 globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex();
875 }
876 }
877
878 #ifdef IS_MPI
879 std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0);
880 MPI::COMM_WORLD.Allreduce(&globalMolMembership[0], &tmpMolMembership[0],
881 nGlobalAtoms,
882 MPI::INT, MPI::SUM);
883
884 info->setGlobalMolMembership(tmpMolMembership);
885 #else
886 info->setGlobalMolMembership(globalMolMembership);
887 #endif
888
889 // nIOPerMol holds the number of integrable objects per molecule
890 // here the molecules are listed by their global indices.
891
892 std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0);
893 for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
894 nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects();
895 }
896
897 #ifdef IS_MPI
898 std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0);
899 MPI::COMM_WORLD.Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0],
900 info->getNGlobalMolecules(), MPI::INT, MPI::SUM);
901 #else
902 std::vector<int> numIntegrableObjectsPerMol = nIOPerMol;
903 #endif
904
905 std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules());
906
907 int startingIndex = 0;
908 for (int i = 0; i < info->getNGlobalMolecules(); i++) {
909 startingIOIndexForMol[i] = startingIndex;
910 startingIndex += numIntegrableObjectsPerMol[i];
911 }
912
913 std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL);
914 for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) {
915 int myGlobalIndex = mol->getGlobalIndex();
916 int globalIO = startingIOIndexForMol[myGlobalIndex];
917 for (StuntDouble* sd = mol->beginIntegrableObject(ioi); sd != NULL;
918 sd = mol->nextIntegrableObject(ioi)) {
919 sd->setGlobalIntegrableObjectIndex(globalIO);
920 IOIndexToIntegrableObject[globalIO] = sd;
921 globalIO++;
922 }
923 }
924
925 info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject);
926
927 }
928
929 void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) {
930 Globals* simParams;
931
932 simParams = info->getSimParams();
933
934 DumpReader reader(info, mdFileName);
935 int nframes = reader.getNFrames();
936
937 if (nframes > 0) {
938 reader.readFrame(nframes - 1);
939 } else {
940 //invalid initial coordinate file
941 sprintf(painCave.errMsg,
942 "Initial configuration file %s should at least contain one frame\n",
943 mdFileName.c_str());
944 painCave.isFatal = 1;
945 simError();
946 }
947 //copy the current snapshot to previous snapshot
948 info->getSnapshotManager()->advance();
949 }
950
951 } //end namespace OpenMD
952
953

Properties

Name Value
svn:keywords Author Id Revision Date