ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-3.0/src/io/DumpWriter.cpp
Revision: 1935
Committed: Wed Jan 12 23:24:55 2005 UTC (19 years, 5 months ago) by tim
File size: 18909 byte(s)
Log Message:
remove useless debug info

File Contents

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