ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/DumpWriter.cpp
Revision: 2425
Committed: Fri Nov 11 15:22:11 2005 UTC (18 years, 8 months ago) by chrisfen
File size: 20214 byte(s)
Log Message:
added in a 5th order polynomial switching function option

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
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
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() {
139
140 #ifdef IS_MPI
141
142 if (worldRank == 0) {
143 #endif // is_mpi
144
145 delete dumpFile_;
146
147 #ifdef IS_MPI
148
149 }
150
151 #endif // is_mpi
152
153 }
154
155 void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) {
156
157 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 << 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
174 //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 << 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
182 os << "\n";
183 }
184
185 void DumpWriter::writeFrame(std::ostream& os) {
186 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 Vector3d frc;
197 Vector3d trq;
198
199 Molecule* mol;
200 StuntDouble* integrableObject;
201 SimInfo::MoleculeIterator mi;
202 Molecule::IntegrableObjectIterator ii;
203
204 int nTotObjects;
205 nTotObjects = info_->getNGlobalIntegrableObjects();
206
207 #ifndef IS_MPI
208
209
210 os << nTotObjects << "\n";
211
212 writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
213
214 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
215
216 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
217 integrableObject = mol->nextIntegrableObject(ii)) {
218
219
220 pos = integrableObject->getPos();
221 vel = integrableObject->getVel();
222
223 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
228 strcpy(writeLine, tempBuffer);
229
230 if (integrableObject->isDirectional()) {
231 q = integrableObject->getQ();
232 ji = integrableObject->getJ();
233
234 sprintf(tempBuffer, "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf",
235 q[0], q[1], q[2], q[3],
236 ji[0], ji[1], ji[2]);
237 strcat(writeLine, tempBuffer);
238 } else {
239 strcat(writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
240 }
241
242 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 os << writeLine;
254
255 }
256 }
257
258 os.flush();
259 #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
298 int * potatoes;
299 int myPotato;
300 int nProc;
301 int which_node;
302 double atomData[19];
303 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
316 if (flag) {
317 MAXTAG = *tagub;
318 } else {
319 MAXTAG = 32767;
320 }
321
322 if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file
323
324 // Node 0 needs a list of the magic potatoes for each processor;
325
326 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
327 potatoes = new int[nProc];
328
329 //write out the comment lines
330 for(int i = 0; i < nProc; i++) {
331 potatoes[i] = 0;
332 }
333
334
335 os << nTotObjects << "\n";
336 writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
337
338 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
339
340 // Get the Node number which has this atom;
341
342 which_node = info_->getMolToProc(i);
343
344 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
349 potatoes[which_node] = 0;
350 MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0,
351 MPI_COMM_WORLD);
352 }
353
354 myPotato = potatoes[which_node];
355
356 //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
361 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
366 potatoes[which_node] = 0;
367 MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node,
368 0, MPI_COMM_WORLD);
369 }
370
371 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR,
372 which_node, myPotato, MPI_COMM_WORLD,
373 &istatus);
374
375 myPotato++;
376
377 MPI_Recv(atomData, 19, MPI_DOUBLE, which_node, myPotato,
378 MPI_COMM_WORLD, &istatus);
379 myPotato++;
380
381 MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen);
382
383 if (msgLen == 13 || msgLen == 19)
384 isDirectional = 1;
385 else
386 isDirectional = 0;
387
388 // If we've survived to here, format the line:
389
390 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
397 strcat(writeLine,
398 "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
399 } else {
400 sprintf(writeLine,
401 "%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 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
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
438 sprintf(writeLine, "\n");
439 os << writeLine;
440
441 } // end for(int l =0)
442
443 potatoes[which_node] = myPotato;
444 } else { //master node has current molecule
445
446 mol = info_->getMoleculeByGlobalIndex(i);
447
448 if (mol == NULL) {
449 sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
450 painCave.isFatal = 1;
451 simError();
452 }
453
454 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
455 integrableObject = mol->nextIntegrableObject(ii)) {
456
457 pos = integrableObject->getPos();
458 vel = integrableObject->getVel();
459
460 atomData[0] = pos[0];
461 atomData[1] = pos[1];
462 atomData[2] = pos[2];
463
464 atomData[3] = vel[0];
465 atomData[4] = vel[1];
466 atomData[5] = vel[2];
467
468 isDirectional = 0;
469
470 if (integrableObject->isDirectional()) {
471 isDirectional = 1;
472
473 q = integrableObject->getQ();
474 ji = integrableObject->getJ();
475
476 for(int j = 0; j < 6; j++) {
477 atomData[j] = atomData[j];
478 }
479
480 atomData[6] = q[0];
481 atomData[7] = q[1];
482 atomData[8] = q[2];
483 atomData[9] = q[3];
484
485 atomData[10] = ji[0];
486 atomData[11] = ji[1];
487 atomData[12] = ji[2];
488 }
489
490 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 // If we've survived to here, format the line:
512
513 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
520 strcat(writeLine,
521 "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0");
522 } else {
523 sprintf(writeLine,
524 "%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 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
541 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
561 sprintf(writeLine, "\n");
562 os << writeLine;
563
564 } //end for(iter = integrableObject.begin())
565 }
566 } //end for(i = 0; i < mpiSim->getNmol())
567
568 os.flush();
569
570 sprintf(checkPointMsg, "Sucessfully took a dump.\n");
571 MPIcheckPoint();
572
573 delete [] potatoes;
574 } else {
575
576 // worldRank != 0, so I'm a remote node.
577
578 // Set my magic potato to 0:
579
580 myPotato = 0;
581
582 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
583
584 // Am I the node which has this integrableObject?
585 int whichNode = info_->getMolToProc(i);
586 if (whichNode == worldRank) {
587 if (myPotato + 1 >= MAXTAG) {
588
589 // 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
593 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
594 &istatus);
595 }
596
597 mol = info_->getMoleculeByGlobalIndex(i);
598
599
600 nCurObj = mol->getNIntegrableObjects();
601
602 MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD);
603 myPotato++;
604
605 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
606 integrableObject = mol->nextIntegrableObject(ii)) {
607
608 if (myPotato + 2 >= MAXTAG) {
609
610 // 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
614 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
615 &istatus);
616 }
617
618 pos = integrableObject->getPos();
619 vel = integrableObject->getVel();
620
621 atomData[0] = pos[0];
622 atomData[1] = pos[1];
623 atomData[2] = pos[2];
624
625 atomData[3] = vel[0];
626 atomData[4] = vel[1];
627 atomData[5] = vel[2];
628
629 isDirectional = 0;
630
631 if (integrableObject->isDirectional()) {
632 isDirectional = 1;
633
634 q = integrableObject->getQ();
635 ji = integrableObject->getJ();
636
637 atomData[6] = q[0];
638 atomData[7] = q[1];
639 atomData[8] = q[2];
640 atomData[9] = q[3];
641
642 atomData[10] = ji[0];
643 atomData[11] = ji[1];
644 atomData[12] = ji[2];
645 }
646
647 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 strncpy(MPIatomTypeString, integrableObject->getType().c_str(), MINIBUFFERSIZE);
671
672 // null terminate the std::string before sending (just in case):
673 MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0';
674
675 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
676 myPotato, MPI_COMM_WORLD);
677
678 myPotato++;
679
680 if (isDirectional && needForceVector_) {
681 MPI_Send(atomData, 19, MPI_DOUBLE, 0, myPotato,
682 MPI_COMM_WORLD);
683 } else if (isDirectional) {
684 MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato,
685 MPI_COMM_WORLD);
686 } else if (needForceVector_) {
687 MPI_Send(atomData, 12, MPI_DOUBLE, 0, myPotato,
688 MPI_COMM_WORLD);
689 } else {
690 MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato,
691 MPI_COMM_WORLD);
692 }
693
694 myPotato++;
695 }
696
697 }
698
699 }
700 sprintf(checkPointMsg, "Sucessfully took a dump.\n");
701 MPIcheckPoint();
702 }
703
704 #endif // is_mpi
705
706 }
707
708 void DumpWriter::writeDump() {
709 writeFrame(*dumpFile_);
710 }
711
712 void DumpWriter::writeEor() {
713 std::ostream* eorStream;
714
715 #ifdef IS_MPI
716 if (worldRank == 0) {
717 #endif // is_mpi
718
719 eorStream = createOStream(eorFilename_);
720
721 #ifdef IS_MPI
722 }
723 #endif // is_mpi
724
725 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 }
737
738
739 void DumpWriter::writeDumpAndEor() {
740 std::vector<std::streambuf*> buffers;
741 std::ostream* eorStream;
742 #ifdef IS_MPI
743 if (worldRank == 0) {
744 #endif // is_mpi
745
746 buffers.push_back(dumpFile_->rdbuf());
747
748 eorStream = createOStream(eorFilename_);
749
750 buffers.push_back(eorStream->rdbuf());
751
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
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
770 }
771
772 std::ostream* DumpWriter::createOStream(const std::string& filename) {
773
774 std::ostream* newOStream;
775 #ifdef HAVE_LIBZ
776 if (needCompression_) {
777 newOStream = new ogzstream(filename.c_str());
778 } else {
779 newOStream = new std::ofstream(filename.c_str());
780 }
781 #else
782 newOStream = new std::ofstream(filename.c_str());
783 #endif
784 return newOStream;
785 }
786
787 }//end namespace oopse