ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-4/src/io/DumpWriter.cpp
Revision: 2761
Committed: Fri May 19 20:45:55 2006 UTC (18 years, 11 months ago) by tim
File size: 21397 byte(s)
Log Message:
Fix a bug of printing empty line in DumpWriter

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