| 1 | gezelter | 403 | /* | 
| 2 |  |  | * Copyright (c) 2005 The University of Notre Dame. All Rights Reserved. | 
| 3 |  |  | * | 
| 4 |  |  | * The University of Notre Dame grants you ("Licensee") a | 
| 5 |  |  | * non-exclusive, royalty free, license to use, modify and | 
| 6 |  |  | * redistribute this software in source and binary code form, provided | 
| 7 |  |  | * that the following conditions are met: | 
| 8 |  |  | * | 
| 9 | gezelter | 1390 | * 1. Redistributions of source code must retain the above copyright | 
| 10 | gezelter | 403 | *    notice, this list of conditions and the following disclaimer. | 
| 11 |  |  | * | 
| 12 | gezelter | 1390 | * 2. Redistributions in binary form must reproduce the above copyright | 
| 13 | gezelter | 403 | *    notice, this list of conditions and the following disclaimer in the | 
| 14 |  |  | *    documentation and/or other materials provided with the | 
| 15 |  |  | *    distribution. | 
| 16 |  |  | * | 
| 17 |  |  | * This software is provided "AS IS," without a warranty of any | 
| 18 |  |  | * kind. All express or implied conditions, representations and | 
| 19 |  |  | * warranties, including any implied warranty of merchantability, | 
| 20 |  |  | * fitness for a particular purpose or non-infringement, are hereby | 
| 21 |  |  | * excluded.  The University of Notre Dame and its licensors shall not | 
| 22 |  |  | * be liable for any damages suffered by licensee as a result of | 
| 23 |  |  | * using, modifying or distributing the software or its | 
| 24 |  |  | * derivatives. In no event will the University of Notre Dame or its | 
| 25 |  |  | * licensors be liable for any lost revenue, profit or data, or for | 
| 26 |  |  | * direct, indirect, special, consequential, incidental or punitive | 
| 27 |  |  | * damages, however caused and regardless of the theory of liability, | 
| 28 |  |  | * arising out of the use of or inability to use software, even if the | 
| 29 |  |  | * University of Notre Dame has been advised of the possibility of | 
| 30 |  |  | * such damages. | 
| 31 | gezelter | 1390 | * | 
| 32 |  |  | * SUPPORT OPEN SCIENCE!  If you use OpenMD or its source code in your | 
| 33 |  |  | * research, please cite the appropriate papers when you publish your | 
| 34 |  |  | * work.  Good starting points are: | 
| 35 |  |  | * | 
| 36 |  |  | * [1]  Meineke, et al., J. Comp. Chem. 26, 252-271 (2005). | 
| 37 |  |  | * [2]  Fennell & Gezelter, J. Chem. Phys. 124, 234104 (2006). | 
| 38 |  |  | * [3]  Sun, Lin & Gezelter, J. Chem. Phys. 128, 24107 (2008). | 
| 39 | gezelter | 1665 | * [4]  Kuang & Gezelter,  J. Chem. Phys. 133, 164101 (2010). | 
| 40 |  |  | * [5]  Vardeman, Stocker & Gezelter, J. Chem. Theory Comput. 7, 834 (2011). | 
| 41 | gezelter | 403 | */ | 
| 42 |  |  |  | 
| 43 |  |  | /** | 
| 44 |  |  | * @file SimCreator.cpp | 
| 45 |  |  | * @author tlin | 
| 46 |  |  | * @date 11/03/2004 | 
| 47 |  |  | * @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 |  |  | #include "UseTheForce/ForceFieldFactory.hpp" | 
| 60 |  |  | #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 | 1710 | #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 | 1627 | #include "mpi.h" | 
| 89 | gezelter | 403 | #include "math/ParallelRandNumGen.hpp" | 
| 90 |  |  | #endif | 
| 91 |  |  |  | 
| 92 | gezelter | 1390 | namespace OpenMD { | 
| 93 | chrisfen | 417 |  | 
| 94 | gezelter | 1613 | 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 | 1613 | 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 | 1613 |  | 
| 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 | 1613 | 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 | 1469 |  | 
| 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 |  |  | int mdOffset; | 
| 260 | gezelter | 1613 | 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 | 1613 |  | 
| 295 |  |  | // found the correct opening string, now try to get the file | 
| 296 |  |  | // format version number. | 
| 297 | tim | 1024 |  | 
| 298 | gezelter | 1613 | 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 | 1613 | Globals* simParams = parseFile(rawMetaDataStream, mdFileName, mdFileVersion, | 
| 371 |  |  | metaDataBlockStart + 1); | 
| 372 | chrisfen | 417 |  | 
| 373 | gezelter | 403 | //create the force field | 
| 374 | gezelter | 1277 | ForceField * ff = ForceFieldFactory::getInstance()->createForceField(simParams->getForceField()); | 
| 375 |  |  |  | 
| 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 | 1710 | //find the storage layout | 
| 426 |  |  |  | 
| 427 |  |  | int storageLayout = computeStorageLayout(info); | 
| 428 |  |  |  | 
| 429 |  |  | cerr << "computed Storage Layout = " << storageLayout << "\n"; | 
| 430 |  |  |  | 
| 431 | gezelter | 945 | //allocate memory for DataStorage(circular reference, need to | 
| 432 |  |  | //break it) | 
| 433 | gezelter | 1710 | info->setSnapshotManager(new SimSnapshotManager(info, storageLayout)); | 
| 434 | gezelter | 403 |  | 
| 435 | gezelter | 945 | //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 | gezelter | 403 | setGlobalIndex(info); | 
| 441 | chrisfen | 417 |  | 
| 442 | gezelter | 1287 | //Although addInteractionPairs is called inside SimInfo's addMolecule | 
| 443 | gezelter | 945 | //method, at that point atoms don't have the global index yet | 
| 444 |  |  | //(their global index are all initialized to -1).  Therefore we | 
| 445 | gezelter | 1287 | //have to call addInteractionPairs explicitly here. A way to work | 
| 446 | gezelter | 945 | //around is that we can determine the beginning global indices of | 
| 447 |  |  | //atoms before they get created. | 
| 448 | gezelter | 403 | SimInfo::MoleculeIterator mi; | 
| 449 |  |  | Molecule* mol; | 
| 450 |  |  | for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { | 
| 451 | gezelter | 1287 | info->addInteractionPairs(mol); | 
| 452 | gezelter | 403 | } | 
| 453 |  |  |  | 
| 454 |  |  | if (loadInitCoords) | 
| 455 | tim | 1024 | loadCoordinates(info, mdFileName); | 
| 456 | gezelter | 403 | return info; | 
| 457 |  |  | } | 
| 458 | chrisfen | 417 |  | 
| 459 | gezelter | 403 | void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) { | 
| 460 | chrisfen | 417 |  | 
| 461 | tim | 749 | //figure out the output file names | 
| 462 | gezelter | 403 | std::string prefix; | 
| 463 | chrisfen | 417 |  | 
| 464 | gezelter | 403 | #ifdef IS_MPI | 
| 465 | chrisfen | 417 |  | 
| 466 | gezelter | 403 | 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 | chrisfen | 417 |  | 
| 475 | gezelter | 403 | info->setFinalConfigFileName(prefix + ".eor"); | 
| 476 |  |  | info->setDumpFileName(prefix + ".dump"); | 
| 477 |  |  | info->setStatFileName(prefix + ".stat"); | 
| 478 | chrisfen | 417 | info->setRestFileName(prefix + ".zang"); | 
| 479 |  |  |  | 
| 480 | gezelter | 403 | #ifdef IS_MPI | 
| 481 | chrisfen | 417 |  | 
| 482 | gezelter | 403 | } | 
| 483 | chrisfen | 417 |  | 
| 484 | gezelter | 403 | #endif | 
| 485 | chrisfen | 417 |  | 
| 486 | gezelter | 403 | } | 
| 487 | chrisfen | 417 |  | 
| 488 | gezelter | 403 | #ifdef IS_MPI | 
| 489 |  |  | void SimCreator::divideMolecules(SimInfo *info) { | 
| 490 | tim | 963 | RealType numerator; | 
| 491 |  |  | RealType denominator; | 
| 492 |  |  | RealType precast; | 
| 493 |  |  | RealType x; | 
| 494 |  |  | RealType y; | 
| 495 |  |  | RealType a; | 
| 496 | gezelter | 403 | 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 |  |  | MPI_Comm_size(MPI_COMM_WORLD, &nProcessors); | 
| 511 | chrisfen | 417 |  | 
| 512 | gezelter | 403 | 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 | gezelter | 1390 | "\tsingle-processor version of OpenMD.\n", nProcessors, nGlobalMols); | 
| 520 | chrisfen | 417 |  | 
| 521 | gezelter | 403 | painCave.isFatal = 1; | 
| 522 |  |  | simError(); | 
| 523 |  |  | } | 
| 524 | chrisfen | 417 |  | 
| 525 | gezelter | 403 | 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 | chrisfen | 417 |  | 
| 535 |  |  |  | 
| 536 | gezelter | 403 | a = 3.0 * nGlobalMols / info->getNGlobalAtoms(); | 
| 537 | chrisfen | 417 |  | 
| 538 | gezelter | 403 | //initialize atomsPerProc | 
| 539 |  |  | atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0); | 
| 540 | chrisfen | 417 |  | 
| 541 | gezelter | 403 | if (worldRank == 0) { | 
| 542 |  |  | numerator = info->getNGlobalAtoms(); | 
| 543 |  |  | denominator = nProcessors; | 
| 544 |  |  | precast = numerator / denominator; | 
| 545 |  |  | nTarget = (int)(precast + 0.5); | 
| 546 | chrisfen | 417 |  | 
| 547 | gezelter | 403 | for(i = 0; i < nGlobalMols; i++) { | 
| 548 |  |  | done = 0; | 
| 549 |  |  | loops = 0; | 
| 550 | chrisfen | 417 |  | 
| 551 | gezelter | 403 | while (!done) { | 
| 552 |  |  | loops++; | 
| 553 | chrisfen | 417 |  | 
| 554 | gezelter | 403 | // Pick a processor at random | 
| 555 | chrisfen | 417 |  | 
| 556 | gezelter | 403 | which_proc = (int) (myRandom->rand() * nProcessors); | 
| 557 | chrisfen | 417 |  | 
| 558 | gezelter | 403 | //get the molecule stamp first | 
| 559 |  |  | int stampId = info->getMoleculeStampId(i); | 
| 560 |  |  | MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId); | 
| 561 | chrisfen | 417 |  | 
| 562 | gezelter | 403 | // 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 | chrisfen | 417 |  | 
| 567 | gezelter | 403 | // 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 | chrisfen | 417 |  | 
| 571 | gezelter | 403 | 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 | chrisfen | 417 |  | 
| 578 | gezelter | 403 | painCave.isFatal = 0; | 
| 579 |  |  | simError(); | 
| 580 | chrisfen | 417 |  | 
| 581 | gezelter | 403 | molToProcMap[i] = which_proc; | 
| 582 |  |  | atomsPerProc[which_proc] += add_atoms; | 
| 583 | chrisfen | 417 |  | 
| 584 | gezelter | 403 | done = 1; | 
| 585 |  |  | continue; | 
| 586 |  |  | } | 
| 587 | chrisfen | 417 |  | 
| 588 | gezelter | 403 | // If we can add this molecule to this processor without sending | 
| 589 |  |  | // it above nTarget, then go ahead and do it: | 
| 590 | chrisfen | 417 |  | 
| 591 | gezelter | 403 | if (new_atoms <= nTarget) { | 
| 592 |  |  | molToProcMap[i] = which_proc; | 
| 593 |  |  | atomsPerProc[which_proc] += add_atoms; | 
| 594 | chrisfen | 417 |  | 
| 595 | gezelter | 403 | done = 1; | 
| 596 |  |  | continue; | 
| 597 |  |  | } | 
| 598 | chrisfen | 417 |  | 
| 599 | gezelter | 403 | // 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 | chrisfen | 417 |  | 
| 603 | gezelter | 403 | // roughly:  x = new_atoms - nTarget | 
| 604 |  |  | //           Pacc(x) = exp(- a * x) | 
| 605 |  |  | // where a = penalty / (average atoms per molecule) | 
| 606 | chrisfen | 417 |  | 
| 607 | tim | 963 | x = (RealType)(new_atoms - nTarget); | 
| 608 | gezelter | 403 | y = myRandom->rand(); | 
| 609 | chrisfen | 417 |  | 
| 610 | gezelter | 403 | if (y < exp(- a * x)) { | 
| 611 |  |  | molToProcMap[i] = which_proc; | 
| 612 |  |  | atomsPerProc[which_proc] += add_atoms; | 
| 613 | chrisfen | 417 |  | 
| 614 | gezelter | 403 | done = 1; | 
| 615 |  |  | continue; | 
| 616 |  |  | } else { | 
| 617 |  |  | continue; | 
| 618 |  |  | } | 
| 619 |  |  | } | 
| 620 |  |  | } | 
| 621 | chrisfen | 417 |  | 
| 622 | gezelter | 403 | delete myRandom; | 
| 623 | chrisfen | 417 |  | 
| 624 | gezelter | 403 | // Spray out this nonsense to all other processors: | 
| 625 | chrisfen | 417 |  | 
| 626 | gezelter | 403 | MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); | 
| 627 |  |  | } else { | 
| 628 | chrisfen | 417 |  | 
| 629 | gezelter | 403 | // Listen to your marching orders from processor 0: | 
| 630 | chrisfen | 417 |  | 
| 631 | gezelter | 403 | MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); | 
| 632 |  |  | } | 
| 633 | chrisfen | 417 |  | 
| 634 | gezelter | 403 | info->setMolToProcMap(molToProcMap); | 
| 635 |  |  | sprintf(checkPointMsg, | 
| 636 |  |  | "Successfully divided the molecules among the processors.\n"); | 
| 637 | gezelter | 1241 | errorCheckPoint(); | 
| 638 | gezelter | 403 | } | 
| 639 | chrisfen | 417 |  | 
| 640 | gezelter | 403 | #endif | 
| 641 | chrisfen | 417 |  | 
| 642 | gezelter | 403 | void SimCreator::createMolecules(SimInfo *info) { | 
| 643 |  |  | MoleculeCreator molCreator; | 
| 644 |  |  | int stampId; | 
| 645 | chrisfen | 417 |  | 
| 646 | gezelter | 403 | for(int i = 0; i < info->getNGlobalMolecules(); i++) { | 
| 647 | chrisfen | 417 |  | 
| 648 | gezelter | 403 | #ifdef IS_MPI | 
| 649 | chrisfen | 417 |  | 
| 650 | gezelter | 403 | if (info->getMolToProc(i) == worldRank) { | 
| 651 |  |  | #endif | 
| 652 | chrisfen | 417 |  | 
| 653 | gezelter | 403 | stampId = info->getMoleculeStampId(i); | 
| 654 | gezelter | 1600 | Molecule * mol = molCreator.createMolecule(info->getForceField(), | 
| 655 |  |  | info->getMoleculeStamp(stampId), | 
| 656 |  |  | stampId, i, | 
| 657 |  |  | info->getLocalIndexManager()); | 
| 658 | chrisfen | 417 |  | 
| 659 | gezelter | 403 | info->addMolecule(mol); | 
| 660 | chrisfen | 417 |  | 
| 661 | gezelter | 403 | #ifdef IS_MPI | 
| 662 | chrisfen | 417 |  | 
| 663 | gezelter | 403 | } | 
| 664 | chrisfen | 417 |  | 
| 665 | gezelter | 403 | #endif | 
| 666 | chrisfen | 417 |  | 
| 667 | gezelter | 403 | } //end for(int i=0) | 
| 668 |  |  | } | 
| 669 | chrisfen | 417 |  | 
| 670 | gezelter | 1710 | int SimCreator::computeStorageLayout(SimInfo* info) { | 
| 671 |  |  |  | 
| 672 | gezelter | 1711 | Globals* simParams = info->getSimParams(); | 
| 673 | gezelter | 1710 | int nRigidBodies = info->getNGlobalRigidBodies(); | 
| 674 |  |  | set<AtomType*> atomTypes = info->getSimulatedAtomTypes(); | 
| 675 |  |  | set<AtomType*>::iterator i; | 
| 676 |  |  | bool hasDirectionalAtoms = false; | 
| 677 |  |  | bool hasFixedCharge = false; | 
| 678 |  |  | bool hasMultipoles = false; | 
| 679 |  |  | bool hasPolarizable = false; | 
| 680 |  |  | bool hasFluctuatingCharge = false; | 
| 681 |  |  | bool hasMetallic = false; | 
| 682 |  |  | int storageLayout = 0; | 
| 683 |  |  | storageLayout |= DataStorage::dslPosition; | 
| 684 |  |  | storageLayout |= DataStorage::dslVelocity; | 
| 685 |  |  | storageLayout |= DataStorage::dslForce; | 
| 686 |  |  |  | 
| 687 |  |  | for (i = atomTypes.begin(); i != atomTypes.end(); ++i) { | 
| 688 |  |  |  | 
| 689 |  |  | DirectionalAdapter da = DirectionalAdapter( (*i) ); | 
| 690 |  |  | MultipoleAdapter ma = MultipoleAdapter( (*i) ); | 
| 691 |  |  | EAMAdapter ea = EAMAdapter( (*i) ); | 
| 692 |  |  | SuttonChenAdapter sca = SuttonChenAdapter( (*i) ); | 
| 693 |  |  | PolarizableAdapter pa = PolarizableAdapter( (*i) ); | 
| 694 |  |  | FixedChargeAdapter fca = FixedChargeAdapter( (*i) ); | 
| 695 |  |  | FluctuatingChargeAdapter fqa = FluctuatingChargeAdapter( (*i) ); | 
| 696 |  |  |  | 
| 697 |  |  | if (da.isDirectional()){ | 
| 698 |  |  | hasDirectionalAtoms = true; | 
| 699 |  |  | } | 
| 700 |  |  | if (ma.isMultipole()){ | 
| 701 |  |  | hasMultipoles = true; | 
| 702 |  |  | } | 
| 703 |  |  | if (ea.isEAM() || sca.isSuttonChen()){ | 
| 704 |  |  | hasMetallic = true; | 
| 705 |  |  | } | 
| 706 |  |  | if ( fca.isFixedCharge() ){ | 
| 707 |  |  | hasFixedCharge = true; | 
| 708 |  |  | } | 
| 709 |  |  | if ( fqa.isFluctuatingCharge() ){ | 
| 710 |  |  | hasFluctuatingCharge = true; | 
| 711 |  |  | } | 
| 712 |  |  | if ( pa.isPolarizable() ){ | 
| 713 |  |  | hasPolarizable = true; | 
| 714 |  |  | } | 
| 715 |  |  | } | 
| 716 |  |  |  | 
| 717 |  |  | if (nRigidBodies > 0 || hasDirectionalAtoms) { | 
| 718 |  |  | storageLayout |= DataStorage::dslAmat; | 
| 719 |  |  | if(storageLayout & DataStorage::dslVelocity) { | 
| 720 |  |  | storageLayout |= DataStorage::dslAngularMomentum; | 
| 721 |  |  | } | 
| 722 |  |  | if (storageLayout & DataStorage::dslForce) { | 
| 723 |  |  | storageLayout |= DataStorage::dslTorque; | 
| 724 |  |  | } | 
| 725 |  |  | } | 
| 726 |  |  | if (hasMultipoles) { | 
| 727 |  |  | storageLayout |= DataStorage::dslElectroFrame; | 
| 728 |  |  | } | 
| 729 |  |  | if (hasFixedCharge || hasFluctuatingCharge) { | 
| 730 |  |  | storageLayout |= DataStorage::dslSkippedCharge; | 
| 731 |  |  | } | 
| 732 |  |  | if (hasMetallic) { | 
| 733 |  |  | storageLayout |= DataStorage::dslDensity; | 
| 734 |  |  | storageLayout |= DataStorage::dslFunctional; | 
| 735 |  |  | storageLayout |= DataStorage::dslFunctionalDerivative; | 
| 736 |  |  | } | 
| 737 |  |  | if (hasPolarizable) { | 
| 738 |  |  | storageLayout |= DataStorage::dslElectricField; | 
| 739 |  |  | } | 
| 740 |  |  | if (hasFluctuatingCharge){ | 
| 741 |  |  | storageLayout |= DataStorage::dslFlucQPosition; | 
| 742 |  |  | if(storageLayout & DataStorage::dslVelocity) { | 
| 743 |  |  | storageLayout |= DataStorage::dslFlucQVelocity; | 
| 744 |  |  | } | 
| 745 |  |  | if (storageLayout & DataStorage::dslForce) { | 
| 746 |  |  | storageLayout |= DataStorage::dslFlucQForce; | 
| 747 |  |  | } | 
| 748 |  |  | } | 
| 749 | gezelter | 1711 |  | 
| 750 |  |  | if (simParams->getOutputParticlePotential()) { | 
| 751 |  |  | storageLayout |= DataStorage::dslParticlePot; | 
| 752 |  |  | } | 
| 753 |  |  |  | 
| 754 | gezelter | 1710 | return storageLayout; | 
| 755 |  |  | } | 
| 756 |  |  |  | 
| 757 | gezelter | 403 | void SimCreator::setGlobalIndex(SimInfo *info) { | 
| 758 |  |  | SimInfo::MoleculeIterator mi; | 
| 759 |  |  | Molecule::AtomIterator ai; | 
| 760 |  |  | Molecule::RigidBodyIterator ri; | 
| 761 |  |  | Molecule::CutoffGroupIterator ci; | 
| 762 | tim | 1024 | Molecule::IntegrableObjectIterator  ioi; | 
| 763 | gezelter | 403 | Molecule * mol; | 
| 764 |  |  | Atom * atom; | 
| 765 |  |  | RigidBody * rb; | 
| 766 |  |  | CutoffGroup * cg; | 
| 767 |  |  | int beginAtomIndex; | 
| 768 |  |  | int beginRigidBodyIndex; | 
| 769 |  |  | int beginCutoffGroupIndex; | 
| 770 |  |  | int nGlobalAtoms = info->getNGlobalAtoms(); | 
| 771 | chrisfen | 417 |  | 
| 772 | gezelter | 403 | beginAtomIndex = 0; | 
| 773 |  |  | beginRigidBodyIndex = 0; | 
| 774 |  |  | beginCutoffGroupIndex = 0; | 
| 775 | gezelter | 1600 |  | 
| 776 |  |  | for(int i = 0; i < info->getNGlobalMolecules(); i++) { | 
| 777 | chrisfen | 417 |  | 
| 778 | gezelter | 1600 | #ifdef IS_MPI | 
| 779 |  |  | if (info->getMolToProc(i) == worldRank) { | 
| 780 |  |  | #endif | 
| 781 |  |  | // stuff to do if I own this molecule | 
| 782 |  |  | mol = info->getMoleculeByGlobalIndex(i); | 
| 783 |  |  |  | 
| 784 |  |  | //local index(index in DataStorge) of atom is important | 
| 785 |  |  | for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { | 
| 786 |  |  | atom->setGlobalIndex(beginAtomIndex++); | 
| 787 |  |  | } | 
| 788 |  |  |  | 
| 789 |  |  | for(rb = mol->beginRigidBody(ri); rb != NULL; | 
| 790 |  |  | rb = mol->nextRigidBody(ri)) { | 
| 791 |  |  | rb->setGlobalIndex(beginRigidBodyIndex++); | 
| 792 |  |  | } | 
| 793 |  |  |  | 
| 794 |  |  | //local index of cutoff group is trivial, it only depends on | 
| 795 |  |  | //the order of travesing | 
| 796 |  |  | for(cg = mol->beginCutoffGroup(ci); cg != NULL; | 
| 797 |  |  | cg = mol->nextCutoffGroup(ci)) { | 
| 798 |  |  | cg->setGlobalIndex(beginCutoffGroupIndex++); | 
| 799 |  |  | } | 
| 800 |  |  |  | 
| 801 |  |  | #ifdef IS_MPI | 
| 802 |  |  | }  else { | 
| 803 |  |  |  | 
| 804 |  |  | // stuff to do if I don't own this molecule | 
| 805 |  |  |  | 
| 806 |  |  | int stampId = info->getMoleculeStampId(i); | 
| 807 |  |  | MoleculeStamp* stamp = info->getMoleculeStamp(stampId); | 
| 808 |  |  |  | 
| 809 |  |  | beginAtomIndex += stamp->getNAtoms(); | 
| 810 |  |  | beginRigidBodyIndex += stamp->getNRigidBodies(); | 
| 811 |  |  | beginCutoffGroupIndex += stamp->getNCutoffGroups() + stamp->getNFreeAtoms(); | 
| 812 | gezelter | 403 | } | 
| 813 | gezelter | 1600 | #endif | 
| 814 |  |  |  | 
| 815 |  |  | } //end for(int i=0) | 
| 816 |  |  |  | 
| 817 | gezelter | 403 | //fill globalGroupMembership | 
| 818 | gezelter | 1600 | std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0); | 
| 819 | gezelter | 403 | for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { | 
| 820 |  |  | for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { | 
| 821 | chrisfen | 417 |  | 
| 822 | gezelter | 403 | for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { | 
| 823 |  |  | globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex(); | 
| 824 |  |  | } | 
| 825 | chrisfen | 417 |  | 
| 826 | gezelter | 403 | } | 
| 827 |  |  | } | 
| 828 | gezelter | 1577 |  | 
| 829 | gezelter | 403 | #ifdef IS_MPI | 
| 830 |  |  | // Since the globalGroupMembership has been zero filled and we've only | 
| 831 |  |  | // poked values into the atoms we know, we can do an Allreduce | 
| 832 |  |  | // to get the full globalGroupMembership array (We think). | 
| 833 |  |  | // This would be prettier if we could use MPI_IN_PLACE like the MPI-2 | 
| 834 |  |  | // docs said we could. | 
| 835 | gezelter | 1313 | std::vector<int> tmpGroupMembership(info->getNGlobalAtoms(), 0); | 
| 836 | gezelter | 403 | MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms, | 
| 837 | gezelter | 1600 | MPI_INT, MPI_SUM, MPI_COMM_WORLD); | 
| 838 | gezelter | 403 | info->setGlobalGroupMembership(tmpGroupMembership); | 
| 839 |  |  | #else | 
| 840 |  |  | info->setGlobalGroupMembership(globalGroupMembership); | 
| 841 |  |  | #endif | 
| 842 | chrisfen | 417 |  | 
| 843 | gezelter | 403 | //fill molMembership | 
| 844 |  |  | std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0); | 
| 845 |  |  |  | 
| 846 |  |  | for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { | 
| 847 |  |  | for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { | 
| 848 |  |  | globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex(); | 
| 849 |  |  | } | 
| 850 |  |  | } | 
| 851 | chrisfen | 417 |  | 
| 852 | gezelter | 403 | #ifdef IS_MPI | 
| 853 | gezelter | 1313 | std::vector<int> tmpMolMembership(info->getNGlobalAtoms(), 0); | 
| 854 | chrisfen | 417 |  | 
| 855 | gezelter | 403 | MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms, | 
| 856 |  |  | MPI_INT, MPI_SUM, MPI_COMM_WORLD); | 
| 857 |  |  |  | 
| 858 |  |  | info->setGlobalMolMembership(tmpMolMembership); | 
| 859 |  |  | #else | 
| 860 |  |  | info->setGlobalMolMembership(globalMolMembership); | 
| 861 |  |  | #endif | 
| 862 | tim | 1024 |  | 
| 863 |  |  | // nIOPerMol holds the number of integrable objects per molecule | 
| 864 |  |  | // here the molecules are listed by their global indices. | 
| 865 |  |  |  | 
| 866 |  |  | std::vector<int> nIOPerMol(info->getNGlobalMolecules(), 0); | 
| 867 |  |  | for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { | 
| 868 |  |  | nIOPerMol[mol->getGlobalIndex()] = mol->getNIntegrableObjects(); | 
| 869 |  |  | } | 
| 870 | chrisfen | 417 |  | 
| 871 | tim | 1024 | #ifdef IS_MPI | 
| 872 |  |  | std::vector<int> numIntegrableObjectsPerMol(info->getNGlobalMolecules(), 0); | 
| 873 |  |  | MPI_Allreduce(&nIOPerMol[0], &numIntegrableObjectsPerMol[0], | 
| 874 |  |  | info->getNGlobalMolecules(), MPI_INT, MPI_SUM, MPI_COMM_WORLD); | 
| 875 |  |  | #else | 
| 876 |  |  | std::vector<int> numIntegrableObjectsPerMol = nIOPerMol; | 
| 877 |  |  | #endif | 
| 878 |  |  |  | 
| 879 | gezelter | 1313 | std::vector<int> startingIOIndexForMol(info->getNGlobalMolecules()); | 
| 880 |  |  |  | 
| 881 |  |  | int startingIndex = 0; | 
| 882 |  |  | for (int i = 0; i < info->getNGlobalMolecules(); i++) { | 
| 883 |  |  | startingIOIndexForMol[i] = startingIndex; | 
| 884 |  |  | startingIndex += numIntegrableObjectsPerMol[i]; | 
| 885 |  |  | } | 
| 886 |  |  |  | 
| 887 |  |  | std::vector<StuntDouble*> IOIndexToIntegrableObject(info->getNGlobalIntegrableObjects(), (StuntDouble*)NULL); | 
| 888 |  |  | for (mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { | 
| 889 | tim | 1024 | int myGlobalIndex = mol->getGlobalIndex(); | 
| 890 |  |  | int globalIO = startingIOIndexForMol[myGlobalIndex]; | 
| 891 |  |  | for (StuntDouble* integrableObject = mol->beginIntegrableObject(ioi); integrableObject != NULL; | 
| 892 |  |  | integrableObject = mol->nextIntegrableObject(ioi)) { | 
| 893 | gezelter | 1313 | integrableObject->setGlobalIntegrableObjectIndex(globalIO); | 
| 894 |  |  | IOIndexToIntegrableObject[globalIO] = integrableObject; | 
| 895 |  |  | globalIO++; | 
| 896 | tim | 1024 | } | 
| 897 |  |  | } | 
| 898 | gezelter | 1597 |  | 
| 899 | gezelter | 1313 | info->setIOIndexToIntegrableObject(IOIndexToIntegrableObject); | 
| 900 |  |  |  | 
| 901 | gezelter | 403 | } | 
| 902 | chrisfen | 417 |  | 
| 903 | tim | 1024 | void SimCreator::loadCoordinates(SimInfo* info, const std::string& mdFileName) { | 
| 904 | gezelter | 403 | Globals* simParams; | 
| 905 | gezelter | 1540 |  | 
| 906 | gezelter | 403 | simParams = info->getSimParams(); | 
| 907 |  |  |  | 
| 908 | tim | 1024 | DumpReader reader(info, mdFileName); | 
| 909 | gezelter | 403 | int nframes = reader.getNFrames(); | 
| 910 | gezelter | 1540 |  | 
| 911 | gezelter | 403 | if (nframes > 0) { | 
| 912 |  |  | reader.readFrame(nframes - 1); | 
| 913 |  |  | } else { | 
| 914 |  |  | //invalid initial coordinate file | 
| 915 | chrisfen | 417 | sprintf(painCave.errMsg, | 
| 916 |  |  | "Initial configuration file %s should at least contain one frame\n", | 
| 917 | tim | 1024 | mdFileName.c_str()); | 
| 918 | gezelter | 403 | painCave.isFatal = 1; | 
| 919 |  |  | simError(); | 
| 920 |  |  | } | 
| 921 |  |  | //copy the current snapshot to previous snapshot | 
| 922 |  |  | info->getSnapshotManager()->advance(); | 
| 923 |  |  | } | 
| 924 | chrisfen | 417 |  | 
| 925 | gezelter | 1390 | } //end namespace OpenMD | 
| 926 | gezelter | 403 |  | 
| 927 |  |  |  |