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