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 |
|
|
#include "io/DumpWriter.hpp" |
43 |
|
|
#include "primitives/Molecule.hpp" |
44 |
|
|
#include "utils/simError.h" |
45 |
tim |
376 |
#include "io/basic_teebuf.hpp" |
46 |
gezelter |
2 |
#ifdef IS_MPI |
47 |
|
|
#include <mpi.h> |
48 |
gezelter |
246 |
#endif //is_mpi |
49 |
gezelter |
2 |
|
50 |
gezelter |
246 |
namespace oopse { |
51 |
gezelter |
2 |
|
52 |
tim |
376 |
DumpWriter::DumpWriter(SimInfo* info) |
53 |
|
|
: info_(info), filename_(info->getDumpFileName()), eorFilename_(info->getFinalConfigFileName()){ |
54 |
|
|
#ifdef IS_MPI |
55 |
|
|
|
56 |
|
|
if (worldRank == 0) { |
57 |
|
|
#endif // is_mpi |
58 |
|
|
|
59 |
|
|
dumpFile_.open(filename_.c_str(), std::ios::out | std::ios::trunc); |
60 |
|
|
|
61 |
|
|
if (!dumpFile_) { |
62 |
|
|
sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n", |
63 |
|
|
filename_.c_str()); |
64 |
|
|
painCave.isFatal = 1; |
65 |
|
|
simError(); |
66 |
|
|
} |
67 |
|
|
|
68 |
|
|
#ifdef IS_MPI |
69 |
|
|
|
70 |
|
|
} |
71 |
|
|
|
72 |
|
|
sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n"); |
73 |
|
|
MPIcheckPoint(); |
74 |
|
|
|
75 |
|
|
#endif // is_mpi |
76 |
|
|
|
77 |
|
|
} |
78 |
|
|
|
79 |
|
|
|
80 |
gezelter |
246 |
DumpWriter::DumpWriter(SimInfo* info, const std::string& filename) |
81 |
|
|
: info_(info), filename_(filename){ |
82 |
|
|
#ifdef IS_MPI |
83 |
gezelter |
2 |
|
84 |
gezelter |
246 |
if (worldRank == 0) { |
85 |
|
|
#endif // is_mpi |
86 |
gezelter |
2 |
|
87 |
tim |
376 |
eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor"; |
88 |
gezelter |
246 |
dumpFile_.open(filename_.c_str(), std::ios::out | std::ios::trunc); |
89 |
gezelter |
2 |
|
90 |
gezelter |
246 |
if (!dumpFile_) { |
91 |
|
|
sprintf(painCave.errMsg, "Could not open \"%s\" for dump output.\n", |
92 |
|
|
filename_.c_str()); |
93 |
|
|
painCave.isFatal = 1; |
94 |
|
|
simError(); |
95 |
|
|
} |
96 |
gezelter |
2 |
|
97 |
|
|
#ifdef IS_MPI |
98 |
|
|
|
99 |
|
|
} |
100 |
|
|
|
101 |
gezelter |
246 |
sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n"); |
102 |
|
|
MPIcheckPoint(); |
103 |
gezelter |
2 |
|
104 |
|
|
#endif // is_mpi |
105 |
gezelter |
246 |
|
106 |
gezelter |
2 |
} |
107 |
|
|
|
108 |
gezelter |
246 |
DumpWriter::~DumpWriter() { |
109 |
gezelter |
2 |
|
110 |
|
|
#ifdef IS_MPI |
111 |
gezelter |
246 |
|
112 |
|
|
if (worldRank == 0) { |
113 |
gezelter |
2 |
#endif // is_mpi |
114 |
|
|
|
115 |
gezelter |
246 |
dumpFile_.close(); |
116 |
gezelter |
2 |
|
117 |
|
|
#ifdef IS_MPI |
118 |
gezelter |
246 |
|
119 |
|
|
} |
120 |
|
|
|
121 |
gezelter |
2 |
#endif // is_mpi |
122 |
gezelter |
246 |
|
123 |
gezelter |
2 |
} |
124 |
|
|
|
125 |
gezelter |
246 |
void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) { |
126 |
gezelter |
2 |
|
127 |
gezelter |
246 |
double currentTime; |
128 |
|
|
Mat3x3d hmat; |
129 |
|
|
double chi; |
130 |
|
|
double integralOfChiDt; |
131 |
|
|
Mat3x3d eta; |
132 |
|
|
|
133 |
|
|
currentTime = s->getTime(); |
134 |
|
|
hmat = s->getHmat(); |
135 |
|
|
chi = s->getChi(); |
136 |
|
|
integralOfChiDt = s->getIntegralOfChiDt(); |
137 |
|
|
eta = s->getEta(); |
138 |
|
|
|
139 |
|
|
os << currentTime << ";\t" |
140 |
|
|
<< hmat(0, 0) << "\t" << hmat(1, 0) << "\t" << hmat(2, 0) << ";\t" |
141 |
|
|
<< hmat(0, 1) << "\t" << hmat(1, 1) << "\t" << hmat(2, 1) << ";\t" |
142 |
|
|
<< hmat(0, 2) << "\t" << hmat(1, 2) << "\t" << hmat(2, 2) << ";\t"; |
143 |
gezelter |
2 |
|
144 |
gezelter |
246 |
//write out additional parameters, such as chi and eta |
145 |
|
|
|
146 |
|
|
os << chi << "\t" << integralOfChiDt << "\t;"; |
147 |
|
|
|
148 |
|
|
os << eta(0, 0) << "\t" << eta(1, 0) << "\t" << eta(2, 0) << ";\t" |
149 |
|
|
<< eta(0, 1) << "\t" << eta(1, 1) << "\t" << eta(2, 1) << ";\t" |
150 |
|
|
<< eta(0, 2) << "\t" << eta(1, 2) << "\t" << eta(2, 2) << ";"; |
151 |
|
|
|
152 |
tim |
376 |
os << "\n"; |
153 |
gezelter |
2 |
} |
154 |
|
|
|
155 |
gezelter |
246 |
void DumpWriter::writeFrame(std::ostream& os) { |
156 |
|
|
const int BUFFERSIZE = 2000; |
157 |
|
|
const int MINIBUFFERSIZE = 100; |
158 |
|
|
|
159 |
|
|
char tempBuffer[BUFFERSIZE]; |
160 |
|
|
char writeLine[BUFFERSIZE]; |
161 |
|
|
|
162 |
|
|
Quat4d q; |
163 |
|
|
Vector3d ji; |
164 |
|
|
Vector3d pos; |
165 |
|
|
Vector3d vel; |
166 |
|
|
|
167 |
|
|
Molecule* mol; |
168 |
|
|
StuntDouble* integrableObject; |
169 |
|
|
SimInfo::MoleculeIterator mi; |
170 |
|
|
Molecule::IntegrableObjectIterator ii; |
171 |
gezelter |
2 |
|
172 |
gezelter |
246 |
int nTotObjects; |
173 |
|
|
nTotObjects = info_->getNGlobalIntegrableObjects(); |
174 |
gezelter |
2 |
|
175 |
gezelter |
246 |
#ifndef IS_MPI |
176 |
gezelter |
2 |
|
177 |
|
|
|
178 |
gezelter |
246 |
os << nTotObjects << "\n"; |
179 |
|
|
|
180 |
|
|
writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot()); |
181 |
gezelter |
2 |
|
182 |
gezelter |
246 |
for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) { |
183 |
gezelter |
2 |
|
184 |
gezelter |
246 |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
185 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
186 |
|
|
|
187 |
gezelter |
2 |
|
188 |
gezelter |
246 |
pos = integrableObject->getPos(); |
189 |
|
|
vel = integrableObject->getVel(); |
190 |
gezelter |
2 |
|
191 |
gezelter |
246 |
sprintf(tempBuffer, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
192 |
|
|
integrableObject->getType().c_str(), |
193 |
|
|
pos[0], pos[1], pos[2], |
194 |
|
|
vel[0], vel[1], vel[2]); |
195 |
gezelter |
2 |
|
196 |
gezelter |
246 |
strcpy(writeLine, tempBuffer); |
197 |
gezelter |
2 |
|
198 |
gezelter |
246 |
if (integrableObject->isDirectional()) { |
199 |
|
|
q = integrableObject->getQ(); |
200 |
|
|
ji = integrableObject->getJ(); |
201 |
gezelter |
2 |
|
202 |
gezelter |
246 |
sprintf(tempBuffer, "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
203 |
|
|
q[0], q[1], q[2], q[3], |
204 |
|
|
ji[0], ji[1], ji[2]); |
205 |
|
|
strcat(writeLine, tempBuffer); |
206 |
|
|
} else { |
207 |
|
|
strcat(writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
208 |
|
|
} |
209 |
gezelter |
2 |
|
210 |
gezelter |
246 |
os << writeLine; |
211 |
gezelter |
2 |
|
212 |
gezelter |
246 |
} |
213 |
gezelter |
2 |
} |
214 |
|
|
|
215 |
tim |
324 |
os.flush(); |
216 |
gezelter |
246 |
#else // is_mpi |
217 |
|
|
/********************************************************************* |
218 |
|
|
* Documentation? You want DOCUMENTATION? |
219 |
|
|
* |
220 |
|
|
* Why all the potatoes below? |
221 |
|
|
* |
222 |
|
|
* To make a long story short, the original version of DumpWriter |
223 |
|
|
* worked in the most inefficient way possible. Node 0 would |
224 |
|
|
* poke each of the node for an individual atom's formatted data |
225 |
|
|
* as node 0 worked its way down the global index. This was particularly |
226 |
|
|
* inefficient since the method blocked all processors at every atom |
227 |
|
|
* (and did it twice!). |
228 |
|
|
* |
229 |
|
|
* An intermediate version of DumpWriter could be described from Node |
230 |
|
|
* zero's perspective as follows: |
231 |
|
|
* |
232 |
|
|
* 1) Have 100 of your friends stand in a circle. |
233 |
|
|
* 2) When you say go, have all of them start tossing potatoes at |
234 |
|
|
* you (one at a time). |
235 |
|
|
* 3) Catch the potatoes. |
236 |
|
|
* |
237 |
|
|
* It was an improvement, but MPI has buffers and caches that could |
238 |
|
|
* best be described in this analogy as "potato nets", so there's no |
239 |
|
|
* need to block the processors atom-by-atom. |
240 |
|
|
* |
241 |
|
|
* This new and improved DumpWriter works in an even more efficient |
242 |
|
|
* way: |
243 |
|
|
* |
244 |
|
|
* 1) Have 100 of your friend stand in a circle. |
245 |
|
|
* 2) When you say go, have them start tossing 5-pound bags of |
246 |
|
|
* potatoes at you. |
247 |
|
|
* 3) Once you've caught a friend's bag of potatoes, |
248 |
|
|
* toss them a spud to let them know they can toss another bag. |
249 |
|
|
* |
250 |
|
|
* How's THAT for documentation? |
251 |
|
|
* |
252 |
|
|
*********************************************************************/ |
253 |
|
|
const int masterNode = 0; |
254 |
gezelter |
2 |
|
255 |
gezelter |
246 |
int * potatoes; |
256 |
|
|
int myPotato; |
257 |
|
|
int nProc; |
258 |
|
|
int which_node; |
259 |
|
|
double atomData[13]; |
260 |
|
|
int isDirectional; |
261 |
|
|
const char * atomTypeString; |
262 |
|
|
char MPIatomTypeString[MINIBUFFERSIZE]; |
263 |
|
|
int msgLen; // the length of message actually recieved at master nodes |
264 |
|
|
int haveError; |
265 |
|
|
MPI_Status istatus; |
266 |
|
|
int nCurObj; |
267 |
|
|
|
268 |
|
|
// code to find maximum tag value |
269 |
|
|
int * tagub; |
270 |
|
|
int flag; |
271 |
|
|
int MAXTAG; |
272 |
|
|
MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag); |
273 |
gezelter |
2 |
|
274 |
gezelter |
246 |
if (flag) { |
275 |
|
|
MAXTAG = *tagub; |
276 |
|
|
} else { |
277 |
|
|
MAXTAG = 32767; |
278 |
|
|
} |
279 |
gezelter |
2 |
|
280 |
gezelter |
246 |
if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file |
281 |
gezelter |
2 |
|
282 |
gezelter |
246 |
// Node 0 needs a list of the magic potatoes for each processor; |
283 |
gezelter |
2 |
|
284 |
gezelter |
246 |
MPI_Comm_size(MPI_COMM_WORLD, &nProc); |
285 |
|
|
potatoes = new int[nProc]; |
286 |
gezelter |
2 |
|
287 |
gezelter |
246 |
//write out the comment lines |
288 |
|
|
for(int i = 0; i < nProc; i++) { |
289 |
|
|
potatoes[i] = 0; |
290 |
|
|
} |
291 |
gezelter |
2 |
|
292 |
|
|
|
293 |
gezelter |
246 |
os << nTotObjects << "\n"; |
294 |
|
|
writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot()); |
295 |
gezelter |
2 |
|
296 |
gezelter |
246 |
for(int i = 0; i < info_->getNGlobalMolecules(); i++) { |
297 |
gezelter |
2 |
|
298 |
gezelter |
246 |
// Get the Node number which has this atom; |
299 |
gezelter |
2 |
|
300 |
gezelter |
246 |
which_node = info_->getMolToProc(i); |
301 |
gezelter |
2 |
|
302 |
gezelter |
246 |
if (which_node != masterNode) { //current molecule is in slave node |
303 |
|
|
if (potatoes[which_node] + 1 >= MAXTAG) { |
304 |
|
|
// The potato was going to exceed the maximum value, |
305 |
|
|
// so wrap this processor potato back to 0: |
306 |
gezelter |
2 |
|
307 |
gezelter |
246 |
potatoes[which_node] = 0; |
308 |
|
|
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0, |
309 |
|
|
MPI_COMM_WORLD); |
310 |
|
|
} |
311 |
gezelter |
2 |
|
312 |
gezelter |
246 |
myPotato = potatoes[which_node]; |
313 |
gezelter |
2 |
|
314 |
gezelter |
246 |
//recieve the number of integrableObject in current molecule |
315 |
|
|
MPI_Recv(&nCurObj, 1, MPI_INT, which_node, myPotato, |
316 |
|
|
MPI_COMM_WORLD, &istatus); |
317 |
|
|
myPotato++; |
318 |
gezelter |
2 |
|
319 |
gezelter |
246 |
for(int l = 0; l < nCurObj; l++) { |
320 |
|
|
if (potatoes[which_node] + 2 >= MAXTAG) { |
321 |
|
|
// The potato was going to exceed the maximum value, |
322 |
|
|
// so wrap this processor potato back to 0: |
323 |
gezelter |
2 |
|
324 |
gezelter |
246 |
potatoes[which_node] = 0; |
325 |
|
|
MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, |
326 |
|
|
0, MPI_COMM_WORLD); |
327 |
|
|
} |
328 |
gezelter |
2 |
|
329 |
gezelter |
246 |
MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, |
330 |
|
|
which_node, myPotato, MPI_COMM_WORLD, |
331 |
|
|
&istatus); |
332 |
gezelter |
2 |
|
333 |
gezelter |
246 |
atomTypeString = MPIatomTypeString; |
334 |
gezelter |
2 |
|
335 |
gezelter |
246 |
myPotato++; |
336 |
gezelter |
2 |
|
337 |
gezelter |
246 |
MPI_Recv(atomData, 13, MPI_DOUBLE, which_node, myPotato, |
338 |
|
|
MPI_COMM_WORLD, &istatus); |
339 |
|
|
myPotato++; |
340 |
gezelter |
2 |
|
341 |
gezelter |
246 |
MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen); |
342 |
gezelter |
2 |
|
343 |
gezelter |
246 |
if (msgLen == 13) |
344 |
|
|
isDirectional = 1; |
345 |
|
|
else |
346 |
|
|
isDirectional = 0; |
347 |
gezelter |
2 |
|
348 |
gezelter |
246 |
// If we've survived to here, format the line: |
349 |
gezelter |
2 |
|
350 |
gezelter |
246 |
if (!isDirectional) { |
351 |
|
|
sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
352 |
|
|
atomTypeString, atomData[0], |
353 |
|
|
atomData[1], atomData[2], |
354 |
|
|
atomData[3], atomData[4], |
355 |
|
|
atomData[5]); |
356 |
gezelter |
2 |
|
357 |
gezelter |
246 |
strcat(writeLine, |
358 |
|
|
"0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
359 |
|
|
} else { |
360 |
|
|
sprintf(writeLine, |
361 |
|
|
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
362 |
|
|
atomTypeString, |
363 |
|
|
atomData[0], |
364 |
|
|
atomData[1], |
365 |
|
|
atomData[2], |
366 |
|
|
atomData[3], |
367 |
|
|
atomData[4], |
368 |
|
|
atomData[5], |
369 |
|
|
atomData[6], |
370 |
|
|
atomData[7], |
371 |
|
|
atomData[8], |
372 |
|
|
atomData[9], |
373 |
|
|
atomData[10], |
374 |
|
|
atomData[11], |
375 |
|
|
atomData[12]); |
376 |
|
|
} |
377 |
gezelter |
2 |
|
378 |
gezelter |
246 |
os << writeLine; |
379 |
gezelter |
2 |
|
380 |
gezelter |
246 |
} // end for(int l =0) |
381 |
gezelter |
2 |
|
382 |
gezelter |
246 |
potatoes[which_node] = myPotato; |
383 |
|
|
} else { //master node has current molecule |
384 |
gezelter |
2 |
|
385 |
gezelter |
246 |
mol = info_->getMoleculeByGlobalIndex(i); |
386 |
gezelter |
2 |
|
387 |
gezelter |
246 |
if (mol == NULL) { |
388 |
|
|
sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank); |
389 |
|
|
painCave.isFatal = 1; |
390 |
|
|
simError(); |
391 |
|
|
} |
392 |
|
|
|
393 |
|
|
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
394 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
395 |
|
|
|
396 |
|
|
atomTypeString = integrableObject->getType().c_str(); |
397 |
gezelter |
2 |
|
398 |
gezelter |
246 |
pos = integrableObject->getPos(); |
399 |
|
|
vel = integrableObject->getVel(); |
400 |
gezelter |
2 |
|
401 |
gezelter |
246 |
atomData[0] = pos[0]; |
402 |
|
|
atomData[1] = pos[1]; |
403 |
|
|
atomData[2] = pos[2]; |
404 |
gezelter |
2 |
|
405 |
gezelter |
246 |
atomData[3] = vel[0]; |
406 |
|
|
atomData[4] = vel[1]; |
407 |
|
|
atomData[5] = vel[2]; |
408 |
gezelter |
2 |
|
409 |
gezelter |
246 |
isDirectional = 0; |
410 |
gezelter |
2 |
|
411 |
gezelter |
246 |
if (integrableObject->isDirectional()) { |
412 |
|
|
isDirectional = 1; |
413 |
gezelter |
2 |
|
414 |
gezelter |
246 |
q = integrableObject->getQ(); |
415 |
|
|
ji = integrableObject->getJ(); |
416 |
gezelter |
2 |
|
417 |
gezelter |
246 |
for(int j = 0; j < 6; j++) { |
418 |
|
|
atomData[j] = atomData[j]; |
419 |
|
|
} |
420 |
gezelter |
2 |
|
421 |
gezelter |
246 |
atomData[6] = q[0]; |
422 |
|
|
atomData[7] = q[1]; |
423 |
|
|
atomData[8] = q[2]; |
424 |
|
|
atomData[9] = q[3]; |
425 |
gezelter |
2 |
|
426 |
gezelter |
246 |
atomData[10] = ji[0]; |
427 |
|
|
atomData[11] = ji[1]; |
428 |
|
|
atomData[12] = ji[2]; |
429 |
|
|
} |
430 |
gezelter |
2 |
|
431 |
gezelter |
246 |
// If we've survived to here, format the line: |
432 |
gezelter |
2 |
|
433 |
gezelter |
246 |
if (!isDirectional) { |
434 |
|
|
sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t", |
435 |
|
|
atomTypeString, atomData[0], |
436 |
|
|
atomData[1], atomData[2], |
437 |
|
|
atomData[3], atomData[4], |
438 |
|
|
atomData[5]); |
439 |
gezelter |
2 |
|
440 |
gezelter |
246 |
strcat(writeLine, |
441 |
|
|
"0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n"); |
442 |
|
|
} else { |
443 |
|
|
sprintf(writeLine, |
444 |
|
|
"%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n", |
445 |
|
|
atomTypeString, |
446 |
|
|
atomData[0], |
447 |
|
|
atomData[1], |
448 |
|
|
atomData[2], |
449 |
|
|
atomData[3], |
450 |
|
|
atomData[4], |
451 |
|
|
atomData[5], |
452 |
|
|
atomData[6], |
453 |
|
|
atomData[7], |
454 |
|
|
atomData[8], |
455 |
|
|
atomData[9], |
456 |
|
|
atomData[10], |
457 |
|
|
atomData[11], |
458 |
|
|
atomData[12]); |
459 |
|
|
} |
460 |
gezelter |
2 |
|
461 |
|
|
|
462 |
gezelter |
246 |
os << writeLine; |
463 |
gezelter |
2 |
|
464 |
gezelter |
246 |
} //end for(iter = integrableObject.begin()) |
465 |
|
|
} |
466 |
|
|
} //end for(i = 0; i < mpiSim->getNmol()) |
467 |
gezelter |
2 |
|
468 |
gezelter |
246 |
os.flush(); |
469 |
tim |
251 |
|
470 |
gezelter |
246 |
sprintf(checkPointMsg, "Sucessfully took a dump.\n"); |
471 |
|
|
MPIcheckPoint(); |
472 |
gezelter |
2 |
|
473 |
gezelter |
246 |
delete [] potatoes; |
474 |
|
|
} else { |
475 |
gezelter |
2 |
|
476 |
gezelter |
246 |
// worldRank != 0, so I'm a remote node. |
477 |
gezelter |
2 |
|
478 |
gezelter |
246 |
// Set my magic potato to 0: |
479 |
gezelter |
2 |
|
480 |
gezelter |
246 |
myPotato = 0; |
481 |
gezelter |
2 |
|
482 |
gezelter |
246 |
for(int i = 0; i < info_->getNGlobalMolecules(); i++) { |
483 |
gezelter |
2 |
|
484 |
gezelter |
246 |
// Am I the node which has this integrableObject? |
485 |
|
|
int whichNode = info_->getMolToProc(i); |
486 |
|
|
if (whichNode == worldRank) { |
487 |
|
|
if (myPotato + 1 >= MAXTAG) { |
488 |
gezelter |
2 |
|
489 |
gezelter |
246 |
// The potato was going to exceed the maximum value, |
490 |
|
|
// so wrap this processor potato back to 0 (and block until |
491 |
|
|
// node 0 says we can go: |
492 |
gezelter |
2 |
|
493 |
gezelter |
246 |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, |
494 |
|
|
&istatus); |
495 |
|
|
} |
496 |
gezelter |
2 |
|
497 |
gezelter |
246 |
mol = info_->getMoleculeByGlobalIndex(i); |
498 |
|
|
|
499 |
gezelter |
2 |
|
500 |
gezelter |
246 |
nCurObj = mol->getNIntegrableObjects(); |
501 |
gezelter |
2 |
|
502 |
gezelter |
246 |
MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD); |
503 |
|
|
myPotato++; |
504 |
gezelter |
2 |
|
505 |
gezelter |
246 |
for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL; |
506 |
|
|
integrableObject = mol->nextIntegrableObject(ii)) { |
507 |
gezelter |
2 |
|
508 |
gezelter |
246 |
if (myPotato + 2 >= MAXTAG) { |
509 |
gezelter |
2 |
|
510 |
gezelter |
246 |
// The potato was going to exceed the maximum value, |
511 |
|
|
// so wrap this processor potato back to 0 (and block until |
512 |
|
|
// node 0 says we can go: |
513 |
gezelter |
2 |
|
514 |
gezelter |
246 |
MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD, |
515 |
|
|
&istatus); |
516 |
|
|
} |
517 |
gezelter |
2 |
|
518 |
gezelter |
246 |
atomTypeString = integrableObject->getType().c_str(); |
519 |
gezelter |
2 |
|
520 |
gezelter |
246 |
pos = integrableObject->getPos(); |
521 |
|
|
vel = integrableObject->getVel(); |
522 |
gezelter |
2 |
|
523 |
gezelter |
246 |
atomData[0] = pos[0]; |
524 |
|
|
atomData[1] = pos[1]; |
525 |
|
|
atomData[2] = pos[2]; |
526 |
gezelter |
2 |
|
527 |
gezelter |
246 |
atomData[3] = vel[0]; |
528 |
|
|
atomData[4] = vel[1]; |
529 |
|
|
atomData[5] = vel[2]; |
530 |
gezelter |
2 |
|
531 |
gezelter |
246 |
isDirectional = 0; |
532 |
gezelter |
2 |
|
533 |
gezelter |
246 |
if (integrableObject->isDirectional()) { |
534 |
|
|
isDirectional = 1; |
535 |
gezelter |
2 |
|
536 |
gezelter |
246 |
q = integrableObject->getQ(); |
537 |
|
|
ji = integrableObject->getJ(); |
538 |
gezelter |
2 |
|
539 |
gezelter |
246 |
atomData[6] = q[0]; |
540 |
|
|
atomData[7] = q[1]; |
541 |
|
|
atomData[8] = q[2]; |
542 |
|
|
atomData[9] = q[3]; |
543 |
|
|
|
544 |
|
|
atomData[10] = ji[0]; |
545 |
|
|
atomData[11] = ji[1]; |
546 |
|
|
atomData[12] = ji[2]; |
547 |
|
|
} |
548 |
|
|
|
549 |
|
|
strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE); |
550 |
|
|
|
551 |
|
|
// null terminate the std::string before sending (just in case): |
552 |
|
|
MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0'; |
553 |
|
|
|
554 |
|
|
MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0, |
555 |
|
|
myPotato, MPI_COMM_WORLD); |
556 |
|
|
|
557 |
|
|
myPotato++; |
558 |
|
|
|
559 |
|
|
if (isDirectional) { |
560 |
|
|
MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato, |
561 |
|
|
MPI_COMM_WORLD); |
562 |
|
|
} else { |
563 |
|
|
MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato, |
564 |
|
|
MPI_COMM_WORLD); |
565 |
|
|
} |
566 |
|
|
|
567 |
|
|
myPotato++; |
568 |
|
|
} |
569 |
|
|
|
570 |
|
|
} |
571 |
|
|
|
572 |
|
|
} |
573 |
|
|
sprintf(checkPointMsg, "Sucessfully took a dump.\n"); |
574 |
|
|
MPIcheckPoint(); |
575 |
|
|
} |
576 |
|
|
|
577 |
|
|
#endif // is_mpi |
578 |
|
|
|
579 |
gezelter |
2 |
} |
580 |
|
|
|
581 |
tim |
376 |
void DumpWriter::writeDump() { |
582 |
|
|
writeFrame(dumpFile_); |
583 |
|
|
|
584 |
|
|
} |
585 |
|
|
|
586 |
|
|
void DumpWriter::writeEor() { |
587 |
|
|
std::ofstream eorStream; |
588 |
|
|
|
589 |
|
|
#ifdef IS_MPI |
590 |
|
|
if (worldRank == 0) { |
591 |
|
|
#endif // is_mpi |
592 |
|
|
|
593 |
|
|
eorStream.open(eorFilename_.c_str()); |
594 |
|
|
if (!eorStream.is_open()) { |
595 |
|
|
sprintf(painCave.errMsg, "DumpWriter : Could not open \"%s\" for writing.\n", |
596 |
|
|
eorFilename_.c_str()); |
597 |
|
|
painCave.isFatal = 1; |
598 |
|
|
simError(); |
599 |
|
|
} |
600 |
|
|
|
601 |
|
|
#ifdef IS_MPI |
602 |
|
|
} |
603 |
|
|
#endif // is_mpi |
604 |
|
|
|
605 |
|
|
writeFrame(eorStream); |
606 |
|
|
} |
607 |
|
|
|
608 |
|
|
|
609 |
|
|
void DumpWriter::writeDumpAndEor() { |
610 |
|
|
std::ofstream eorStream; |
611 |
|
|
std::vector<std::streambuf*> buffers; |
612 |
|
|
#ifdef IS_MPI |
613 |
|
|
if (worldRank == 0) { |
614 |
|
|
#endif // is_mpi |
615 |
|
|
|
616 |
|
|
buffers.push_back(dumpFile_.rdbuf()); |
617 |
|
|
|
618 |
|
|
eorStream.open(eorFilename_.c_str()); |
619 |
|
|
if (!eorStream.is_open()) { |
620 |
|
|
sprintf(painCave.errMsg, "DumpWriter : Could not open \"%s\" for writing.\n", |
621 |
|
|
eorFilename_.c_str()); |
622 |
|
|
painCave.isFatal = 1; |
623 |
|
|
simError(); |
624 |
|
|
} |
625 |
|
|
|
626 |
|
|
buffers.push_back(eorStream.rdbuf()); |
627 |
|
|
|
628 |
|
|
#ifdef IS_MPI |
629 |
|
|
} |
630 |
|
|
#endif // is_mpi |
631 |
|
|
|
632 |
|
|
TeeBuf tbuf(buffers.begin(), buffers.end()); |
633 |
|
|
std::ostream os(&tbuf); |
634 |
|
|
|
635 |
|
|
writeFrame(os); |
636 |
|
|
|
637 |
|
|
} |
638 |
|
|
|
639 |
|
|
|
640 |
|
|
|
641 |
gezelter |
246 |
}//end namespace oopse |