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 |
|
|
* 1. Acknowledgement of the program authors must be made in any |
10 |
|
|
* publication of scientific results based in part on use of the |
11 |
|
|
* program. An acceptable form of acknowledgement is citation of |
12 |
|
|
* the article in which the program was described (Matthew |
13 |
|
|
* A. Meineke, Charles F. Vardeman II, Teng Lin, Christopher |
14 |
|
|
* J. Fennell and J. Daniel Gezelter, "OOPSE: An Object-Oriented |
15 |
|
|
* Parallel Simulation Engine for Molecular Dynamics," |
16 |
|
|
* J. Comput. Chem. 26, pp. 252-271 (2005)) |
17 |
|
|
* |
18 |
|
|
* 2. Redistributions of source code must retain the above copyright |
19 |
|
|
* notice, this list of conditions and the following disclaimer. |
20 |
|
|
* |
21 |
|
|
* 3. Redistributions in binary form must reproduce the above copyright |
22 |
|
|
* notice, this list of conditions and the following disclaimer in the |
23 |
|
|
* documentation and/or other materials provided with the |
24 |
|
|
* distribution. |
25 |
|
|
* |
26 |
|
|
* This software is provided "AS IS," without a warranty of any |
27 |
|
|
* kind. All express or implied conditions, representations and |
28 |
|
|
* warranties, including any implied warranty of merchantability, |
29 |
|
|
* fitness for a particular purpose or non-infringement, are hereby |
30 |
|
|
* excluded. The University of Notre Dame and its licensors shall not |
31 |
|
|
* be liable for any damages suffered by licensee as a result of |
32 |
|
|
* using, modifying or distributing the software or its |
33 |
|
|
* derivatives. In no event will the University of Notre Dame or its |
34 |
|
|
* licensors be liable for any lost revenue, profit or data, or for |
35 |
|
|
* direct, indirect, special, consequential, incidental or punitive |
36 |
|
|
* damages, however caused and regardless of the theory of liability, |
37 |
|
|
* arising out of the use of or inability to use software, even if the |
38 |
|
|
* University of Notre Dame has been advised of the possibility of |
39 |
|
|
* such damages. |
40 |
|
|
*/ |
41 |
|
|
|
42 |
|
|
/** |
43 |
|
|
* @file SimCreator.cpp |
44 |
|
|
* @author tlin |
45 |
|
|
* @date 11/03/2004 |
46 |
|
|
* @time 13:51am |
47 |
|
|
* @version 1.0 |
48 |
|
|
*/ |
49 |
tim |
845 |
#include <exception> |
50 |
tim |
770 |
#include <iostream> |
51 |
|
|
#include <sstream> |
52 |
|
|
#include <string> |
53 |
|
|
|
54 |
gezelter |
403 |
#include "brains/MoleculeCreator.hpp" |
55 |
|
|
#include "brains/SimCreator.hpp" |
56 |
|
|
#include "brains/SimSnapshotManager.hpp" |
57 |
|
|
#include "io/DumpReader.hpp" |
58 |
|
|
#include "UseTheForce/ForceFieldFactory.hpp" |
59 |
|
|
#include "utils/simError.h" |
60 |
|
|
#include "utils/StringUtils.hpp" |
61 |
|
|
#include "math/SeqRandNumGen.hpp" |
62 |
tim |
770 |
#include "mdParser/MDLexer.hpp" |
63 |
|
|
#include "mdParser/MDParser.hpp" |
64 |
|
|
#include "mdParser/MDTreeParser.hpp" |
65 |
|
|
#include "mdParser/SimplePreprocessor.hpp" |
66 |
tim |
816 |
#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 |
tim |
770 |
|
73 |
tim |
816 |
#include "antlr/MismatchedCharException.hpp" |
74 |
|
|
#include "antlr/MismatchedTokenException.hpp" |
75 |
|
|
#include "antlr/NoViableAltForCharException.hpp" |
76 |
|
|
#include "antlr/NoViableAltException.hpp" |
77 |
tim |
770 |
|
78 |
gezelter |
403 |
#ifdef IS_MPI |
79 |
|
|
#include "math/ParallelRandNumGen.hpp" |
80 |
|
|
#endif |
81 |
|
|
|
82 |
|
|
namespace oopse { |
83 |
chrisfen |
417 |
|
84 |
tim |
770 |
Globals* SimCreator::parseFile(const std::string mdFileName){ |
85 |
|
|
Globals* simParams = NULL; |
86 |
|
|
try { |
87 |
|
|
|
88 |
|
|
// Create a preprocessor that preprocesses md file into an ostringstream |
89 |
|
|
std::stringstream ppStream; |
90 |
|
|
#ifdef IS_MPI |
91 |
|
|
int streamSize; |
92 |
|
|
const int masterNode = 0; |
93 |
|
|
int commStatus; |
94 |
|
|
if (worldRank == masterNode) { |
95 |
|
|
#endif |
96 |
|
|
|
97 |
|
|
SimplePreprocessor preprocessor; |
98 |
|
|
preprocessor.preprocess(mdFileName, ppStream); |
99 |
|
|
|
100 |
|
|
#ifdef IS_MPI |
101 |
|
|
//brocasting the stream size |
102 |
|
|
streamSize = ppStream.str().size() +1; |
103 |
|
|
commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD); |
104 |
|
|
|
105 |
tim |
832 |
commStatus = MPI_Bcast(static_cast<void*>(const_cast<char*>(ppStream.str().c_str())), streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
106 |
tim |
770 |
|
107 |
|
|
|
108 |
|
|
} else { |
109 |
|
|
//get stream size |
110 |
|
|
commStatus = MPI_Bcast(&streamSize, 1, MPI_LONG, masterNode, MPI_COMM_WORLD); |
111 |
|
|
|
112 |
|
|
char* buf = new char[streamSize]; |
113 |
|
|
assert(buf); |
114 |
|
|
|
115 |
|
|
//receive file content |
116 |
|
|
commStatus = MPI_Bcast(buf, streamSize, MPI_CHAR, masterNode, MPI_COMM_WORLD); |
117 |
|
|
|
118 |
|
|
ppStream.str(buf); |
119 |
|
|
delete buf; |
120 |
|
|
|
121 |
|
|
} |
122 |
|
|
#endif |
123 |
|
|
// Create a scanner that reads from the input stream |
124 |
|
|
MDLexer lexer(ppStream); |
125 |
|
|
lexer.setFilename(mdFileName); |
126 |
|
|
lexer.initDeferredLineCount(); |
127 |
chrisfen |
417 |
|
128 |
tim |
770 |
// Create a parser that reads from the scanner |
129 |
|
|
MDParser parser(lexer); |
130 |
|
|
parser.setFilename(mdFileName); |
131 |
|
|
|
132 |
|
|
// Create an observer that synchorizes file name change |
133 |
|
|
FilenameObserver observer; |
134 |
|
|
observer.setLexer(&lexer); |
135 |
|
|
observer.setParser(&parser); |
136 |
|
|
lexer.setObserver(&observer); |
137 |
chrisfen |
417 |
|
138 |
tim |
770 |
antlr::ASTFactory factory; |
139 |
|
|
parser.initializeASTFactory(factory); |
140 |
|
|
parser.setASTFactory(&factory); |
141 |
|
|
parser.mdfile(); |
142 |
|
|
|
143 |
|
|
// Create a tree parser that reads information into Globals |
144 |
|
|
MDTreeParser treeParser; |
145 |
|
|
treeParser.initializeASTFactory(factory); |
146 |
|
|
treeParser.setASTFactory(&factory); |
147 |
|
|
simParams = treeParser.walkTree(parser.getAST()); |
148 |
|
|
|
149 |
|
|
} |
150 |
tim |
845 |
|
151 |
tim |
816 |
|
152 |
|
|
catch(antlr::MismatchedCharException& e) { |
153 |
tim |
845 |
sprintf(painCave.errMsg, |
154 |
|
|
"parser exception: %s %s:%d:%d\n", |
155 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
156 |
|
|
painCave.isFatal = 1; |
157 |
|
|
simError(); |
158 |
tim |
816 |
} |
159 |
|
|
catch(antlr::MismatchedTokenException &e) { |
160 |
tim |
845 |
sprintf(painCave.errMsg, |
161 |
|
|
"parser exception: %s %s:%d:%d\n", |
162 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
163 |
|
|
painCave.isFatal = 1; |
164 |
|
|
simError(); |
165 |
tim |
816 |
} |
166 |
|
|
catch(antlr::NoViableAltForCharException &e) { |
167 |
tim |
845 |
sprintf(painCave.errMsg, |
168 |
|
|
"parser exception: %s %s:%d:%d\n", |
169 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
170 |
|
|
painCave.isFatal = 1; |
171 |
|
|
simError(); |
172 |
tim |
816 |
} |
173 |
|
|
catch(antlr::NoViableAltException &e) { |
174 |
tim |
845 |
sprintf(painCave.errMsg, |
175 |
|
|
"parser exception: %s %s:%d:%d\n", |
176 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
177 |
|
|
painCave.isFatal = 1; |
178 |
|
|
simError(); |
179 |
tim |
816 |
} |
180 |
tim |
845 |
|
181 |
tim |
816 |
catch(antlr::TokenStreamRecognitionException& e) { |
182 |
tim |
845 |
sprintf(painCave.errMsg, |
183 |
|
|
"parser exception: %s %s:%d:%d\n", |
184 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
185 |
|
|
painCave.isFatal = 1; |
186 |
|
|
simError(); |
187 |
tim |
816 |
} |
188 |
tim |
845 |
|
189 |
tim |
816 |
catch(antlr::TokenStreamIOException& e) { |
190 |
tim |
845 |
sprintf(painCave.errMsg, |
191 |
|
|
"parser exception: %s\n", |
192 |
|
|
e.getMessage().c_str()); |
193 |
|
|
painCave.isFatal = 1; |
194 |
|
|
simError(); |
195 |
tim |
816 |
} |
196 |
tim |
845 |
|
197 |
tim |
816 |
catch(antlr::TokenStreamException& e) { |
198 |
tim |
845 |
sprintf(painCave.errMsg, |
199 |
|
|
"parser exception: %s\n", |
200 |
|
|
e.getMessage().c_str()); |
201 |
|
|
painCave.isFatal = 1; |
202 |
|
|
simError(); |
203 |
tim |
816 |
} |
204 |
|
|
catch (antlr::RecognitionException& e) { |
205 |
tim |
845 |
sprintf(painCave.errMsg, |
206 |
|
|
"parser exception: %s %s:%d:%d\n", |
207 |
|
|
e.getMessage().c_str(),e.getFilename().c_str(), e.getLine(), e.getColumn()); |
208 |
|
|
painCave.isFatal = 1; |
209 |
|
|
simError(); |
210 |
tim |
816 |
} |
211 |
|
|
catch (antlr::CharStreamException& e) { |
212 |
tim |
845 |
sprintf(painCave.errMsg, |
213 |
|
|
"parser exception: %s\n", |
214 |
|
|
e.getMessage().c_str()); |
215 |
|
|
painCave.isFatal = 1; |
216 |
|
|
simError(); |
217 |
tim |
816 |
} |
218 |
tim |
845 |
catch (OOPSEException& e) { |
219 |
|
|
sprintf(painCave.errMsg, |
220 |
|
|
"%s\n", |
221 |
|
|
e.getMessage().c_str()); |
222 |
|
|
painCave.isFatal = 1; |
223 |
|
|
simError(); |
224 |
|
|
} |
225 |
|
|
catch (std::exception& e) { |
226 |
|
|
sprintf(painCave.errMsg, |
227 |
|
|
"parser exception: %s\n", |
228 |
|
|
e.what()); |
229 |
|
|
painCave.isFatal = 1; |
230 |
|
|
simError(); |
231 |
|
|
} |
232 |
tim |
770 |
|
233 |
|
|
return simParams; |
234 |
gezelter |
403 |
} |
235 |
chrisfen |
417 |
|
236 |
chrisfen |
514 |
SimInfo* SimCreator::createSim(const std::string & mdFileName, |
237 |
|
|
bool loadInitCoords) { |
238 |
tim |
770 |
|
239 |
gezelter |
403 |
//parse meta-data file |
240 |
tim |
770 |
Globals* simParams = parseFile(mdFileName); |
241 |
chrisfen |
417 |
|
242 |
gezelter |
403 |
//create the force field |
243 |
chrisfen |
514 |
ForceField * ff = ForceFieldFactory::getInstance() |
244 |
|
|
->createForceField(simParams->getForceField()); |
245 |
gezelter |
403 |
|
246 |
|
|
if (ff == NULL) { |
247 |
chrisfen |
514 |
sprintf(painCave.errMsg, |
248 |
|
|
"ForceField Factory can not create %s force field\n", |
249 |
tim |
665 |
simParams->getForceField().c_str()); |
250 |
gezelter |
403 |
painCave.isFatal = 1; |
251 |
|
|
simError(); |
252 |
|
|
} |
253 |
chrisfen |
417 |
|
254 |
gezelter |
403 |
if (simParams->haveForceFieldFileName()) { |
255 |
|
|
ff->setForceFieldFileName(simParams->getForceFieldFileName()); |
256 |
|
|
} |
257 |
|
|
|
258 |
|
|
std::string forcefieldFileName; |
259 |
|
|
forcefieldFileName = ff->getForceFieldFileName(); |
260 |
chrisfen |
417 |
|
261 |
gezelter |
403 |
if (simParams->haveForceFieldVariant()) { |
262 |
|
|
//If the force field has variant, the variant force field name will be |
263 |
|
|
//Base.variant.frc. For exampel EAM.u6.frc |
264 |
chrisfen |
417 |
|
265 |
gezelter |
403 |
std::string variant = simParams->getForceFieldVariant(); |
266 |
chrisfen |
417 |
|
267 |
gezelter |
403 |
std::string::size_type pos = forcefieldFileName.rfind(".frc"); |
268 |
|
|
variant = "." + variant; |
269 |
|
|
if (pos != std::string::npos) { |
270 |
|
|
forcefieldFileName.insert(pos, variant); |
271 |
|
|
} else { |
272 |
|
|
//If the default force field file name does not containt .frc suffix, just append the .variant |
273 |
|
|
forcefieldFileName.append(variant); |
274 |
|
|
} |
275 |
|
|
} |
276 |
|
|
|
277 |
|
|
ff->parse(forcefieldFileName); |
278 |
chuckv |
823 |
ff->setFortranForceOptions(); |
279 |
gezelter |
403 |
//create SimInfo |
280 |
tim |
770 |
SimInfo * info = new SimInfo(ff, simParams); |
281 |
tim |
490 |
|
282 |
gezelter |
403 |
//gather parameters (SimCreator only retrieves part of the parameters) |
283 |
|
|
gatherParameters(info, mdFileName); |
284 |
chrisfen |
417 |
|
285 |
gezelter |
403 |
//divide the molecules and determine the global index of molecules |
286 |
|
|
#ifdef IS_MPI |
287 |
|
|
divideMolecules(info); |
288 |
|
|
#endif |
289 |
chrisfen |
417 |
|
290 |
gezelter |
403 |
//create the molecules |
291 |
|
|
createMolecules(info); |
292 |
chrisfen |
417 |
|
293 |
|
|
|
294 |
gezelter |
403 |
//allocate memory for DataStorage(circular reference, need to break it) |
295 |
|
|
info->setSnapshotManager(new SimSnapshotManager(info)); |
296 |
|
|
|
297 |
|
|
//set the global index of atoms, rigidbodies and cutoffgroups (only need to be set once, the |
298 |
|
|
//global index will never change again). Local indices of atoms and rigidbodies are already set by |
299 |
|
|
//MoleculeCreator class which actually delegates the responsibility to LocalIndexManager. |
300 |
|
|
setGlobalIndex(info); |
301 |
chrisfen |
417 |
|
302 |
gezelter |
403 |
//Alought addExculdePairs is called inside SimInfo's addMolecule method, at that point |
303 |
|
|
//atoms don't have the global index yet (their global index are all initialized to -1). |
304 |
|
|
//Therefore we have to call addExcludePairs explicitly here. A way to work around is that |
305 |
|
|
//we can determine the beginning global indices of atoms before they get created. |
306 |
|
|
SimInfo::MoleculeIterator mi; |
307 |
|
|
Molecule* mol; |
308 |
|
|
for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
309 |
|
|
info->addExcludePairs(mol); |
310 |
|
|
} |
311 |
|
|
|
312 |
|
|
if (loadInitCoords) |
313 |
|
|
loadCoordinates(info); |
314 |
|
|
|
315 |
|
|
return info; |
316 |
|
|
} |
317 |
chrisfen |
417 |
|
318 |
gezelter |
403 |
void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) { |
319 |
chrisfen |
417 |
|
320 |
tim |
749 |
//figure out the output file names |
321 |
gezelter |
403 |
std::string prefix; |
322 |
chrisfen |
417 |
|
323 |
gezelter |
403 |
#ifdef IS_MPI |
324 |
chrisfen |
417 |
|
325 |
gezelter |
403 |
if (worldRank == 0) { |
326 |
|
|
#endif // is_mpi |
327 |
|
|
Globals * simParams = info->getSimParams(); |
328 |
|
|
if (simParams->haveFinalConfig()) { |
329 |
|
|
prefix = getPrefix(simParams->getFinalConfig()); |
330 |
|
|
} else { |
331 |
|
|
prefix = getPrefix(mdfile); |
332 |
|
|
} |
333 |
chrisfen |
417 |
|
334 |
gezelter |
403 |
info->setFinalConfigFileName(prefix + ".eor"); |
335 |
|
|
info->setDumpFileName(prefix + ".dump"); |
336 |
|
|
info->setStatFileName(prefix + ".stat"); |
337 |
chrisfen |
417 |
info->setRestFileName(prefix + ".zang"); |
338 |
|
|
|
339 |
gezelter |
403 |
#ifdef IS_MPI |
340 |
chrisfen |
417 |
|
341 |
gezelter |
403 |
} |
342 |
chrisfen |
417 |
|
343 |
gezelter |
403 |
#endif |
344 |
chrisfen |
417 |
|
345 |
gezelter |
403 |
} |
346 |
chrisfen |
417 |
|
347 |
gezelter |
403 |
#ifdef IS_MPI |
348 |
|
|
void SimCreator::divideMolecules(SimInfo *info) { |
349 |
|
|
double numerator; |
350 |
|
|
double denominator; |
351 |
|
|
double precast; |
352 |
|
|
double x; |
353 |
|
|
double y; |
354 |
|
|
double a; |
355 |
|
|
int old_atoms; |
356 |
|
|
int add_atoms; |
357 |
|
|
int new_atoms; |
358 |
|
|
int nTarget; |
359 |
|
|
int done; |
360 |
|
|
int i; |
361 |
|
|
int j; |
362 |
|
|
int loops; |
363 |
|
|
int which_proc; |
364 |
|
|
int nProcessors; |
365 |
|
|
std::vector<int> atomsPerProc; |
366 |
|
|
int nGlobalMols = info->getNGlobalMolecules(); |
367 |
|
|
std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition: |
368 |
|
|
|
369 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &nProcessors); |
370 |
chrisfen |
417 |
|
371 |
gezelter |
403 |
if (nProcessors > nGlobalMols) { |
372 |
|
|
sprintf(painCave.errMsg, |
373 |
|
|
"nProcessors (%d) > nMol (%d)\n" |
374 |
|
|
"\tThe number of processors is larger than\n" |
375 |
|
|
"\tthe number of molecules. This will not result in a \n" |
376 |
|
|
"\tusable division of atoms for force decomposition.\n" |
377 |
|
|
"\tEither try a smaller number of processors, or run the\n" |
378 |
|
|
"\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols); |
379 |
chrisfen |
417 |
|
380 |
gezelter |
403 |
painCave.isFatal = 1; |
381 |
|
|
simError(); |
382 |
|
|
} |
383 |
chrisfen |
417 |
|
384 |
gezelter |
403 |
int seedValue; |
385 |
|
|
Globals * simParams = info->getSimParams(); |
386 |
|
|
SeqRandNumGen* myRandom; //divide labor does not need Parallel random number generator |
387 |
|
|
if (simParams->haveSeed()) { |
388 |
|
|
seedValue = simParams->getSeed(); |
389 |
|
|
myRandom = new SeqRandNumGen(seedValue); |
390 |
|
|
}else { |
391 |
|
|
myRandom = new SeqRandNumGen(); |
392 |
|
|
} |
393 |
chrisfen |
417 |
|
394 |
|
|
|
395 |
gezelter |
403 |
a = 3.0 * nGlobalMols / info->getNGlobalAtoms(); |
396 |
chrisfen |
417 |
|
397 |
gezelter |
403 |
//initialize atomsPerProc |
398 |
|
|
atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0); |
399 |
chrisfen |
417 |
|
400 |
gezelter |
403 |
if (worldRank == 0) { |
401 |
|
|
numerator = info->getNGlobalAtoms(); |
402 |
|
|
denominator = nProcessors; |
403 |
|
|
precast = numerator / denominator; |
404 |
|
|
nTarget = (int)(precast + 0.5); |
405 |
chrisfen |
417 |
|
406 |
gezelter |
403 |
for(i = 0; i < nGlobalMols; i++) { |
407 |
|
|
done = 0; |
408 |
|
|
loops = 0; |
409 |
chrisfen |
417 |
|
410 |
gezelter |
403 |
while (!done) { |
411 |
|
|
loops++; |
412 |
chrisfen |
417 |
|
413 |
gezelter |
403 |
// Pick a processor at random |
414 |
chrisfen |
417 |
|
415 |
gezelter |
403 |
which_proc = (int) (myRandom->rand() * nProcessors); |
416 |
chrisfen |
417 |
|
417 |
gezelter |
403 |
//get the molecule stamp first |
418 |
|
|
int stampId = info->getMoleculeStampId(i); |
419 |
|
|
MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId); |
420 |
chrisfen |
417 |
|
421 |
gezelter |
403 |
// How many atoms does this processor have so far? |
422 |
|
|
old_atoms = atomsPerProc[which_proc]; |
423 |
|
|
add_atoms = moleculeStamp->getNAtoms(); |
424 |
|
|
new_atoms = old_atoms + add_atoms; |
425 |
chrisfen |
417 |
|
426 |
gezelter |
403 |
// If we've been through this loop too many times, we need |
427 |
|
|
// to just give up and assign the molecule to this processor |
428 |
|
|
// and be done with it. |
429 |
chrisfen |
417 |
|
430 |
gezelter |
403 |
if (loops > 100) { |
431 |
|
|
sprintf(painCave.errMsg, |
432 |
|
|
"I've tried 100 times to assign molecule %d to a " |
433 |
|
|
" processor, but can't find a good spot.\n" |
434 |
|
|
"I'm assigning it at random to processor %d.\n", |
435 |
|
|
i, which_proc); |
436 |
chrisfen |
417 |
|
437 |
gezelter |
403 |
painCave.isFatal = 0; |
438 |
|
|
simError(); |
439 |
chrisfen |
417 |
|
440 |
gezelter |
403 |
molToProcMap[i] = which_proc; |
441 |
|
|
atomsPerProc[which_proc] += add_atoms; |
442 |
chrisfen |
417 |
|
443 |
gezelter |
403 |
done = 1; |
444 |
|
|
continue; |
445 |
|
|
} |
446 |
chrisfen |
417 |
|
447 |
gezelter |
403 |
// If we can add this molecule to this processor without sending |
448 |
|
|
// it above nTarget, then go ahead and do it: |
449 |
chrisfen |
417 |
|
450 |
gezelter |
403 |
if (new_atoms <= nTarget) { |
451 |
|
|
molToProcMap[i] = which_proc; |
452 |
|
|
atomsPerProc[which_proc] += add_atoms; |
453 |
chrisfen |
417 |
|
454 |
gezelter |
403 |
done = 1; |
455 |
|
|
continue; |
456 |
|
|
} |
457 |
chrisfen |
417 |
|
458 |
gezelter |
403 |
// The only situation left is when new_atoms > nTarget. We |
459 |
|
|
// want to accept this with some probability that dies off the |
460 |
|
|
// farther we are from nTarget |
461 |
chrisfen |
417 |
|
462 |
gezelter |
403 |
// roughly: x = new_atoms - nTarget |
463 |
|
|
// Pacc(x) = exp(- a * x) |
464 |
|
|
// where a = penalty / (average atoms per molecule) |
465 |
chrisfen |
417 |
|
466 |
gezelter |
403 |
x = (double)(new_atoms - nTarget); |
467 |
|
|
y = myRandom->rand(); |
468 |
chrisfen |
417 |
|
469 |
gezelter |
403 |
if (y < exp(- a * x)) { |
470 |
|
|
molToProcMap[i] = which_proc; |
471 |
|
|
atomsPerProc[which_proc] += add_atoms; |
472 |
chrisfen |
417 |
|
473 |
gezelter |
403 |
done = 1; |
474 |
|
|
continue; |
475 |
|
|
} else { |
476 |
|
|
continue; |
477 |
|
|
} |
478 |
|
|
} |
479 |
|
|
} |
480 |
chrisfen |
417 |
|
481 |
gezelter |
403 |
delete myRandom; |
482 |
chrisfen |
417 |
|
483 |
gezelter |
403 |
// Spray out this nonsense to all other processors: |
484 |
chrisfen |
417 |
|
485 |
gezelter |
403 |
MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); |
486 |
|
|
} else { |
487 |
chrisfen |
417 |
|
488 |
gezelter |
403 |
// Listen to your marching orders from processor 0: |
489 |
chrisfen |
417 |
|
490 |
gezelter |
403 |
MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); |
491 |
|
|
} |
492 |
chrisfen |
417 |
|
493 |
gezelter |
403 |
info->setMolToProcMap(molToProcMap); |
494 |
|
|
sprintf(checkPointMsg, |
495 |
|
|
"Successfully divided the molecules among the processors.\n"); |
496 |
|
|
MPIcheckPoint(); |
497 |
|
|
} |
498 |
chrisfen |
417 |
|
499 |
gezelter |
403 |
#endif |
500 |
chrisfen |
417 |
|
501 |
gezelter |
403 |
void SimCreator::createMolecules(SimInfo *info) { |
502 |
|
|
MoleculeCreator molCreator; |
503 |
|
|
int stampId; |
504 |
chrisfen |
417 |
|
505 |
gezelter |
403 |
for(int i = 0; i < info->getNGlobalMolecules(); i++) { |
506 |
chrisfen |
417 |
|
507 |
gezelter |
403 |
#ifdef IS_MPI |
508 |
chrisfen |
417 |
|
509 |
gezelter |
403 |
if (info->getMolToProc(i) == worldRank) { |
510 |
|
|
#endif |
511 |
chrisfen |
417 |
|
512 |
gezelter |
403 |
stampId = info->getMoleculeStampId(i); |
513 |
|
|
Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId), |
514 |
|
|
stampId, i, info->getLocalIndexManager()); |
515 |
chrisfen |
417 |
|
516 |
gezelter |
403 |
info->addMolecule(mol); |
517 |
chrisfen |
417 |
|
518 |
gezelter |
403 |
#ifdef IS_MPI |
519 |
chrisfen |
417 |
|
520 |
gezelter |
403 |
} |
521 |
chrisfen |
417 |
|
522 |
gezelter |
403 |
#endif |
523 |
chrisfen |
417 |
|
524 |
gezelter |
403 |
} //end for(int i=0) |
525 |
|
|
} |
526 |
chrisfen |
417 |
|
527 |
gezelter |
403 |
void SimCreator::setGlobalIndex(SimInfo *info) { |
528 |
|
|
SimInfo::MoleculeIterator mi; |
529 |
|
|
Molecule::AtomIterator ai; |
530 |
|
|
Molecule::RigidBodyIterator ri; |
531 |
|
|
Molecule::CutoffGroupIterator ci; |
532 |
|
|
Molecule * mol; |
533 |
|
|
Atom * atom; |
534 |
|
|
RigidBody * rb; |
535 |
|
|
CutoffGroup * cg; |
536 |
|
|
int beginAtomIndex; |
537 |
|
|
int beginRigidBodyIndex; |
538 |
|
|
int beginCutoffGroupIndex; |
539 |
|
|
int nGlobalAtoms = info->getNGlobalAtoms(); |
540 |
|
|
|
541 |
|
|
#ifndef IS_MPI |
542 |
chrisfen |
417 |
|
543 |
gezelter |
403 |
beginAtomIndex = 0; |
544 |
|
|
beginRigidBodyIndex = 0; |
545 |
|
|
beginCutoffGroupIndex = 0; |
546 |
chrisfen |
417 |
|
547 |
gezelter |
403 |
#else |
548 |
chrisfen |
417 |
|
549 |
gezelter |
403 |
int nproc; |
550 |
|
|
int myNode; |
551 |
chrisfen |
417 |
|
552 |
gezelter |
403 |
myNode = worldRank; |
553 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &nproc); |
554 |
chrisfen |
417 |
|
555 |
gezelter |
403 |
std::vector < int > tmpAtomsInProc(nproc, 0); |
556 |
|
|
std::vector < int > tmpRigidBodiesInProc(nproc, 0); |
557 |
|
|
std::vector < int > tmpCutoffGroupsInProc(nproc, 0); |
558 |
|
|
std::vector < int > NumAtomsInProc(nproc, 0); |
559 |
|
|
std::vector < int > NumRigidBodiesInProc(nproc, 0); |
560 |
|
|
std::vector < int > NumCutoffGroupsInProc(nproc, 0); |
561 |
chrisfen |
417 |
|
562 |
gezelter |
403 |
tmpAtomsInProc[myNode] = info->getNAtoms(); |
563 |
|
|
tmpRigidBodiesInProc[myNode] = info->getNRigidBodies(); |
564 |
|
|
tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups(); |
565 |
chrisfen |
417 |
|
566 |
gezelter |
403 |
//do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups |
567 |
|
|
MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT, |
568 |
|
|
MPI_SUM, MPI_COMM_WORLD); |
569 |
|
|
MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc, |
570 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
571 |
|
|
MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc, |
572 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
573 |
chrisfen |
417 |
|
574 |
gezelter |
403 |
beginAtomIndex = 0; |
575 |
|
|
beginRigidBodyIndex = 0; |
576 |
|
|
beginCutoffGroupIndex = 0; |
577 |
chrisfen |
417 |
|
578 |
gezelter |
403 |
for(int i = 0; i < myNode; i++) { |
579 |
|
|
beginAtomIndex += NumAtomsInProc[i]; |
580 |
|
|
beginRigidBodyIndex += NumRigidBodiesInProc[i]; |
581 |
|
|
beginCutoffGroupIndex += NumCutoffGroupsInProc[i]; |
582 |
|
|
} |
583 |
chrisfen |
417 |
|
584 |
gezelter |
403 |
#endif |
585 |
chrisfen |
417 |
|
586 |
gezelter |
403 |
//rigidbody's index begins right after atom's |
587 |
|
|
beginRigidBodyIndex += info->getNGlobalAtoms(); |
588 |
chrisfen |
417 |
|
589 |
gezelter |
403 |
for(mol = info->beginMolecule(mi); mol != NULL; |
590 |
|
|
mol = info->nextMolecule(mi)) { |
591 |
chrisfen |
417 |
|
592 |
gezelter |
403 |
//local index(index in DataStorge) of atom is important |
593 |
|
|
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
594 |
|
|
atom->setGlobalIndex(beginAtomIndex++); |
595 |
|
|
} |
596 |
chrisfen |
417 |
|
597 |
gezelter |
403 |
for(rb = mol->beginRigidBody(ri); rb != NULL; |
598 |
|
|
rb = mol->nextRigidBody(ri)) { |
599 |
|
|
rb->setGlobalIndex(beginRigidBodyIndex++); |
600 |
|
|
} |
601 |
chrisfen |
417 |
|
602 |
gezelter |
403 |
//local index of cutoff group is trivial, it only depends on the order of travesing |
603 |
|
|
for(cg = mol->beginCutoffGroup(ci); cg != NULL; |
604 |
|
|
cg = mol->nextCutoffGroup(ci)) { |
605 |
|
|
cg->setGlobalIndex(beginCutoffGroupIndex++); |
606 |
|
|
} |
607 |
|
|
} |
608 |
chrisfen |
417 |
|
609 |
gezelter |
403 |
//fill globalGroupMembership |
610 |
|
|
std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0); |
611 |
|
|
for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
612 |
|
|
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
613 |
chrisfen |
417 |
|
614 |
gezelter |
403 |
for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { |
615 |
|
|
globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex(); |
616 |
|
|
} |
617 |
chrisfen |
417 |
|
618 |
gezelter |
403 |
} |
619 |
|
|
} |
620 |
chrisfen |
417 |
|
621 |
gezelter |
403 |
#ifdef IS_MPI |
622 |
|
|
// Since the globalGroupMembership has been zero filled and we've only |
623 |
|
|
// poked values into the atoms we know, we can do an Allreduce |
624 |
|
|
// to get the full globalGroupMembership array (We think). |
625 |
|
|
// This would be prettier if we could use MPI_IN_PLACE like the MPI-2 |
626 |
|
|
// docs said we could. |
627 |
|
|
std::vector<int> tmpGroupMembership(nGlobalAtoms, 0); |
628 |
|
|
MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms, |
629 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
630 |
|
|
info->setGlobalGroupMembership(tmpGroupMembership); |
631 |
|
|
#else |
632 |
|
|
info->setGlobalGroupMembership(globalGroupMembership); |
633 |
|
|
#endif |
634 |
chrisfen |
417 |
|
635 |
gezelter |
403 |
//fill molMembership |
636 |
|
|
std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0); |
637 |
|
|
|
638 |
|
|
for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
639 |
chrisfen |
417 |
|
640 |
gezelter |
403 |
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
641 |
|
|
globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex(); |
642 |
|
|
} |
643 |
|
|
} |
644 |
chrisfen |
417 |
|
645 |
gezelter |
403 |
#ifdef IS_MPI |
646 |
|
|
std::vector<int> tmpMolMembership(nGlobalAtoms, 0); |
647 |
chrisfen |
417 |
|
648 |
gezelter |
403 |
MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms, |
649 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
650 |
|
|
|
651 |
|
|
info->setGlobalMolMembership(tmpMolMembership); |
652 |
|
|
#else |
653 |
|
|
info->setGlobalMolMembership(globalMolMembership); |
654 |
|
|
#endif |
655 |
chrisfen |
417 |
|
656 |
gezelter |
403 |
} |
657 |
chrisfen |
417 |
|
658 |
gezelter |
403 |
void SimCreator::loadCoordinates(SimInfo* info) { |
659 |
|
|
Globals* simParams; |
660 |
|
|
simParams = info->getSimParams(); |
661 |
|
|
|
662 |
|
|
if (!simParams->haveInitialConfig()) { |
663 |
|
|
sprintf(painCave.errMsg, |
664 |
|
|
"Cannot intialize a simulation without an initial configuration file.\n"); |
665 |
|
|
painCave.isFatal = 1;; |
666 |
|
|
simError(); |
667 |
|
|
} |
668 |
chrisfen |
417 |
|
669 |
gezelter |
403 |
DumpReader reader(info, simParams->getInitialConfig()); |
670 |
|
|
int nframes = reader.getNFrames(); |
671 |
chrisfen |
417 |
|
672 |
gezelter |
403 |
if (nframes > 0) { |
673 |
|
|
reader.readFrame(nframes - 1); |
674 |
|
|
} else { |
675 |
|
|
//invalid initial coordinate file |
676 |
chrisfen |
417 |
sprintf(painCave.errMsg, |
677 |
|
|
"Initial configuration file %s should at least contain one frame\n", |
678 |
tim |
665 |
simParams->getInitialConfig().c_str()); |
679 |
gezelter |
403 |
painCave.isFatal = 1; |
680 |
|
|
simError(); |
681 |
|
|
} |
682 |
chrisfen |
417 |
|
683 |
gezelter |
403 |
//copy the current snapshot to previous snapshot |
684 |
|
|
info->getSnapshotManager()->advance(); |
685 |
|
|
} |
686 |
chrisfen |
417 |
|
687 |
gezelter |
403 |
} //end namespace oopse |
688 |
|
|
|
689 |
|
|
|