ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 910
Committed: Thu Jan 8 18:05:37 2004 UTC (20 years, 6 months ago) by gezelter
File size: 19466 byte(s)
Log Message:
added strncpy to DumpWriter

File Contents

# Content
1 #define _FILE_OFFSET_BITS 64
2
3 #include <string.h>
4 #include <iostream>
5 #include <fstream>
6
7 #ifdef IS_MPI
8 #include <mpi.h>
9 #include "mpiSimulation.hpp"
10
11 namespace dWrite{
12 void DieDieDie( void );
13 }
14
15 using namespace dWrite;
16 #endif //is_mpi
17
18 #include "ReadWrite.hpp"
19 #include "simError.h"
20
21 DumpWriter::DumpWriter( SimInfo* the_entry_plug ){
22
23 entry_plug = the_entry_plug;
24
25 #ifdef IS_MPI
26 if(worldRank == 0 ){
27 #endif // is_mpi
28
29 strcpy( outName, entry_plug->sampleName );
30
31 outFile.open(outName, ios::out | ios::trunc );
32
33 if( !outFile ){
34
35 sprintf( painCave.errMsg,
36 "Could not open \"%s\" for dump output.\n",
37 outName);
38 painCave.isFatal = 1;
39 simError();
40 }
41
42 //outFile.setf( ios::scientific );
43
44 #ifdef IS_MPI
45 }
46
47 sprintf( checkPointMsg,
48 "Sucessfully opened output file for dumping.\n");
49 MPIcheckPoint();
50 #endif // is_mpi
51 }
52
53 DumpWriter::~DumpWriter( ){
54
55 #ifdef IS_MPI
56 if(worldRank == 0 ){
57 #endif // is_mpi
58
59 outFile.close();
60
61 #ifdef IS_MPI
62 }
63 #endif // is_mpi
64 }
65
66 void DumpWriter::writeDump( double currentTime ){
67
68 const int BUFFERSIZE = 2000;
69 const int MINIBUFFERSIZE = 10;
70
71 char tempBuffer[BUFFERSIZE];
72 char writeLine[BUFFERSIZE];
73
74 int i;
75 #ifdef IS_MPI
76 int j, which_node, done, which_atom, local_index;
77 double atomTransData[6];
78 double atomOrientData[7];
79 int isDirectional;
80 char* atomTypeString;
81 char MPIatomTypeString[MINIBUFFERSIZE];
82 int me;
83 int atomTypeTag;
84 int atomIsDirectionalTag;
85 int atomTransDataTag;
86 int atomOrientDataTag;
87 #else //is_mpi
88 int nAtoms = entry_plug->n_atoms;
89 #endif //is_mpi
90
91 double q[4];
92 DirectionalAtom* dAtom;
93 Atom** atoms = entry_plug->atoms;
94 double pos[3], vel[3];
95
96 // write current frame to the eor file
97
98 this->writeFinal( currentTime );
99
100 #ifndef IS_MPI
101
102 outFile << nAtoms << "\n";
103
104 outFile << currentTime << ";\t"
105 << entry_plug->Hmat[0][0] << "\t"
106 << entry_plug->Hmat[1][0] << "\t"
107 << entry_plug->Hmat[2][0] << ";\t"
108
109 << entry_plug->Hmat[0][1] << "\t"
110 << entry_plug->Hmat[1][1] << "\t"
111 << entry_plug->Hmat[2][1] << ";\t"
112
113 << entry_plug->Hmat[0][2] << "\t"
114 << entry_plug->Hmat[1][2] << "\t"
115 << entry_plug->Hmat[2][2] << ";";
116 //write out additional parameters, such as chi and eta
117 outFile << entry_plug->the_integrator->getAdditionalParameters();
118 outFile << endl;
119
120 for( i=0; i<nAtoms; i++ ){
121
122 atoms[i]->getPos(pos);
123 atoms[i]->getVel(vel);
124
125 sprintf( tempBuffer,
126 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
127 atoms[i]->getType(),
128 pos[0],
129 pos[1],
130 pos[2],
131 vel[0],
132 vel[1],
133 vel[2]);
134 strcpy( writeLine, tempBuffer );
135
136 if( atoms[i]->isDirectional() ){
137
138 dAtom = (DirectionalAtom *)atoms[i];
139 dAtom->getQ( q );
140
141 sprintf( tempBuffer,
142 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
143 q[0],
144 q[1],
145 q[2],
146 q[3],
147 dAtom->getJx(),
148 dAtom->getJy(),
149 dAtom->getJz());
150 strcat( writeLine, tempBuffer );
151 }
152 else
153 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
154
155 outFile << writeLine;
156 }
157 outFile.flush();
158
159 #else // is_mpi
160
161 // first thing first, suspend fatalities.
162 painCave.isEventLoop = 1;
163
164 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
165 int haveError;
166
167 MPI_Status istatus;
168 int *AtomToProcMap = mpiSim->getAtomToProcMap();
169
170 // write out header and node 0's coordinates
171
172 if( worldRank == 0 ){
173 outFile << mpiSim->getTotAtoms() << "\n";
174
175 outFile << currentTime << ";\t"
176 << entry_plug->Hmat[0][0] << "\t"
177 << entry_plug->Hmat[1][0] << "\t"
178 << entry_plug->Hmat[2][0] << ";\t"
179
180 << entry_plug->Hmat[0][1] << "\t"
181 << entry_plug->Hmat[1][1] << "\t"
182 << entry_plug->Hmat[2][1] << ";\t"
183
184 << entry_plug->Hmat[0][2] << "\t"
185 << entry_plug->Hmat[1][2] << "\t"
186 << entry_plug->Hmat[2][2] << ";";
187
188 outFile << entry_plug->the_integrator->getAdditionalParameters();
189 outFile << endl;
190 outFile.flush();
191 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
192 // Get the Node number which has this atom;
193
194 which_node = AtomToProcMap[i];
195
196 if (which_node != 0) {
197
198 atomTypeTag = 4*i;
199 atomIsDirectionalTag = 4*i + 1;
200 atomTransDataTag = 4*i + 2;
201 atomOrientDataTag = 4*i + 3;
202
203 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
204 atomTypeTag, MPI_COMM_WORLD, &istatus);
205
206 strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
207
208 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
209 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
210
211 MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
212 atomTransDataTag, MPI_COMM_WORLD, &istatus);
213
214 if (isDirectional) {
215
216 MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
217 atomOrientDataTag, MPI_COMM_WORLD, &istatus);
218
219 }
220
221 } else {
222
223 haveError = 0;
224 which_atom = i;
225 local_index=-1;
226
227 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
228 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
229 }
230
231 if (local_index != -1) {
232
233 atomTypeString = atoms[local_index]->getType();
234
235 atoms[local_index]->getPos(pos);
236 atoms[local_index]->getVel(vel);
237
238 atomTransData[0] = pos[0];
239 atomTransData[1] = pos[1];
240 atomTransData[2] = pos[2];
241
242 atomTransData[3] = vel[0];
243 atomTransData[4] = vel[1];
244 atomTransData[5] = vel[2];
245
246 isDirectional = 0;
247
248 if( atoms[local_index]->isDirectional() ){
249
250 isDirectional = 1;
251
252 dAtom = (DirectionalAtom *)atoms[local_index];
253 dAtom->getQ( q );
254
255 atomOrientData[0] = q[0];
256 atomOrientData[1] = q[1];
257 atomOrientData[2] = q[2];
258 atomOrientData[3] = q[3];
259
260 atomOrientData[4] = dAtom->getJx();
261 atomOrientData[5] = dAtom->getJy();
262 atomOrientData[6] = dAtom->getJz();
263 }
264
265 } else {
266 sprintf(painCave.errMsg,
267 "Atom %d not found on processor %d\n",
268 i, worldRank );
269 haveError= 1;
270 simError();
271 }
272
273 if(haveError) DieDieDie();
274
275 // If we've survived to here, format the line:
276
277 sprintf( tempBuffer,
278 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
279 atomTypeString,
280 atomTransData[0],
281 atomTransData[1],
282 atomTransData[2],
283 atomTransData[3],
284 atomTransData[4],
285 atomTransData[5]);
286
287 strcpy( writeLine, tempBuffer );
288
289 if (isDirectional) {
290
291 sprintf( tempBuffer,
292 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
293 atomOrientData[0],
294 atomOrientData[1],
295 atomOrientData[2],
296 atomOrientData[3],
297 atomOrientData[4],
298 atomOrientData[5],
299 atomOrientData[6]);
300 strcat( writeLine, tempBuffer );
301
302 } else {
303 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
304 }
305
306 outFile << writeLine;
307 outFile.flush();
308 }
309 }
310
311 outFile.flush();
312 sprintf( checkPointMsg,
313 "Sucessfully took a dump.\n");
314 MPIcheckPoint();
315
316 } else {
317
318 // worldRank != 0, so I'm a remote node.
319
320 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
321
322 // Am I the node which has this atom?
323
324 if (AtomToProcMap[i] == worldRank) {
325
326 local_index=-1;
327 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
328 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
329 }
330 if (local_index != -1) {
331
332 atomTypeString = atoms[local_index]->getType();
333
334 atoms[local_index]->getPos(pos);
335 atoms[local_index]->getVel(vel);
336
337 atomTransData[0] = pos[0];
338 atomTransData[1] = pos[1];
339 atomTransData[2] = pos[2];
340
341 atomTransData[3] = vel[0];
342 atomTransData[4] = vel[1];
343 atomTransData[5] = vel[2];
344
345 isDirectional = 0;
346
347 if( atoms[local_index]->isDirectional() ){
348
349 isDirectional = 1;
350
351 dAtom = (DirectionalAtom *)atoms[local_index];
352 dAtom->getQ( q );
353
354 atomOrientData[0] = q[0];
355 atomOrientData[1] = q[1];
356 atomOrientData[2] = q[2];
357 atomOrientData[3] = q[3];
358
359 atomOrientData[4] = dAtom->getJx();
360 atomOrientData[5] = dAtom->getJy();
361 atomOrientData[6] = dAtom->getJz();
362 }
363
364 } else {
365 sprintf(painCave.errMsg,
366 "Atom %d not found on processor %d\n",
367 i, worldRank );
368 haveError= 1;
369 simError();
370 }
371
372 // I've survived this far, so send off the data!
373
374 atomTypeTag = 4*i;
375 atomIsDirectionalTag = 4*i + 1;
376 atomTransDataTag = 4*i + 2;
377 atomOrientDataTag = 4*i + 3;
378
379
380 strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
381
382 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
383 atomTypeTag, MPI_COMM_WORLD);
384
385 MPI_Send(&isDirectional, 1, MPI_INT, 0,
386 atomIsDirectionalTag, MPI_COMM_WORLD);
387
388 MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
389 atomTransDataTag, MPI_COMM_WORLD);
390
391 if (isDirectional) {
392
393 MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
394 atomOrientDataTag, MPI_COMM_WORLD);
395
396 }
397
398 }
399 }
400
401 sprintf( checkPointMsg,
402 "Sucessfully took a dump.\n");
403 MPIcheckPoint();
404
405 }
406
407 painCave.isEventLoop = 0;
408
409 #endif // is_mpi
410 }
411
412 void DumpWriter::writeFinal(double finalTime){
413
414 char finalName[500];
415 ofstream finalOut;
416
417 const int BUFFERSIZE = 2000;
418 const int MINIBUFFERSIZE = 10;
419 char tempBuffer[BUFFERSIZE];
420 char writeLine[BUFFERSIZE];
421
422 double q[4];
423 DirectionalAtom* dAtom;
424 Atom** atoms = entry_plug->atoms;
425 int i;
426 #ifdef IS_MPI
427 int j, which_node, done, which_atom, local_index;
428 double atomTransData[6];
429 double atomOrientData[7];
430 int isDirectional;
431 char* atomTypeString;
432 char MPIatomTypeString[MINIBUFFERSIZE];
433 int atomTypeTag;
434 int atomIsDirectionalTag;
435 int atomTransDataTag;
436 int atomOrientDataTag;
437 #else //is_mpi
438 int nAtoms = entry_plug->n_atoms;
439 #endif //is_mpi
440
441 double pos[3], vel[3];
442
443 #ifdef IS_MPI
444 if(worldRank == 0 ){
445 #endif // is_mpi
446
447 strcpy( finalName, entry_plug->finalName );
448
449 finalOut.open( finalName, ios::out | ios::trunc );
450 if( !finalOut ){
451 sprintf( painCave.errMsg,
452 "Could not open \"%s\" for final dump output.\n",
453 finalName );
454 painCave.isFatal = 1;
455 simError();
456 }
457
458 // finalOut.setf( ios::scientific );
459
460 #ifdef IS_MPI
461 }
462
463 sprintf(checkPointMsg,"Opened file for final configuration\n");
464 MPIcheckPoint();
465
466 #endif //is_mpi
467
468
469 #ifndef IS_MPI
470
471 finalOut << nAtoms << "\n";
472
473 finalOut << finalTime << ";\t"
474 << entry_plug->Hmat[0][0] << "\t"
475 << entry_plug->Hmat[1][0] << "\t"
476 << entry_plug->Hmat[2][0] << ";\t"
477
478 << entry_plug->Hmat[0][1] << "\t"
479 << entry_plug->Hmat[1][1] << "\t"
480 << entry_plug->Hmat[2][1] << ";\t"
481
482 << entry_plug->Hmat[0][2] << "\t"
483 << entry_plug->Hmat[1][2] << "\t"
484 << entry_plug->Hmat[2][2] << ";";
485
486 //write out additional parameters, such as chi and eta
487 finalOut << entry_plug->the_integrator->getAdditionalParameters();
488 finalOut << endl;
489
490 for( i=0; i<nAtoms; i++ ){
491
492 atoms[i]->getPos(pos);
493 atoms[i]->getVel(vel);
494
495 sprintf( tempBuffer,
496 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
497 atoms[i]->getType(),
498 pos[0],
499 pos[1],
500 pos[2],
501 vel[0],
502 vel[1],
503 vel[2]);
504 strcpy( writeLine, tempBuffer );
505
506 if( atoms[i]->isDirectional() ){
507
508 dAtom = (DirectionalAtom *)atoms[i];
509 dAtom->getQ( q );
510
511 sprintf( tempBuffer,
512 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
513 q[0],
514 q[1],
515 q[2],
516 q[3],
517 dAtom->getJx(),
518 dAtom->getJy(),
519 dAtom->getJz());
520 strcat( writeLine, tempBuffer );
521 }
522 else
523 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
524
525 finalOut << writeLine;
526 }
527 finalOut.flush();
528 finalOut.close();
529
530 #else // is_mpi
531
532 // first thing first, suspend fatalities.
533 painCave.isEventLoop = 1;
534
535 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
536 int haveError;
537
538 MPI_Status istatus;
539 int *AtomToProcMap = mpiSim->getAtomToProcMap();
540
541 // write out header and node 0's coordinates
542
543 if( worldRank == 0 ){
544 finalOut << mpiSim->getTotAtoms() << "\n";
545
546 finalOut << finalTime << ";\t"
547 << entry_plug->Hmat[0][0] << "\t"
548 << entry_plug->Hmat[1][0] << "\t"
549 << entry_plug->Hmat[2][0] << ";\t"
550
551 << entry_plug->Hmat[0][1] << "\t"
552 << entry_plug->Hmat[1][1] << "\t"
553 << entry_plug->Hmat[2][1] << ";\t"
554
555 << entry_plug->Hmat[0][2] << "\t"
556 << entry_plug->Hmat[1][2] << "\t"
557 << entry_plug->Hmat[2][2] << ";";
558
559 finalOut << entry_plug->the_integrator->getAdditionalParameters();
560 finalOut << endl;
561 finalOut.flush();
562 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
563 // Get the Node number which has this atom;
564
565 which_node = AtomToProcMap[i];
566
567 if (which_node != 0) {
568
569 atomTypeTag = 4*i;
570 atomIsDirectionalTag = 4*i + 1;
571 atomTransDataTag = 4*i + 2;
572 atomOrientDataTag = 4*i + 3;
573
574 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
575 atomTypeTag, MPI_COMM_WORLD, &istatus);
576
577 strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
578
579 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
580 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
581
582 MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
583 atomTransDataTag, MPI_COMM_WORLD, &istatus);
584
585 if (isDirectional) {
586
587 MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
588 atomOrientDataTag, MPI_COMM_WORLD, &istatus);
589
590 }
591
592 } else {
593
594 haveError = 0;
595 which_atom = i;
596 local_index=-1;
597
598 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
599 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
600 }
601
602 if (local_index != -1) {
603
604 atomTypeString = atoms[local_index]->getType();
605
606 atoms[local_index]->getPos(pos);
607 atoms[local_index]->getVel(vel);
608
609 atomTransData[0] = pos[0];
610 atomTransData[1] = pos[1];
611 atomTransData[2] = pos[2];
612
613 atomTransData[3] = vel[0];
614 atomTransData[4] = vel[1];
615 atomTransData[5] = vel[2];
616
617 isDirectional = 0;
618
619 if( atoms[local_index]->isDirectional() ){
620
621 isDirectional = 1;
622
623 dAtom = (DirectionalAtom *)atoms[local_index];
624 dAtom->getQ( q );
625
626 atomOrientData[0] = q[0];
627 atomOrientData[1] = q[1];
628 atomOrientData[2] = q[2];
629 atomOrientData[3] = q[3];
630
631 atomOrientData[4] = dAtom->getJx();
632 atomOrientData[5] = dAtom->getJy();
633 atomOrientData[6] = dAtom->getJz();
634 }
635
636 } else {
637 sprintf(painCave.errMsg,
638 "Atom %d not found on processor %d\n",
639 i, worldRank );
640 haveError= 1;
641 simError();
642 }
643
644 if(haveError) DieDieDie();
645
646 // If we've survived to here, format the line:
647
648 sprintf( tempBuffer,
649 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
650 atomTypeString,
651 atomTransData[0],
652 atomTransData[1],
653 atomTransData[2],
654 atomTransData[3],
655 atomTransData[4],
656 atomTransData[5]);
657
658 strcpy( writeLine, tempBuffer );
659
660 if (isDirectional) {
661
662 sprintf( tempBuffer,
663 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
664 atomOrientData[0],
665 atomOrientData[1],
666 atomOrientData[2],
667 atomOrientData[3],
668 atomOrientData[4],
669 atomOrientData[5],
670 atomOrientData[6]);
671 strcat( writeLine, tempBuffer );
672
673 } else {
674 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
675 }
676
677 finalOut << writeLine;
678 finalOut.flush();
679 }
680 }
681
682 finalOut.flush();
683 sprintf( checkPointMsg,
684 "Sucessfully took a dump.\n");
685 MPIcheckPoint();
686
687 } else {
688
689 // worldRank != 0, so I'm a remote node.
690
691 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
692
693 // Am I the node which has this atom?
694
695 if (AtomToProcMap[i] == worldRank) {
696
697 local_index=-1;
698 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
699 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
700 }
701 if (local_index != -1) {
702
703 atomTypeString = atoms[local_index]->getType();
704
705 atoms[local_index]->getPos(pos);
706 atoms[local_index]->getVel(vel);
707
708 atomTransData[0] = pos[0];
709 atomTransData[1] = pos[1];
710 atomTransData[2] = pos[2];
711
712 atomTransData[3] = vel[0];
713 atomTransData[4] = vel[1];
714 atomTransData[5] = vel[2];
715
716 isDirectional = 0;
717
718 if( atoms[local_index]->isDirectional() ){
719
720 isDirectional = 1;
721
722 dAtom = (DirectionalAtom *)atoms[local_index];
723 dAtom->getQ( q );
724
725 atomOrientData[0] = q[0];
726 atomOrientData[1] = q[1];
727 atomOrientData[2] = q[2];
728 atomOrientData[3] = q[3];
729
730 atomOrientData[4] = dAtom->getJx();
731 atomOrientData[5] = dAtom->getJy();
732 atomOrientData[6] = dAtom->getJz();
733 }
734
735 } else {
736 sprintf(painCave.errMsg,
737 "Atom %d not found on processor %d\n",
738 i, worldRank );
739 haveError= 1;
740 simError();
741 }
742
743 // I've survived this far, so send off the data!
744
745 atomTypeTag = 4*i;
746 atomIsDirectionalTag = 4*i + 1;
747 atomTransDataTag = 4*i + 2;
748 atomOrientDataTag = 4*i + 3;
749
750 strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
751
752 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
753 atomTypeTag, MPI_COMM_WORLD);
754
755 MPI_Send(&isDirectional, 1, MPI_INT, 0,
756 atomIsDirectionalTag, MPI_COMM_WORLD);
757
758 MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
759 atomTransDataTag, MPI_COMM_WORLD);
760
761 if (isDirectional) {
762
763 MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
764 atomOrientDataTag, MPI_COMM_WORLD);
765
766 }
767
768 }
769 }
770
771 sprintf( checkPointMsg,
772 "Sucessfully wrote final file.\n");
773 MPIcheckPoint();
774
775 }
776
777 painCave.isEventLoop = 0;
778
779 if( worldRank == 0 ) finalOut.close();
780 #endif // is_mpi
781 }
782
783
784
785 #ifdef IS_MPI
786
787 // a couple of functions to let us escape the write loop
788
789 void dWrite::DieDieDie( void ){
790
791 MPI_Finalize();
792 exit (0);
793 }
794
795 #endif //is_mpi