ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE-2.0/src/io/DumpWriter.cpp
Revision: 2079
Committed: Thu Mar 3 14:40:20 2005 UTC (19 years, 4 months ago) by tim
File size: 20729 byte(s)
Log Message:
avoid using const char*(hope can fixed the missing atom type problem)

File Contents

# User Rev Content
1 gezelter 1930 /*
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 2060 #include "io/basic_teebuf.hpp"
46 gezelter 1490 #ifdef IS_MPI
47     #include <mpi.h>
48 gezelter 1930 #endif //is_mpi
49 gezelter 1490
50 gezelter 1930 namespace oopse {
51 gezelter 1490
52 tim 2060 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 1930 DumpWriter::DumpWriter(SimInfo* info, const std::string& filename)
81     : info_(info), filename_(filename){
82     #ifdef IS_MPI
83 gezelter 1490
84 gezelter 1930 if (worldRank == 0) {
85     #endif // is_mpi
86 gezelter 1490
87 tim 2060 eorFilename_ = filename_.substr(0, filename_.rfind(".")) + ".eor";
88 gezelter 1930 dumpFile_.open(filename_.c_str(), std::ios::out | std::ios::trunc);
89 gezelter 1490
90 gezelter 1930 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 1490
97     #ifdef IS_MPI
98    
99     }
100    
101 gezelter 1930 sprintf(checkPointMsg, "Sucessfully opened output file for dumping.\n");
102     MPIcheckPoint();
103 gezelter 1490
104     #endif // is_mpi
105 gezelter 1930
106 gezelter 1490 }
107    
108 gezelter 1930 DumpWriter::~DumpWriter() {
109 gezelter 1490
110     #ifdef IS_MPI
111 gezelter 1930
112     if (worldRank == 0) {
113 gezelter 1490 #endif // is_mpi
114    
115 gezelter 1930 dumpFile_.close();
116 gezelter 1490
117     #ifdef IS_MPI
118 gezelter 1930
119     }
120    
121 gezelter 1490 #endif // is_mpi
122 gezelter 1930
123 gezelter 1490 }
124    
125 gezelter 1930 void DumpWriter::writeCommentLine(std::ostream& os, Snapshot* s) {
126 gezelter 1490
127 gezelter 1930 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 1490
144 gezelter 1930 //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 2060 os << "\n";
153 gezelter 1490 }
154    
155 gezelter 1930 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 1490
172 gezelter 1930 int nTotObjects;
173     nTotObjects = info_->getNGlobalIntegrableObjects();
174 gezelter 1490
175 gezelter 1930 #ifndef IS_MPI
176 gezelter 1490
177    
178 gezelter 1930 os << nTotObjects << "\n";
179    
180     writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
181 gezelter 1490
182 gezelter 1930 for (mol = info_->beginMolecule(mi); mol != NULL; mol = info_->nextMolecule(mi)) {
183 gezelter 1490
184 gezelter 1930 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
185     integrableObject = mol->nextIntegrableObject(ii)) {
186    
187 gezelter 1490
188 gezelter 1930 pos = integrableObject->getPos();
189     vel = integrableObject->getVel();
190 gezelter 1490
191 gezelter 1930 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 1490
196 gezelter 1930 strcpy(writeLine, tempBuffer);
197 gezelter 1490
198 gezelter 1930 if (integrableObject->isDirectional()) {
199     q = integrableObject->getQ();
200     ji = integrableObject->getJ();
201 gezelter 1490
202 gezelter 1930 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 1490
210 gezelter 1930 os << writeLine;
211 gezelter 1490
212 gezelter 1930 }
213 gezelter 1490 }
214    
215 tim 2008 os.flush();
216 gezelter 1930 #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 1490
255 gezelter 1930 int * potatoes;
256     int myPotato;
257     int nProc;
258     int which_node;
259     double atomData[13];
260     int isDirectional;
261     char MPIatomTypeString[MINIBUFFERSIZE];
262     int msgLen; // the length of message actually recieved at master nodes
263     int haveError;
264     MPI_Status istatus;
265     int nCurObj;
266    
267     // code to find maximum tag value
268     int * tagub;
269     int flag;
270     int MAXTAG;
271     MPI_Attr_get(MPI_COMM_WORLD, MPI_TAG_UB, &tagub, &flag);
272 gezelter 1490
273 gezelter 1930 if (flag) {
274     MAXTAG = *tagub;
275     } else {
276     MAXTAG = 32767;
277     }
278 gezelter 1490
279 gezelter 1930 if (worldRank == masterNode) { //master node (node 0) is responsible for writing the dump file
280 gezelter 1490
281 gezelter 1930 // Node 0 needs a list of the magic potatoes for each processor;
282 gezelter 1490
283 gezelter 1930 MPI_Comm_size(MPI_COMM_WORLD, &nProc);
284     potatoes = new int[nProc];
285 gezelter 1490
286 gezelter 1930 //write out the comment lines
287     for(int i = 0; i < nProc; i++) {
288     potatoes[i] = 0;
289     }
290 gezelter 1490
291    
292 gezelter 1930 os << nTotObjects << "\n";
293     writeCommentLine(os, info_->getSnapshotManager()->getCurrentSnapshot());
294 gezelter 1490
295 gezelter 1930 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
296 gezelter 1490
297 gezelter 1930 // Get the Node number which has this atom;
298 gezelter 1490
299 gezelter 1930 which_node = info_->getMolToProc(i);
300 gezelter 1490
301 gezelter 1930 if (which_node != masterNode) { //current molecule is in slave node
302     if (potatoes[which_node] + 1 >= MAXTAG) {
303     // The potato was going to exceed the maximum value,
304     // so wrap this processor potato back to 0:
305 gezelter 1490
306 gezelter 1930 potatoes[which_node] = 0;
307     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node, 0,
308     MPI_COMM_WORLD);
309     }
310 gezelter 1490
311 gezelter 1930 myPotato = potatoes[which_node];
312 gezelter 1490
313 gezelter 1930 //recieve the number of integrableObject in current molecule
314     MPI_Recv(&nCurObj, 1, MPI_INT, which_node, myPotato,
315     MPI_COMM_WORLD, &istatus);
316     myPotato++;
317 gezelter 1490
318 gezelter 1930 for(int l = 0; l < nCurObj; l++) {
319     if (potatoes[which_node] + 2 >= MAXTAG) {
320     // The potato was going to exceed the maximum value,
321     // so wrap this processor potato back to 0:
322 gezelter 1490
323 gezelter 1930 potatoes[which_node] = 0;
324     MPI_Send(&potatoes[which_node], 1, MPI_INT, which_node,
325     0, MPI_COMM_WORLD);
326     }
327 gezelter 1490
328 gezelter 1930 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR,
329     which_node, myPotato, MPI_COMM_WORLD,
330     &istatus);
331 gezelter 1490
332 gezelter 1930 myPotato++;
333 gezelter 1490
334 gezelter 1930 MPI_Recv(atomData, 13, MPI_DOUBLE, which_node, myPotato,
335     MPI_COMM_WORLD, &istatus);
336     myPotato++;
337 gezelter 1490
338 gezelter 1930 MPI_Get_count(&istatus, MPI_DOUBLE, &msgLen);
339 gezelter 1490
340 gezelter 1930 if (msgLen == 13)
341     isDirectional = 1;
342     else
343     isDirectional = 0;
344 gezelter 1490
345 gezelter 1930 // If we've survived to here, format the line:
346 gezelter 1490
347 gezelter 1930 if (!isDirectional) {
348     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
349 tim 2079 MPIatomTypeString, atomData[0],
350 gezelter 1930 atomData[1], atomData[2],
351     atomData[3], atomData[4],
352     atomData[5]);
353 gezelter 1490
354 gezelter 1930 strcat(writeLine,
355     "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n");
356     } else {
357     sprintf(writeLine,
358     "%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",
359 tim 2079 MPIatomTypeString,
360 gezelter 1930 atomData[0],
361     atomData[1],
362     atomData[2],
363     atomData[3],
364     atomData[4],
365     atomData[5],
366     atomData[6],
367     atomData[7],
368     atomData[8],
369     atomData[9],
370     atomData[10],
371     atomData[11],
372     atomData[12]);
373     }
374 gezelter 1490
375 gezelter 1930 os << writeLine;
376 gezelter 1490
377 gezelter 1930 } // end for(int l =0)
378 gezelter 1490
379 gezelter 1930 potatoes[which_node] = myPotato;
380     } else { //master node has current molecule
381 gezelter 1490
382 gezelter 1930 mol = info_->getMoleculeByGlobalIndex(i);
383 gezelter 1490
384 gezelter 1930 if (mol == NULL) {
385     sprintf(painCave.errMsg, "Molecule not found on node %d!", worldRank);
386     painCave.isFatal = 1;
387     simError();
388     }
389    
390     for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
391 tim 2079 integrableObject = mol->nextIntegrableObject(ii)) {
392 gezelter 1490
393 gezelter 1930 pos = integrableObject->getPos();
394     vel = integrableObject->getVel();
395 gezelter 1490
396 gezelter 1930 atomData[0] = pos[0];
397     atomData[1] = pos[1];
398     atomData[2] = pos[2];
399 gezelter 1490
400 gezelter 1930 atomData[3] = vel[0];
401     atomData[4] = vel[1];
402     atomData[5] = vel[2];
403 gezelter 1490
404 gezelter 1930 isDirectional = 0;
405 gezelter 1490
406 gezelter 1930 if (integrableObject->isDirectional()) {
407     isDirectional = 1;
408 gezelter 1490
409 gezelter 1930 q = integrableObject->getQ();
410     ji = integrableObject->getJ();
411 gezelter 1490
412 gezelter 1930 for(int j = 0; j < 6; j++) {
413     atomData[j] = atomData[j];
414     }
415 gezelter 1490
416 gezelter 1930 atomData[6] = q[0];
417     atomData[7] = q[1];
418     atomData[8] = q[2];
419     atomData[9] = q[3];
420 gezelter 1490
421 gezelter 1930 atomData[10] = ji[0];
422     atomData[11] = ji[1];
423     atomData[12] = ji[2];
424     }
425 gezelter 1490
426 gezelter 1930 // If we've survived to here, format the line:
427 gezelter 1490
428 gezelter 1930 if (!isDirectional) {
429     sprintf(writeLine, "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
430 tim 2079 integrableObject->getType().c_str(), atomData[0],
431 gezelter 1930 atomData[1], atomData[2],
432     atomData[3], atomData[4],
433     atomData[5]);
434 gezelter 1490
435 gezelter 1930 strcat(writeLine,
436     "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n");
437     } else {
438     sprintf(writeLine,
439     "%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",
440 tim 2079 integrableObject->getType().c_str(),
441 gezelter 1930 atomData[0],
442     atomData[1],
443     atomData[2],
444     atomData[3],
445     atomData[4],
446     atomData[5],
447     atomData[6],
448     atomData[7],
449     atomData[8],
450     atomData[9],
451     atomData[10],
452     atomData[11],
453     atomData[12]);
454     }
455 gezelter 1490
456    
457 gezelter 1930 os << writeLine;
458 gezelter 1490
459 gezelter 1930 } //end for(iter = integrableObject.begin())
460     }
461     } //end for(i = 0; i < mpiSim->getNmol())
462 gezelter 1490
463 gezelter 1930 os.flush();
464 tim 1935
465 gezelter 1930 sprintf(checkPointMsg, "Sucessfully took a dump.\n");
466     MPIcheckPoint();
467 gezelter 1490
468 gezelter 1930 delete [] potatoes;
469     } else {
470 gezelter 1490
471 gezelter 1930 // worldRank != 0, so I'm a remote node.
472 gezelter 1490
473 gezelter 1930 // Set my magic potato to 0:
474 gezelter 1490
475 gezelter 1930 myPotato = 0;
476 gezelter 1490
477 gezelter 1930 for(int i = 0; i < info_->getNGlobalMolecules(); i++) {
478 gezelter 1490
479 gezelter 1930 // Am I the node which has this integrableObject?
480     int whichNode = info_->getMolToProc(i);
481     if (whichNode == worldRank) {
482     if (myPotato + 1 >= MAXTAG) {
483 gezelter 1490
484 gezelter 1930 // The potato was going to exceed the maximum value,
485     // so wrap this processor potato back to 0 (and block until
486     // node 0 says we can go:
487 gezelter 1490
488 gezelter 1930 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
489     &istatus);
490     }
491 gezelter 1490
492 gezelter 1930 mol = info_->getMoleculeByGlobalIndex(i);
493    
494 gezelter 1490
495 gezelter 1930 nCurObj = mol->getNIntegrableObjects();
496 gezelter 1490
497 gezelter 1930 MPI_Send(&nCurObj, 1, MPI_INT, 0, myPotato, MPI_COMM_WORLD);
498     myPotato++;
499 gezelter 1490
500 gezelter 1930 for (integrableObject = mol->beginIntegrableObject(ii); integrableObject != NULL;
501     integrableObject = mol->nextIntegrableObject(ii)) {
502 gezelter 1490
503 gezelter 1930 if (myPotato + 2 >= MAXTAG) {
504 gezelter 1490
505 gezelter 1930 // The potato was going to exceed the maximum value,
506     // so wrap this processor potato back to 0 (and block until
507     // node 0 says we can go:
508 gezelter 1490
509 gezelter 1930 MPI_Recv(&myPotato, 1, MPI_INT, 0, 0, MPI_COMM_WORLD,
510     &istatus);
511     }
512 gezelter 1490
513 gezelter 1930 pos = integrableObject->getPos();
514     vel = integrableObject->getVel();
515 gezelter 1490
516 gezelter 1930 atomData[0] = pos[0];
517     atomData[1] = pos[1];
518     atomData[2] = pos[2];
519 gezelter 1490
520 gezelter 1930 atomData[3] = vel[0];
521     atomData[4] = vel[1];
522     atomData[5] = vel[2];
523 gezelter 1490
524 gezelter 1930 isDirectional = 0;
525 gezelter 1490
526 gezelter 1930 if (integrableObject->isDirectional()) {
527     isDirectional = 1;
528 gezelter 1490
529 gezelter 1930 q = integrableObject->getQ();
530     ji = integrableObject->getJ();
531 gezelter 1490
532 gezelter 1930 atomData[6] = q[0];
533     atomData[7] = q[1];
534     atomData[8] = q[2];
535     atomData[9] = q[3];
536    
537     atomData[10] = ji[0];
538     atomData[11] = ji[1];
539     atomData[12] = ji[2];
540     }
541    
542 tim 2079 strncpy(MPIatomTypeString, integrableObject->getType().c_str(), MINIBUFFERSIZE);
543 gezelter 1930
544     // null terminate the std::string before sending (just in case):
545     MPIatomTypeString[MINIBUFFERSIZE - 1] = '\0';
546    
547     MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
548     myPotato, MPI_COMM_WORLD);
549    
550     myPotato++;
551    
552     if (isDirectional) {
553     MPI_Send(atomData, 13, MPI_DOUBLE, 0, myPotato,
554     MPI_COMM_WORLD);
555     } else {
556     MPI_Send(atomData, 6, MPI_DOUBLE, 0, myPotato,
557     MPI_COMM_WORLD);
558     }
559    
560     myPotato++;
561     }
562    
563     }
564    
565     }
566     sprintf(checkPointMsg, "Sucessfully took a dump.\n");
567     MPIcheckPoint();
568     }
569    
570     #endif // is_mpi
571    
572 gezelter 1490 }
573    
574 tim 2060 void DumpWriter::writeDump() {
575     writeFrame(dumpFile_);
576    
577     }
578    
579     void DumpWriter::writeEor() {
580     std::ofstream eorStream;
581    
582     #ifdef IS_MPI
583     if (worldRank == 0) {
584     #endif // is_mpi
585    
586     eorStream.open(eorFilename_.c_str());
587     if (!eorStream.is_open()) {
588     sprintf(painCave.errMsg, "DumpWriter : Could not open \"%s\" for writing.\n",
589     eorFilename_.c_str());
590     painCave.isFatal = 1;
591     simError();
592     }
593    
594     #ifdef IS_MPI
595     }
596     #endif // is_mpi
597    
598     writeFrame(eorStream);
599     }
600    
601    
602     void DumpWriter::writeDumpAndEor() {
603     std::ofstream eorStream;
604     std::vector<std::streambuf*> buffers;
605     #ifdef IS_MPI
606     if (worldRank == 0) {
607     #endif // is_mpi
608    
609     buffers.push_back(dumpFile_.rdbuf());
610    
611     eorStream.open(eorFilename_.c_str());
612     if (!eorStream.is_open()) {
613     sprintf(painCave.errMsg, "DumpWriter : Could not open \"%s\" for writing.\n",
614     eorFilename_.c_str());
615     painCave.isFatal = 1;
616     simError();
617     }
618    
619     buffers.push_back(eorStream.rdbuf());
620    
621     #ifdef IS_MPI
622     }
623     #endif // is_mpi
624    
625     TeeBuf tbuf(buffers.begin(), buffers.end());
626     std::ostream os(&tbuf);
627    
628     writeFrame(os);
629    
630     }
631    
632    
633    
634 gezelter 1930 }//end namespace oopse