1 |
gezelter |
246 |
/* |
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 |
|
|
|
50 |
|
|
#include <sprng.h> |
51 |
|
|
|
52 |
|
|
#include "brains/MoleculeCreator.hpp" |
53 |
|
|
#include "brains/SimCreator.hpp" |
54 |
|
|
#include "brains/SimSnapshotManager.hpp" |
55 |
|
|
#include "io/DumpReader.hpp" |
56 |
|
|
#include "io/parse_me.h" |
57 |
|
|
#include "UseTheForce/ForceFieldFactory.hpp" |
58 |
|
|
#include "utils/simError.h" |
59 |
|
|
#include "utils/StringUtils.hpp" |
60 |
|
|
#ifdef IS_MPI |
61 |
|
|
#include "io/mpiBASS.h" |
62 |
|
|
#include "math/randomSPRNG.hpp" |
63 |
|
|
#endif |
64 |
|
|
|
65 |
|
|
namespace oopse { |
66 |
|
|
|
67 |
|
|
void SimCreator::parseFile(const std::string mdFileName, MakeStamps* stamps, Globals* simParams){ |
68 |
|
|
|
69 |
|
|
#ifdef IS_MPI |
70 |
|
|
|
71 |
|
|
if (worldRank == 0) { |
72 |
|
|
#endif // is_mpi |
73 |
|
|
|
74 |
|
|
simParams->initalize(); |
75 |
|
|
set_interface_stamps(stamps, simParams); |
76 |
|
|
|
77 |
|
|
#ifdef IS_MPI |
78 |
|
|
|
79 |
|
|
mpiEventInit(); |
80 |
|
|
|
81 |
|
|
#endif |
82 |
|
|
|
83 |
|
|
yacc_BASS(mdFileName.c_str()); |
84 |
|
|
|
85 |
|
|
#ifdef IS_MPI |
86 |
|
|
|
87 |
|
|
throwMPIEvent(NULL); |
88 |
|
|
} else { |
89 |
|
|
set_interface_stamps(stamps, simParams); |
90 |
|
|
mpiEventInit(); |
91 |
|
|
MPIcheckPoint(); |
92 |
|
|
mpiEventLoop(); |
93 |
|
|
} |
94 |
|
|
|
95 |
|
|
#endif |
96 |
|
|
|
97 |
|
|
} |
98 |
|
|
|
99 |
|
|
SimInfo* SimCreator::createSim(const std::string & mdFileName, bool loadInitCoords) { |
100 |
|
|
|
101 |
|
|
MakeStamps * stamps = new MakeStamps(); |
102 |
|
|
|
103 |
|
|
Globals * simParams = new Globals(); |
104 |
|
|
|
105 |
|
|
//parse meta-data file |
106 |
|
|
parseFile(mdFileName, stamps, simParams); |
107 |
|
|
|
108 |
|
|
//create the force field |
109 |
|
|
ForceField * ff = ForceFieldFactory::getInstance()->createForceField( |
110 |
|
|
simParams->getForceField()); |
111 |
|
|
|
112 |
|
|
if (ff == NULL) { |
113 |
|
|
sprintf(painCave.errMsg, "ForceField Factory can not create %s force field\n", |
114 |
|
|
simParams->getForceField()); |
115 |
|
|
painCave.isFatal = 1; |
116 |
|
|
simError(); |
117 |
|
|
} |
118 |
|
|
|
119 |
tim |
273 |
if (simParams->haveForceFieldFileName()) { |
120 |
|
|
ff->setForceFieldFileName(simParams->getForceFieldFileName()); |
121 |
|
|
} |
122 |
|
|
|
123 |
gezelter |
246 |
std::string forcefieldFileName; |
124 |
|
|
forcefieldFileName = ff->getForceFieldFileName(); |
125 |
|
|
|
126 |
|
|
if (simParams->haveForceFieldVariant()) { |
127 |
|
|
//If the force field has variant, the variant force field name will be |
128 |
|
|
//Base.variant.frc. For exampel EAM.u6.frc |
129 |
|
|
|
130 |
|
|
std::string variant = simParams->getForceFieldVariant(); |
131 |
|
|
|
132 |
|
|
std::string::size_type pos = forcefieldFileName.rfind(".frc"); |
133 |
|
|
variant = "." + variant; |
134 |
|
|
if (pos != std::string::npos) { |
135 |
|
|
forcefieldFileName.insert(pos, variant); |
136 |
|
|
} else { |
137 |
|
|
//If the default force field file name does not containt .frc suffix, just append the .variant |
138 |
|
|
forcefieldFileName.append(variant); |
139 |
|
|
} |
140 |
|
|
} |
141 |
|
|
|
142 |
|
|
ff->parse(forcefieldFileName); |
143 |
|
|
|
144 |
|
|
//extract the molecule stamps |
145 |
|
|
std::vector < std::pair<MoleculeStamp *, int> > moleculeStampPairs; |
146 |
|
|
compList(stamps, simParams, moleculeStampPairs); |
147 |
|
|
|
148 |
|
|
//create SimInfo |
149 |
|
|
SimInfo * info = new SimInfo(moleculeStampPairs, ff, simParams); |
150 |
|
|
|
151 |
|
|
//gather parameters (SimCreator only retrieves part of the parameters) |
152 |
|
|
gatherParameters(info, mdFileName); |
153 |
|
|
|
154 |
|
|
//divide the molecules and determine the global index of molecules |
155 |
|
|
#ifdef IS_MPI |
156 |
|
|
divideMolecules(info); |
157 |
|
|
#endif |
158 |
|
|
|
159 |
|
|
//create the molecules |
160 |
|
|
createMolecules(info); |
161 |
|
|
|
162 |
|
|
|
163 |
|
|
//allocate memory for DataStorage(circular reference, need to break it) |
164 |
|
|
info->setSnapshotManager(new SimSnapshotManager(info)); |
165 |
|
|
|
166 |
|
|
//set the global index of atoms, rigidbodies and cutoffgroups (only need to be set once, the |
167 |
|
|
//global index will never change again). Local indices of atoms and rigidbodies are already set by |
168 |
|
|
//MoleculeCreator class which actually delegates the responsibility to LocalIndexManager. |
169 |
|
|
setGlobalIndex(info); |
170 |
|
|
|
171 |
|
|
//Alought addExculdePairs is called inside SimInfo's addMolecule method, at that point |
172 |
|
|
//atoms don't have the global index yet (their global index are all initialized to -1). |
173 |
|
|
//Therefore we have to call addExcludePairs explicitly here. A way to work around is that |
174 |
|
|
//we can determine the beginning global indices of atoms before they get created. |
175 |
|
|
SimInfo::MoleculeIterator mi; |
176 |
|
|
Molecule* mol; |
177 |
|
|
for (mol= info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
178 |
|
|
info->addExcludePairs(mol); |
179 |
|
|
} |
180 |
|
|
|
181 |
|
|
|
182 |
|
|
//load initial coordinates, some extra information are pushed into SimInfo's property map ( such as |
183 |
|
|
//eta, chi for NPT integrator) |
184 |
|
|
if (loadInitCoords) |
185 |
|
|
loadCoordinates(info); |
186 |
|
|
|
187 |
|
|
return info; |
188 |
|
|
} |
189 |
|
|
|
190 |
|
|
void SimCreator::gatherParameters(SimInfo *info, const std::string& mdfile) { |
191 |
|
|
|
192 |
|
|
//setup seed for random number generator |
193 |
|
|
int seedValue; |
194 |
|
|
Globals * simParams = info->getSimParams(); |
195 |
|
|
|
196 |
|
|
if (simParams->haveSeed()) { |
197 |
|
|
seedValue = simParams->getSeed(); |
198 |
|
|
|
199 |
|
|
if (seedValue < 100000000 ) { |
200 |
|
|
sprintf(painCave.errMsg, |
201 |
|
|
"Seed for sprng library should contain at least 9 digits\n" |
202 |
|
|
"OOPSE will generate a seed for user\n"); |
203 |
|
|
|
204 |
|
|
painCave.isFatal = 0; |
205 |
|
|
simError(); |
206 |
|
|
|
207 |
|
|
//using seed generated by system instead of invalid seed set by user |
208 |
|
|
|
209 |
|
|
#ifndef IS_MPI |
210 |
|
|
|
211 |
|
|
seedValue = make_sprng_seed(); |
212 |
|
|
|
213 |
|
|
#else |
214 |
|
|
|
215 |
|
|
if (worldRank == 0) { |
216 |
|
|
seedValue = make_sprng_seed(); |
217 |
|
|
} |
218 |
|
|
|
219 |
|
|
MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD); |
220 |
|
|
|
221 |
|
|
#endif |
222 |
|
|
|
223 |
|
|
} //end if (seedValue /1000000000 == 0) |
224 |
|
|
} else { |
225 |
|
|
|
226 |
|
|
#ifndef IS_MPI |
227 |
|
|
|
228 |
|
|
seedValue = make_sprng_seed(); |
229 |
|
|
|
230 |
|
|
#else |
231 |
|
|
|
232 |
|
|
if (worldRank == 0) { |
233 |
|
|
seedValue = make_sprng_seed(); |
234 |
|
|
} |
235 |
|
|
|
236 |
|
|
MPI_Bcast(&seedValue, 1, MPI_INT, 0, MPI_COMM_WORLD); |
237 |
|
|
|
238 |
|
|
#endif |
239 |
|
|
|
240 |
|
|
} //end of simParams->haveSeed() |
241 |
|
|
|
242 |
|
|
info->setSeed(seedValue); |
243 |
|
|
|
244 |
|
|
|
245 |
|
|
//figure out the ouput file names |
246 |
|
|
std::string prefix; |
247 |
|
|
|
248 |
|
|
#ifdef IS_MPI |
249 |
|
|
|
250 |
|
|
if (worldRank == 0) { |
251 |
|
|
#endif // is_mpi |
252 |
|
|
|
253 |
|
|
if (simParams->haveFinalConfig()) { |
254 |
|
|
prefix = getPrefix(simParams->getFinalConfig()); |
255 |
|
|
} else { |
256 |
|
|
prefix = getPrefix(mdfile); |
257 |
|
|
} |
258 |
|
|
|
259 |
|
|
info->setFinalConfigFileName(prefix + ".eor"); |
260 |
|
|
info->setDumpFileName(prefix + ".dump"); |
261 |
|
|
info->setStatFileName(prefix + ".stat"); |
262 |
|
|
|
263 |
|
|
#ifdef IS_MPI |
264 |
|
|
|
265 |
|
|
} |
266 |
|
|
|
267 |
|
|
#endif |
268 |
|
|
|
269 |
|
|
} |
270 |
|
|
|
271 |
|
|
#ifdef IS_MPI |
272 |
|
|
void SimCreator::divideMolecules(SimInfo *info) { |
273 |
|
|
double numerator; |
274 |
|
|
double denominator; |
275 |
|
|
double precast; |
276 |
|
|
double x; |
277 |
|
|
double y; |
278 |
|
|
double a; |
279 |
|
|
int old_atoms; |
280 |
|
|
int add_atoms; |
281 |
|
|
int new_atoms; |
282 |
|
|
int nTarget; |
283 |
|
|
int done; |
284 |
|
|
int i; |
285 |
|
|
int j; |
286 |
|
|
int loops; |
287 |
|
|
int which_proc; |
288 |
|
|
int nProcessors; |
289 |
|
|
std::vector<int> atomsPerProc; |
290 |
|
|
randomSPRNG myRandom(info->getSeed()); |
291 |
|
|
int nGlobalMols = info->getNGlobalMolecules(); |
292 |
|
|
std::vector<int> molToProcMap(nGlobalMols, -1); // default to an error condition: |
293 |
|
|
|
294 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &nProcessors); |
295 |
|
|
|
296 |
|
|
if (nProcessors > nGlobalMols) { |
297 |
|
|
sprintf(painCave.errMsg, |
298 |
|
|
"nProcessors (%d) > nMol (%d)\n" |
299 |
|
|
"\tThe number of processors is larger than\n" |
300 |
|
|
"\tthe number of molecules. This will not result in a \n" |
301 |
|
|
"\tusable division of atoms for force decomposition.\n" |
302 |
|
|
"\tEither try a smaller number of processors, or run the\n" |
303 |
|
|
"\tsingle-processor version of OOPSE.\n", nProcessors, nGlobalMols); |
304 |
|
|
|
305 |
|
|
painCave.isFatal = 1; |
306 |
|
|
simError(); |
307 |
|
|
} |
308 |
|
|
|
309 |
|
|
a = 3.0 * nGlobalMols / info->getNGlobalAtoms(); |
310 |
|
|
|
311 |
|
|
//initialize atomsPerProc |
312 |
|
|
atomsPerProc.insert(atomsPerProc.end(), nProcessors, 0); |
313 |
|
|
|
314 |
|
|
if (worldRank == 0) { |
315 |
|
|
numerator = info->getNGlobalAtoms(); |
316 |
|
|
denominator = nProcessors; |
317 |
|
|
precast = numerator / denominator; |
318 |
|
|
nTarget = (int)(precast + 0.5); |
319 |
|
|
|
320 |
|
|
for(i = 0; i < nGlobalMols; i++) { |
321 |
|
|
done = 0; |
322 |
|
|
loops = 0; |
323 |
|
|
|
324 |
|
|
while (!done) { |
325 |
|
|
loops++; |
326 |
|
|
|
327 |
|
|
// Pick a processor at random |
328 |
|
|
|
329 |
|
|
which_proc = (int) (myRandom.getRandom() * nProcessors); |
330 |
|
|
|
331 |
|
|
//get the molecule stamp first |
332 |
|
|
int stampId = info->getMoleculeStampId(i); |
333 |
|
|
MoleculeStamp * moleculeStamp = info->getMoleculeStamp(stampId); |
334 |
|
|
|
335 |
|
|
// How many atoms does this processor have so far? |
336 |
|
|
old_atoms = atomsPerProc[which_proc]; |
337 |
|
|
add_atoms = moleculeStamp->getNAtoms(); |
338 |
|
|
new_atoms = old_atoms + add_atoms; |
339 |
|
|
|
340 |
|
|
// If we've been through this loop too many times, we need |
341 |
|
|
// to just give up and assign the molecule to this processor |
342 |
|
|
// and be done with it. |
343 |
|
|
|
344 |
|
|
if (loops > 100) { |
345 |
|
|
sprintf(painCave.errMsg, |
346 |
|
|
"I've tried 100 times to assign molecule %d to a " |
347 |
|
|
" processor, but can't find a good spot.\n" |
348 |
|
|
"I'm assigning it at random to processor %d.\n", |
349 |
|
|
i, which_proc); |
350 |
|
|
|
351 |
|
|
painCave.isFatal = 0; |
352 |
|
|
simError(); |
353 |
|
|
|
354 |
|
|
molToProcMap[i] = which_proc; |
355 |
|
|
atomsPerProc[which_proc] += add_atoms; |
356 |
|
|
|
357 |
|
|
done = 1; |
358 |
|
|
continue; |
359 |
|
|
} |
360 |
|
|
|
361 |
|
|
// If we can add this molecule to this processor without sending |
362 |
|
|
// it above nTarget, then go ahead and do it: |
363 |
|
|
|
364 |
|
|
if (new_atoms <= nTarget) { |
365 |
|
|
molToProcMap[i] = which_proc; |
366 |
|
|
atomsPerProc[which_proc] += add_atoms; |
367 |
|
|
|
368 |
|
|
done = 1; |
369 |
|
|
continue; |
370 |
|
|
} |
371 |
|
|
|
372 |
|
|
// The only situation left is when new_atoms > nTarget. We |
373 |
|
|
// want to accept this with some probability that dies off the |
374 |
|
|
// farther we are from nTarget |
375 |
|
|
|
376 |
|
|
// roughly: x = new_atoms - nTarget |
377 |
|
|
// Pacc(x) = exp(- a * x) |
378 |
|
|
// where a = penalty / (average atoms per molecule) |
379 |
|
|
|
380 |
|
|
x = (double)(new_atoms - nTarget); |
381 |
|
|
y = myRandom.getRandom(); |
382 |
|
|
|
383 |
|
|
if (y < exp(- a * x)) { |
384 |
|
|
molToProcMap[i] = which_proc; |
385 |
|
|
atomsPerProc[which_proc] += add_atoms; |
386 |
|
|
|
387 |
|
|
done = 1; |
388 |
|
|
continue; |
389 |
|
|
} else { |
390 |
|
|
continue; |
391 |
|
|
} |
392 |
|
|
} |
393 |
|
|
} |
394 |
|
|
|
395 |
|
|
// Spray out this nonsense to all other processors: |
396 |
|
|
|
397 |
|
|
MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); |
398 |
|
|
} else { |
399 |
|
|
|
400 |
|
|
// Listen to your marching orders from processor 0: |
401 |
|
|
|
402 |
|
|
MPI_Bcast(&molToProcMap[0], nGlobalMols, MPI_INT, 0, MPI_COMM_WORLD); |
403 |
|
|
} |
404 |
|
|
|
405 |
|
|
info->setMolToProcMap(molToProcMap); |
406 |
|
|
sprintf(checkPointMsg, |
407 |
|
|
"Successfully divided the molecules among the processors.\n"); |
408 |
|
|
MPIcheckPoint(); |
409 |
|
|
} |
410 |
|
|
|
411 |
|
|
#endif |
412 |
|
|
|
413 |
|
|
void SimCreator::createMolecules(SimInfo *info) { |
414 |
|
|
MoleculeCreator molCreator; |
415 |
|
|
int stampId; |
416 |
|
|
|
417 |
|
|
for(int i = 0; i < info->getNGlobalMolecules(); i++) { |
418 |
|
|
|
419 |
|
|
#ifdef IS_MPI |
420 |
|
|
|
421 |
|
|
if (info->getMolToProc(i) == worldRank) { |
422 |
|
|
#endif |
423 |
|
|
|
424 |
|
|
stampId = info->getMoleculeStampId(i); |
425 |
|
|
Molecule * mol = molCreator.createMolecule(info->getForceField(), info->getMoleculeStamp(stampId), |
426 |
|
|
stampId, i, info->getLocalIndexManager()); |
427 |
|
|
|
428 |
|
|
info->addMolecule(mol); |
429 |
|
|
|
430 |
|
|
#ifdef IS_MPI |
431 |
|
|
|
432 |
|
|
} |
433 |
|
|
|
434 |
|
|
#endif |
435 |
|
|
|
436 |
|
|
} //end for(int i=0) |
437 |
|
|
} |
438 |
|
|
|
439 |
|
|
void SimCreator::compList(MakeStamps *stamps, Globals* simParams, |
440 |
|
|
std::vector < std::pair<MoleculeStamp *, int> > &moleculeStampPairs) { |
441 |
|
|
int i; |
442 |
|
|
char * id; |
443 |
|
|
MoleculeStamp * currentStamp; |
444 |
|
|
Component** the_components = simParams->getComponents(); |
445 |
|
|
int n_components = simParams->getNComponents(); |
446 |
|
|
|
447 |
|
|
if (!simParams->haveNMol()) { |
448 |
|
|
// we don't have the total number of molecules, so we assume it is |
449 |
|
|
// given in each component |
450 |
|
|
|
451 |
|
|
for(i = 0; i < n_components; i++) { |
452 |
|
|
if (!the_components[i]->haveNMol()) { |
453 |
|
|
// we have a problem |
454 |
|
|
sprintf(painCave.errMsg, |
455 |
|
|
"SimCreator Error. No global NMol or component NMol given.\n" |
456 |
|
|
"\tCannot calculate the number of atoms.\n"); |
457 |
|
|
|
458 |
|
|
painCave.isFatal = 1; |
459 |
|
|
simError(); |
460 |
|
|
} |
461 |
|
|
|
462 |
|
|
id = the_components[i]->getType(); |
463 |
|
|
currentStamp = (stamps->extractMolStamp(id))->getStamp(); |
464 |
|
|
|
465 |
|
|
if (currentStamp == NULL) { |
466 |
|
|
sprintf(painCave.errMsg, |
467 |
|
|
"SimCreator error: Component \"%s\" was not found in the " |
468 |
|
|
"list of declared molecules\n", id); |
469 |
|
|
|
470 |
|
|
painCave.isFatal = 1; |
471 |
|
|
simError(); |
472 |
|
|
} |
473 |
|
|
|
474 |
|
|
moleculeStampPairs.push_back( |
475 |
|
|
std::make_pair(currentStamp, the_components[i]->getNMol())); |
476 |
|
|
} //end for (i = 0; i < n_components; i++) |
477 |
|
|
} else { |
478 |
|
|
sprintf(painCave.errMsg, "SimSetup error.\n" |
479 |
|
|
"\tSorry, the ability to specify total" |
480 |
|
|
" nMols and then give molfractions in the components\n" |
481 |
|
|
"\tis not currently supported." |
482 |
|
|
" Please give nMol in the components.\n"); |
483 |
|
|
|
484 |
|
|
painCave.isFatal = 1; |
485 |
|
|
simError(); |
486 |
|
|
} |
487 |
|
|
|
488 |
|
|
#ifdef IS_MPI |
489 |
|
|
|
490 |
|
|
strcpy(checkPointMsg, "Component stamps successfully extracted\n"); |
491 |
|
|
MPIcheckPoint(); |
492 |
|
|
|
493 |
|
|
#endif // is_mpi |
494 |
|
|
|
495 |
|
|
} |
496 |
|
|
|
497 |
|
|
void SimCreator::setGlobalIndex(SimInfo *info) { |
498 |
|
|
SimInfo::MoleculeIterator mi; |
499 |
|
|
Molecule::AtomIterator ai; |
500 |
|
|
Molecule::RigidBodyIterator ri; |
501 |
|
|
Molecule::CutoffGroupIterator ci; |
502 |
|
|
Molecule * mol; |
503 |
|
|
Atom * atom; |
504 |
|
|
RigidBody * rb; |
505 |
|
|
CutoffGroup * cg; |
506 |
|
|
int beginAtomIndex; |
507 |
|
|
int beginRigidBodyIndex; |
508 |
|
|
int beginCutoffGroupIndex; |
509 |
|
|
int nGlobalAtoms = info->getNGlobalAtoms(); |
510 |
|
|
|
511 |
|
|
#ifndef IS_MPI |
512 |
|
|
|
513 |
|
|
beginAtomIndex = 0; |
514 |
|
|
beginRigidBodyIndex = 0; |
515 |
|
|
beginCutoffGroupIndex = 0; |
516 |
|
|
|
517 |
|
|
#else |
518 |
|
|
|
519 |
|
|
int nproc; |
520 |
|
|
int myNode; |
521 |
|
|
|
522 |
|
|
myNode = worldRank; |
523 |
|
|
MPI_Comm_size(MPI_COMM_WORLD, &nproc); |
524 |
|
|
|
525 |
|
|
std::vector < int > tmpAtomsInProc(nproc, 0); |
526 |
|
|
std::vector < int > tmpRigidBodiesInProc(nproc, 0); |
527 |
|
|
std::vector < int > tmpCutoffGroupsInProc(nproc, 0); |
528 |
|
|
std::vector < int > NumAtomsInProc(nproc, 0); |
529 |
|
|
std::vector < int > NumRigidBodiesInProc(nproc, 0); |
530 |
|
|
std::vector < int > NumCutoffGroupsInProc(nproc, 0); |
531 |
|
|
|
532 |
|
|
tmpAtomsInProc[myNode] = info->getNAtoms(); |
533 |
|
|
tmpRigidBodiesInProc[myNode] = info->getNRigidBodies(); |
534 |
|
|
tmpCutoffGroupsInProc[myNode] = info->getNCutoffGroups(); |
535 |
|
|
|
536 |
|
|
//do MPI_ALLREDUCE to exchange the total number of atoms, rigidbodies and cutoff groups |
537 |
|
|
MPI_Allreduce(&tmpAtomsInProc[0], &NumAtomsInProc[0], nproc, MPI_INT, |
538 |
|
|
MPI_SUM, MPI_COMM_WORLD); |
539 |
|
|
MPI_Allreduce(&tmpRigidBodiesInProc[0], &NumRigidBodiesInProc[0], nproc, |
540 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
541 |
|
|
MPI_Allreduce(&tmpCutoffGroupsInProc[0], &NumCutoffGroupsInProc[0], nproc, |
542 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
543 |
|
|
|
544 |
|
|
beginAtomIndex = 0; |
545 |
|
|
beginRigidBodyIndex = 0; |
546 |
|
|
beginCutoffGroupIndex = 0; |
547 |
|
|
|
548 |
|
|
for(int i = 0; i < myNode; i++) { |
549 |
|
|
beginAtomIndex += NumAtomsInProc[i]; |
550 |
|
|
beginRigidBodyIndex += NumRigidBodiesInProc[i]; |
551 |
|
|
beginCutoffGroupIndex += NumCutoffGroupsInProc[i]; |
552 |
|
|
} |
553 |
|
|
|
554 |
tim |
285 |
//rigidbody's index begins right after atom's |
555 |
|
|
beginRigidBodyIndex += info->getNGlobalAtoms(); |
556 |
gezelter |
246 |
#endif |
557 |
|
|
|
558 |
|
|
for(mol = info->beginMolecule(mi); mol != NULL; |
559 |
|
|
mol = info->nextMolecule(mi)) { |
560 |
|
|
|
561 |
|
|
//local index(index in DataStorge) of atom is important |
562 |
|
|
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
563 |
|
|
atom->setGlobalIndex(beginAtomIndex++); |
564 |
|
|
} |
565 |
|
|
|
566 |
|
|
for(rb = mol->beginRigidBody(ri); rb != NULL; |
567 |
|
|
rb = mol->nextRigidBody(ri)) { |
568 |
|
|
rb->setGlobalIndex(beginRigidBodyIndex++); |
569 |
|
|
} |
570 |
|
|
|
571 |
|
|
//local index of cutoff group is trivial, it only depends on the order of travesing |
572 |
|
|
for(cg = mol->beginCutoffGroup(ci); cg != NULL; |
573 |
|
|
cg = mol->nextCutoffGroup(ci)) { |
574 |
|
|
cg->setGlobalIndex(beginCutoffGroupIndex++); |
575 |
|
|
} |
576 |
|
|
} |
577 |
|
|
|
578 |
|
|
//fill globalGroupMembership |
579 |
|
|
std::vector<int> globalGroupMembership(info->getNGlobalAtoms(), 0); |
580 |
|
|
for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
581 |
|
|
for (cg = mol->beginCutoffGroup(ci); cg != NULL; cg = mol->nextCutoffGroup(ci)) { |
582 |
|
|
|
583 |
|
|
for(atom = cg->beginAtom(ai); atom != NULL; atom = cg->nextAtom(ai)) { |
584 |
|
|
globalGroupMembership[atom->getGlobalIndex()] = cg->getGlobalIndex(); |
585 |
|
|
} |
586 |
|
|
|
587 |
|
|
} |
588 |
|
|
} |
589 |
|
|
|
590 |
|
|
#ifdef IS_MPI |
591 |
|
|
// Since the globalGroupMembership has been zero filled and we've only |
592 |
|
|
// poked values into the atoms we know, we can do an Allreduce |
593 |
|
|
// to get the full globalGroupMembership array (We think). |
594 |
|
|
// This would be prettier if we could use MPI_IN_PLACE like the MPI-2 |
595 |
|
|
// docs said we could. |
596 |
|
|
std::vector<int> tmpGroupMembership(nGlobalAtoms, 0); |
597 |
|
|
MPI_Allreduce(&globalGroupMembership[0], &tmpGroupMembership[0], nGlobalAtoms, |
598 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
599 |
|
|
info->setGlobalGroupMembership(tmpGroupMembership); |
600 |
|
|
#else |
601 |
|
|
info->setGlobalGroupMembership(globalGroupMembership); |
602 |
|
|
#endif |
603 |
|
|
|
604 |
|
|
//fill molMembership |
605 |
|
|
std::vector<int> globalMolMembership(info->getNGlobalAtoms(), 0); |
606 |
|
|
|
607 |
|
|
for(mol = info->beginMolecule(mi); mol != NULL; mol = info->nextMolecule(mi)) { |
608 |
|
|
|
609 |
|
|
for(atom = mol->beginAtom(ai); atom != NULL; atom = mol->nextAtom(ai)) { |
610 |
|
|
globalMolMembership[atom->getGlobalIndex()] = mol->getGlobalIndex(); |
611 |
|
|
} |
612 |
|
|
} |
613 |
|
|
|
614 |
|
|
#ifdef IS_MPI |
615 |
|
|
std::vector<int> tmpMolMembership(nGlobalAtoms, 0); |
616 |
|
|
|
617 |
|
|
MPI_Allreduce(&globalMolMembership[0], &tmpMolMembership[0], nGlobalAtoms, |
618 |
|
|
MPI_INT, MPI_SUM, MPI_COMM_WORLD); |
619 |
|
|
|
620 |
|
|
info->setGlobalMolMembership(tmpMolMembership); |
621 |
|
|
#else |
622 |
|
|
info->setGlobalMolMembership(globalMolMembership); |
623 |
|
|
#endif |
624 |
|
|
|
625 |
|
|
} |
626 |
|
|
|
627 |
|
|
void SimCreator::loadCoordinates(SimInfo* info) { |
628 |
|
|
Globals* simParams; |
629 |
|
|
simParams = info->getSimParams(); |
630 |
|
|
|
631 |
|
|
if (!simParams->haveInitialConfig()) { |
632 |
|
|
sprintf(painCave.errMsg, |
633 |
|
|
"Cannot intialize a simulation without an initial configuration file.\n"); |
634 |
|
|
painCave.isFatal = 1;; |
635 |
|
|
simError(); |
636 |
|
|
} |
637 |
|
|
|
638 |
|
|
DumpReader reader(info, simParams->getInitialConfig()); |
639 |
|
|
int nframes = reader.getNFrames(); |
640 |
|
|
|
641 |
|
|
if (nframes > 0) { |
642 |
|
|
reader.readFrame(nframes - 1); |
643 |
|
|
} else { |
644 |
|
|
//invalid initial coordinate file |
645 |
|
|
sprintf(painCave.errMsg, "Initial configuration file %s should at least contain one frame\n", |
646 |
|
|
simParams->getInitialConfig()); |
647 |
|
|
painCave.isFatal = 1; |
648 |
|
|
simError(); |
649 |
|
|
} |
650 |
|
|
|
651 |
|
|
//copy the current snapshot to previous snapshot |
652 |
|
|
info->getSnapshotManager()->advance(); |
653 |
|
|
} |
654 |
|
|
|
655 |
|
|
} //end namespace oopse |
656 |
|
|
|
657 |
|
|
|