ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 913
Committed: Thu Jan 8 22:25:52 2004 UTC (20 years, 6 months ago) by chuckv
File size: 19877 byte(s)
Log Message:
A work in progress...

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