ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 912
Committed: Thu Jan 8 18:59:36 2004 UTC (20 years, 6 months ago) by gezelter
File size: 19717 byte(s)
Log Message:
null terminate some strings just in case

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 = 100;
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 // Null terminate the atomTypeString just in case:
209
210 atomTypeString[strlen(atomTypeString) - 1] = '\0';
211
212 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
213 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
214
215 MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
216 atomTransDataTag, MPI_COMM_WORLD, &istatus);
217
218 if (isDirectional) {
219
220 MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
221 atomOrientDataTag, MPI_COMM_WORLD, &istatus);
222
223 }
224
225 } else {
226
227 haveError = 0;
228 which_atom = i;
229 local_index=-1;
230
231 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
232 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
233 }
234
235 if (local_index != -1) {
236
237 atomTypeString = atoms[local_index]->getType();
238
239 atoms[local_index]->getPos(pos);
240 atoms[local_index]->getVel(vel);
241
242 atomTransData[0] = pos[0];
243 atomTransData[1] = pos[1];
244 atomTransData[2] = pos[2];
245
246 atomTransData[3] = vel[0];
247 atomTransData[4] = vel[1];
248 atomTransData[5] = vel[2];
249
250 isDirectional = 0;
251
252 if( atoms[local_index]->isDirectional() ){
253
254 isDirectional = 1;
255
256 dAtom = (DirectionalAtom *)atoms[local_index];
257 dAtom->getQ( q );
258
259 atomOrientData[0] = q[0];
260 atomOrientData[1] = q[1];
261 atomOrientData[2] = q[2];
262 atomOrientData[3] = q[3];
263
264 atomOrientData[4] = dAtom->getJx();
265 atomOrientData[5] = dAtom->getJy();
266 atomOrientData[6] = dAtom->getJz();
267 }
268
269 } else {
270 sprintf(painCave.errMsg,
271 "Atom %d not found on processor %d\n",
272 i, worldRank );
273 haveError= 1;
274 simError();
275 }
276
277 if(haveError) DieDieDie();
278
279 // If we've survived to here, format the line:
280
281 sprintf( tempBuffer,
282 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
283 atomTypeString,
284 atomTransData[0],
285 atomTransData[1],
286 atomTransData[2],
287 atomTransData[3],
288 atomTransData[4],
289 atomTransData[5]);
290
291 strcpy( writeLine, tempBuffer );
292
293 if (isDirectional) {
294
295 sprintf( tempBuffer,
296 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
297 atomOrientData[0],
298 atomOrientData[1],
299 atomOrientData[2],
300 atomOrientData[3],
301 atomOrientData[4],
302 atomOrientData[5],
303 atomOrientData[6]);
304 strcat( writeLine, tempBuffer );
305
306 } else {
307 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
308 }
309
310 outFile << writeLine;
311 outFile.flush();
312 }
313 }
314
315 outFile.flush();
316 sprintf( checkPointMsg,
317 "Sucessfully took a dump.\n");
318 MPIcheckPoint();
319
320 } else {
321
322 // worldRank != 0, so I'm a remote node.
323
324 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
325
326 // Am I the node which has this atom?
327
328 if (AtomToProcMap[i] == worldRank) {
329
330 local_index=-1;
331 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
332 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
333 }
334 if (local_index != -1) {
335
336 atomTypeString = atoms[local_index]->getType();
337
338 atoms[local_index]->getPos(pos);
339 atoms[local_index]->getVel(vel);
340
341 atomTransData[0] = pos[0];
342 atomTransData[1] = pos[1];
343 atomTransData[2] = pos[2];
344
345 atomTransData[3] = vel[0];
346 atomTransData[4] = vel[1];
347 atomTransData[5] = vel[2];
348
349 isDirectional = 0;
350
351 if( atoms[local_index]->isDirectional() ){
352
353 isDirectional = 1;
354
355 dAtom = (DirectionalAtom *)atoms[local_index];
356 dAtom->getQ( q );
357
358 atomOrientData[0] = q[0];
359 atomOrientData[1] = q[1];
360 atomOrientData[2] = q[2];
361 atomOrientData[3] = q[3];
362
363 atomOrientData[4] = dAtom->getJx();
364 atomOrientData[5] = dAtom->getJy();
365 atomOrientData[6] = dAtom->getJz();
366 }
367
368 } else {
369 sprintf(painCave.errMsg,
370 "Atom %d not found on processor %d\n",
371 i, worldRank );
372 haveError= 1;
373 simError();
374 }
375
376 // I've survived this far, so send off the data!
377
378 atomTypeTag = 4*i;
379 atomIsDirectionalTag = 4*i + 1;
380 atomTransDataTag = 4*i + 2;
381 atomOrientDataTag = 4*i + 3;
382
383
384 strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
385
386 // null terminate the string before sending (just in case):
387 MPIatomTypeString[MINIBUFFERSIZE-1] = '\0';
388
389 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
390 atomTypeTag, MPI_COMM_WORLD);
391
392 MPI_Send(&isDirectional, 1, MPI_INT, 0,
393 atomIsDirectionalTag, MPI_COMM_WORLD);
394
395 MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
396 atomTransDataTag, MPI_COMM_WORLD);
397
398 if (isDirectional) {
399
400 MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
401 atomOrientDataTag, MPI_COMM_WORLD);
402
403 }
404
405 }
406 }
407
408 sprintf( checkPointMsg,
409 "Sucessfully took a dump.\n");
410 MPIcheckPoint();
411
412 }
413
414 painCave.isEventLoop = 0;
415
416 #endif // is_mpi
417 }
418
419 void DumpWriter::writeFinal(double finalTime){
420
421 char finalName[500];
422 ofstream finalOut;
423
424 const int BUFFERSIZE = 2000;
425 const int MINIBUFFERSIZE = 100;
426 char tempBuffer[BUFFERSIZE];
427 char writeLine[BUFFERSIZE];
428
429 double q[4];
430 DirectionalAtom* dAtom;
431 Atom** atoms = entry_plug->atoms;
432 int i;
433 #ifdef IS_MPI
434 int j, which_node, done, which_atom, local_index;
435 double atomTransData[6];
436 double atomOrientData[7];
437 int isDirectional;
438 char* atomTypeString;
439 char MPIatomTypeString[MINIBUFFERSIZE];
440 int atomTypeTag;
441 int atomIsDirectionalTag;
442 int atomTransDataTag;
443 int atomOrientDataTag;
444 #else //is_mpi
445 int nAtoms = entry_plug->n_atoms;
446 #endif //is_mpi
447
448 double pos[3], vel[3];
449
450 #ifdef IS_MPI
451 if(worldRank == 0 ){
452 #endif // is_mpi
453
454 strcpy( finalName, entry_plug->finalName );
455
456 finalOut.open( finalName, ios::out | ios::trunc );
457 if( !finalOut ){
458 sprintf( painCave.errMsg,
459 "Could not open \"%s\" for final dump output.\n",
460 finalName );
461 painCave.isFatal = 1;
462 simError();
463 }
464
465 // finalOut.setf( ios::scientific );
466
467 #ifdef IS_MPI
468 }
469
470 sprintf(checkPointMsg,"Opened file for final configuration\n");
471 MPIcheckPoint();
472
473 #endif //is_mpi
474
475
476 #ifndef IS_MPI
477
478 finalOut << nAtoms << "\n";
479
480 finalOut << finalTime << ";\t"
481 << entry_plug->Hmat[0][0] << "\t"
482 << entry_plug->Hmat[1][0] << "\t"
483 << entry_plug->Hmat[2][0] << ";\t"
484
485 << entry_plug->Hmat[0][1] << "\t"
486 << entry_plug->Hmat[1][1] << "\t"
487 << entry_plug->Hmat[2][1] << ";\t"
488
489 << entry_plug->Hmat[0][2] << "\t"
490 << entry_plug->Hmat[1][2] << "\t"
491 << entry_plug->Hmat[2][2] << ";";
492
493 //write out additional parameters, such as chi and eta
494 finalOut << entry_plug->the_integrator->getAdditionalParameters();
495 finalOut << endl;
496
497 for( i=0; i<nAtoms; i++ ){
498
499 atoms[i]->getPos(pos);
500 atoms[i]->getVel(vel);
501
502 sprintf( tempBuffer,
503 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
504 atoms[i]->getType(),
505 pos[0],
506 pos[1],
507 pos[2],
508 vel[0],
509 vel[1],
510 vel[2]);
511 strcpy( writeLine, tempBuffer );
512
513 if( atoms[i]->isDirectional() ){
514
515 dAtom = (DirectionalAtom *)atoms[i];
516 dAtom->getQ( q );
517
518 sprintf( tempBuffer,
519 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
520 q[0],
521 q[1],
522 q[2],
523 q[3],
524 dAtom->getJx(),
525 dAtom->getJy(),
526 dAtom->getJz());
527 strcat( writeLine, tempBuffer );
528 }
529 else
530 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
531
532 finalOut << writeLine;
533 }
534 finalOut.flush();
535 finalOut.close();
536
537 #else // is_mpi
538
539 // first thing first, suspend fatalities.
540 painCave.isEventLoop = 1;
541
542 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
543 int haveError;
544
545 MPI_Status istatus;
546 int *AtomToProcMap = mpiSim->getAtomToProcMap();
547
548 // write out header and node 0's coordinates
549
550 if( worldRank == 0 ){
551 finalOut << mpiSim->getTotAtoms() << "\n";
552
553 finalOut << finalTime << ";\t"
554 << entry_plug->Hmat[0][0] << "\t"
555 << entry_plug->Hmat[1][0] << "\t"
556 << entry_plug->Hmat[2][0] << ";\t"
557
558 << entry_plug->Hmat[0][1] << "\t"
559 << entry_plug->Hmat[1][1] << "\t"
560 << entry_plug->Hmat[2][1] << ";\t"
561
562 << entry_plug->Hmat[0][2] << "\t"
563 << entry_plug->Hmat[1][2] << "\t"
564 << entry_plug->Hmat[2][2] << ";";
565
566 finalOut << entry_plug->the_integrator->getAdditionalParameters();
567 finalOut << endl;
568 finalOut.flush();
569 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
570 // Get the Node number which has this atom;
571
572 which_node = AtomToProcMap[i];
573
574 if (which_node != 0) {
575
576 atomTypeTag = 4*i;
577 atomIsDirectionalTag = 4*i + 1;
578 atomTransDataTag = 4*i + 2;
579 atomOrientDataTag = 4*i + 3;
580
581 MPI_Recv(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, which_node,
582 atomTypeTag, MPI_COMM_WORLD, &istatus);
583
584 strncpy(atomTypeString, MPIatomTypeString, MINIBUFFERSIZE);
585
586 MPI_Recv(&isDirectional, 1, MPI_INT, which_node,
587 atomIsDirectionalTag, MPI_COMM_WORLD, &istatus);
588
589 MPI_Recv(atomTransData, 6, MPI_DOUBLE, which_node,
590 atomTransDataTag, MPI_COMM_WORLD, &istatus);
591
592 if (isDirectional) {
593
594 MPI_Recv(atomOrientData, 7, MPI_DOUBLE, which_node,
595 atomOrientDataTag, MPI_COMM_WORLD, &istatus);
596
597 }
598
599 } else {
600
601 haveError = 0;
602 which_atom = i;
603 local_index=-1;
604
605 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
606 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
607 }
608
609 if (local_index != -1) {
610
611 atomTypeString = atoms[local_index]->getType();
612
613 atoms[local_index]->getPos(pos);
614 atoms[local_index]->getVel(vel);
615
616 atomTransData[0] = pos[0];
617 atomTransData[1] = pos[1];
618 atomTransData[2] = pos[2];
619
620 atomTransData[3] = vel[0];
621 atomTransData[4] = vel[1];
622 atomTransData[5] = vel[2];
623
624 isDirectional = 0;
625
626 if( atoms[local_index]->isDirectional() ){
627
628 isDirectional = 1;
629
630 dAtom = (DirectionalAtom *)atoms[local_index];
631 dAtom->getQ( q );
632
633 atomOrientData[0] = q[0];
634 atomOrientData[1] = q[1];
635 atomOrientData[2] = q[2];
636 atomOrientData[3] = q[3];
637
638 atomOrientData[4] = dAtom->getJx();
639 atomOrientData[5] = dAtom->getJy();
640 atomOrientData[6] = dAtom->getJz();
641 }
642
643 } else {
644 sprintf(painCave.errMsg,
645 "Atom %d not found on processor %d\n",
646 i, worldRank );
647 haveError= 1;
648 simError();
649 }
650
651 if(haveError) DieDieDie();
652
653 // If we've survived to here, format the line:
654
655 sprintf( tempBuffer,
656 "%s\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t",
657 atomTypeString,
658 atomTransData[0],
659 atomTransData[1],
660 atomTransData[2],
661 atomTransData[3],
662 atomTransData[4],
663 atomTransData[5]);
664
665 strcpy( writeLine, tempBuffer );
666
667 if (isDirectional) {
668
669 sprintf( tempBuffer,
670 "%lf\t%lf\t%lf\t%lf\t%lf\t%lf\t%lf\n",
671 atomOrientData[0],
672 atomOrientData[1],
673 atomOrientData[2],
674 atomOrientData[3],
675 atomOrientData[4],
676 atomOrientData[5],
677 atomOrientData[6]);
678 strcat( writeLine, tempBuffer );
679
680 } else {
681 strcat( writeLine, "0.0\t0.0\t0.0\t0.0\t0.0\t0.0\t0.0\n" );
682 }
683
684 finalOut << writeLine;
685 finalOut.flush();
686 }
687 }
688
689 finalOut.flush();
690 sprintf( checkPointMsg,
691 "Sucessfully took a dump.\n");
692 MPIcheckPoint();
693
694 } else {
695
696 // worldRank != 0, so I'm a remote node.
697
698 for (i = 0 ; i < mpiSim->getTotAtoms(); i++ ) {
699
700 // Am I the node which has this atom?
701
702 if (AtomToProcMap[i] == worldRank) {
703
704 local_index=-1;
705 for (j=0; (j<mpiSim->getMyNlocal()) && (local_index < 0); j++) {
706 if (atoms[j]->getGlobalIndex() == which_atom) local_index = j;
707 }
708 if (local_index != -1) {
709
710 atomTypeString = atoms[local_index]->getType();
711
712 atoms[local_index]->getPos(pos);
713 atoms[local_index]->getVel(vel);
714
715 atomTransData[0] = pos[0];
716 atomTransData[1] = pos[1];
717 atomTransData[2] = pos[2];
718
719 atomTransData[3] = vel[0];
720 atomTransData[4] = vel[1];
721 atomTransData[5] = vel[2];
722
723 isDirectional = 0;
724
725 if( atoms[local_index]->isDirectional() ){
726
727 isDirectional = 1;
728
729 dAtom = (DirectionalAtom *)atoms[local_index];
730 dAtom->getQ( q );
731
732 atomOrientData[0] = q[0];
733 atomOrientData[1] = q[1];
734 atomOrientData[2] = q[2];
735 atomOrientData[3] = q[3];
736
737 atomOrientData[4] = dAtom->getJx();
738 atomOrientData[5] = dAtom->getJy();
739 atomOrientData[6] = dAtom->getJz();
740 }
741
742 } else {
743 sprintf(painCave.errMsg,
744 "Atom %d not found on processor %d\n",
745 i, worldRank );
746 haveError= 1;
747 simError();
748 }
749
750 // I've survived this far, so send off the data!
751
752 atomTypeTag = 4*i;
753 atomIsDirectionalTag = 4*i + 1;
754 atomTransDataTag = 4*i + 2;
755 atomOrientDataTag = 4*i + 3;
756
757 strncpy(MPIatomTypeString, atomTypeString, MINIBUFFERSIZE);
758
759 MPI_Send(MPIatomTypeString, MINIBUFFERSIZE, MPI_CHAR, 0,
760 atomTypeTag, MPI_COMM_WORLD);
761
762 MPI_Send(&isDirectional, 1, MPI_INT, 0,
763 atomIsDirectionalTag, MPI_COMM_WORLD);
764
765 MPI_Send(atomTransData, 6, MPI_DOUBLE, 0,
766 atomTransDataTag, MPI_COMM_WORLD);
767
768 if (isDirectional) {
769
770 MPI_Send(atomOrientData, 7, MPI_DOUBLE, 0,
771 atomOrientDataTag, MPI_COMM_WORLD);
772
773 }
774
775 }
776 }
777
778 sprintf( checkPointMsg,
779 "Sucessfully wrote final file.\n");
780 MPIcheckPoint();
781
782 }
783
784 painCave.isEventLoop = 0;
785
786 if( worldRank == 0 ) finalOut.close();
787 #endif // is_mpi
788 }
789
790
791
792 #ifdef IS_MPI
793
794 // a couple of functions to let us escape the write loop
795
796 void dWrite::DieDieDie( void ){
797
798 MPI_Finalize();
799 exit (0);
800 }
801
802 #endif //is_mpi