ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpWriter.cpp
Revision: 723
Committed: Tue Aug 26 20:12:51 2003 UTC (20 years, 10 months ago) by mmeineke
File size: 15074 byte(s)
Log Message:
added define statemewnt to Statwriter and Dumpwriter to handle files larger than 2 gb.

commented out some print statements in Zconstraint

File Contents

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