ViewVC Help
View File | Revision Log | Show Annotations | View Changeset | Root Listing
root/group/trunk/OOPSE/libmdtools/DumpReader.cpp
Revision: 647
Committed: Tue Jul 22 21:49:38 2003 UTC (20 years, 11 months ago) by mmeineke
File size: 19028 byte(s)
Log Message:
added the scan function to the DumpReader. It should now save the start of each frame in a vector.

File Contents

# Content
1 #define _FILE_OFFSET_BITS 64
2
3 #include <iostream>
4 #include <cmath>
5 #include <vector>
6
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include <string.h>
10 #include <unistd.h>
11 #include <sys/types.h>
12 #include <sys/stat.h>
13
14 #include "ReadWrite.hpp"
15 #include "simError.h"
16
17 #ifdef IS_MPI
18 #include <mpi.h>
19 #include "mpiSimulation.hpp"
20 #define TAKE_THIS_TAG_CHAR 0
21 #define TAKE_THIS_TAG_INT 1
22 #endif // is_mpi
23
24 namespace dumpRead{
25
26 #ifdef IS_MPI
27 void nodeZeroError( void );
28 void anonymousNodeDie( void );
29 #endif // is_mpi
30
31 class FilePos{
32
33 public:
34 FilePos(){ myPos = NULL; }
35 FilePos( fpos_t* thePos ) { myPos = thePos; }
36 ~FilePos(){ if( myPos != NULL ) delete myPos; }
37
38 FilePos &operator=(fpos_t *thePos){ myPos = thePos; return *this; }
39
40 void setPos( fpos_t *thePos ){ myPos = thePos; }
41 fpos_t *getPos( void ){ return myPos; }
42
43 private:
44
45 fpos_t *myPos
46
47 };
48
49 bool operator<(FilePos a, FilePos b){
50 return (a.getPos())->__pos < (b.getPos())->__pos; }
51
52 bool operator==(FilePos a, FilePos b){
53 return (a.getPos())->__pos == (b.getPos())->__pos; }
54
55 vector<FilePos> frameStart;
56 }
57
58 using namespace dumpRead;
59
60
61 DumpReader :: DumpReader( char *in_name ){
62
63 isScanned = false;
64
65 #ifdef IS_MPI
66 if (worldRank == 0) {
67 #endif
68
69 inFile = fopen(in_name, "r");
70 if(inFile == NULL){
71 sprintf(painCave.errMsg,
72 "Cannot open file: %s\n", in_name);
73 painCave.isFatal = 1;
74 simError();
75 }
76
77 strcpy( inName, in_name);
78 #ifdef IS_MPI
79 }
80 strcpy( checkPointMsg, "Dump file opened for reading successfully." );
81 MPIcheckPoint();
82 #endif
83 return;
84 }
85
86 DumpReader :: ~DumpReader( ){
87 #ifdef IS_MPI
88 if (worldRank == 0) {
89 #endif
90 int error;
91 error = fclose( inFile );
92 if( error ){
93 sprintf( painCave.errMsg,
94 "Error closing %s\n", inName );
95 simError();
96 }
97 #ifdef IS_MPI
98 }
99 strcpy( checkPointMsg, "Dump file closed successfully." );
100 MPIcheckPoint();
101 #endif
102
103 return;
104 }
105
106 int DumpReader::getNframes( void ){
107
108 if( !isScanned ) scanFile();
109 return nFrames;
110 }
111
112 void DumpReader::scanFile( void ){
113
114 int vectorSize;
115 int i, j, k;
116 int lineNum = 0;
117 char readBuffer[2000];
118 char* foo;
119 fpos_t *currPos;
120
121 #ifdef IS_MPI
122 if( worldRank == 0 ){
123 #endif // is_mpi
124
125 rewind( inFile );
126
127 currPos = new fpos_t;
128 fgetpos( inFile, currPos );
129 fgets( readBuffer, sizeof( readBuffer ), in_file );
130 lineNum++;
131 if( feof( in_file ) ){
132 sprintf( painCave.errMsg,
133 "File \"%s\" ended unexpectedly at line %d\n",
134 inName,
135 lineNum );
136 painCave.isFatal = 1;
137 simError();
138 }
139
140 while( !feof( in_file ) ){
141
142 frameStart.push_back( FilePos( currPos ) );
143 i = atoi(readBuffer);
144
145 fgets( readBuffer, sizeof( readBuffer ), in_file );
146 lineNum++;
147 if( feof( in_file ) ){
148 sprintf( painCave.errMsg,
149 "File \"%s\" ended unexpectedly at line %d\n",
150 inName,
151 lineNum );
152 painCave.isFatal = 1;
153 simError();
154 }
155
156 if(outTime){
157 foo = strtok( readBuffer, " ,;\t" );
158 time = atof( foo );
159 }
160
161 for(j=0; j<i; j++){
162
163 fgets( readBuffer, sizeof( readBuffer ), in_file );
164 lineNum++;
165 if( feof( in_file ) ){
166 sprintf( painCave.errMsg,
167 "File \"%s\" ended unexpectedly at line %d,"
168 " with atom %d\n",
169 inName,
170 lineNum,
171 j );
172 painCave.isFatal = 1;
173 simError();
174 }
175
176 }
177
178 currPos = new fpos_t;
179 fgetpos( inFile, currPos );
180 fgets( readBuffer, sizeof( readBuffer ), in_file );
181 lineNum++;
182 }
183
184 delete currPos;
185 vectorSize = frameStart.size();
186 rewind( inFile );
187
188 #ifdef IS_MPI
189 }
190 strcpy( checkPointMsg, "Successfully scanned DumpFile\n" );
191 MPIcheckPoint();
192 #endif // is_mpi
193
194 nFrames = vectorSize;
195 isScanned = true;
196 }
197
198
199
200 void DumpReader :: getFrame( SimInfo* the_simnfo, int whichFrame){
201
202 int i, j, done, which_node, which_atom; // loop counter
203
204 const int BUFFERSIZE = 2000; // size of the read buffer
205 int n_atoms; // the number of atoms
206 char read_buffer[BUFFERSIZE]; //the line buffer for reading
207 #ifdef IS_MPI
208 char send_buffer[BUFFERSIZE];
209 #endif
210
211 char *eof_test; // ptr to see when we reach the end of the file
212 char *parseErr;
213 int procIndex;
214 double boxMat[9];
215 double theBoxMat3[3][3];
216
217 simnfo = the_simnfo;
218
219 #ifndef IS_MPI
220 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
221 if( eof_test == NULL ){
222 sprintf( painCave.errMsg,
223 "DumpReader error: error reading 1st line of \"%s\"\n",
224 inName );
225 painCave.isFatal = 1;
226 simError();
227 }
228
229 n_atoms = atoi( read_buffer );
230
231 Atom **atoms = simnfo->atoms;
232 DirectionalAtom* dAtom;
233
234 if( n_atoms != simnfo->n_atoms ){
235 sprintf( painCave.errMsg,
236 "DumpReader error. %s n_atoms, %d, "
237 "does not match the BASS file's n_atoms, %d.\n",
238 inName, n_atoms, simnfo->n_atoms );
239 painCave.isFatal = 1;
240 simError();
241 }
242
243 //read the box mat from the comment line
244
245 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
246 if(eof_test == NULL){
247 sprintf( painCave.errMsg,
248 "error in reading commment in %s\n", inName);
249 painCave.isFatal = 1;
250 simError();
251 }
252
253 parseErr = parseCommentLine( read_buffer, time, boxMat );
254 if( parseErr != NULL ){
255 strcpy( painCave.errMsg, parseErr );
256 painCave.isFatal = 1;
257 simError();
258 }
259
260 simnfo->setTime( time );
261
262 for(i=0;i<3;i++)
263 for(j=0;j<3;j++) theBoxMat3[i][j] = boxMat[3*j+i];
264
265
266 simnfo->setBoxM( theBoxMat3 );
267
268
269 for( i=0; i < n_atoms; i++){
270
271 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
272 if(eof_test == NULL){
273 sprintf(painCave.errMsg,
274 "error in reading file %s\n"
275 "natoms = %d; index = %d\n"
276 "error reading the line from the file.\n",
277 inName, n_atoms, i );
278 painCave.isFatal = 1;
279 simError();
280 }
281
282
283 parseErr = parseDumpLine( read_buffer, i );
284 if( parseErr != NULL ){
285 strcpy( painCave.errMsg, parseErr );
286 painCave.isFatal = 1;
287 simError();
288 }
289 }
290
291
292 // MPI Section of code..........
293 #else //IS_MPI
294
295 // first thing first, suspend fatalities.
296 painCave.isEventLoop = 1;
297
298 int myStatus; // 1 = wakeup & success; 0 = error; -1 = AllDone
299 int haveError;
300
301 MPI_Status istatus;
302 int *AtomToProcMap = mpiSim->getAtomToProcMap();
303
304
305 haveError = 0;
306 if (worldRank == 0) {
307
308 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
309 if( eof_test == NULL ){
310 sprintf( painCave.errMsg,
311 "Error reading 1st line of %d \n ",inName);
312 haveError = 1;
313 simError();
314 }
315
316 n_atoms = atoi( read_buffer );
317
318 Atom **atoms = simnfo->atoms;
319 DirectionalAtom* dAtom;
320
321 // Check to see that the number of atoms in the intial configuration file is the
322 // same as declared in simBass.
323
324 if( n_atoms != mpiSim->getTotAtoms() ){
325 sprintf( painCave.errMsg,
326 "Initialize from File error. %s n_atoms, %d, "
327 "does not match the BASS file's n_atoms, %d.\n",
328 inName, n_atoms, simnfo->n_atoms );
329 haveError= 1;
330 simError();
331 }
332
333 //read the time and boxMat from the comment line
334
335 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
336 if(eof_test == NULL){
337 sprintf( painCave.errMsg,
338 "error in reading commment in %s\n", inName);
339 haveError = 1;
340 simError();
341 }
342
343 parseErr = parseCommentLine( read_buffer, time, boxMat );
344 if( parseErr != NULL ){
345 strcpy( painCave.errMsg, parseErr );
346 haveError = 1;
347 simError();
348 }
349
350 MPI_Bcast(time, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD );
351
352 MPI_Bcast(boxMat, 9, MPI_DOUBLE, 0, MPI_COMM_WORLD );
353
354 if(haveError) nodeZeroError();
355
356 for (i=0 ; i < mpiSim->getTotAtoms(); i++) {
357
358 eof_test = fgets(read_buffer, sizeof(read_buffer), inFile);
359 if(eof_test == NULL){
360 sprintf(painCave.errMsg,
361 "error in reading file %s\n"
362 "natoms = %d; index = %d\n"
363 "error reading the line from the file.\n",
364 inName, n_atoms, i );
365 haveError= 1;
366 simError();
367 }
368
369 if(haveError) nodeZeroError();
370
371 // Get the Node number which wants this atom:
372 which_node = AtomToProcMap[i];
373 if (which_node == 0) {
374 parseErr = parseDumpLine( read_buffer, i );
375 if( parseErr != NULL ){
376 strcpy( painCave.errMsg, parseErr );
377 haveError = 1;
378 simError();
379 }
380 if(haveError) nodeZeroError();
381 }
382
383 else {
384
385 myStatus = 1;
386 MPI_Send(&myStatus, 1, MPI_INT, which_node,
387 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
388 MPI_Send(read_buffer, BUFFERSIZE, MPI_CHAR, which_node,
389 TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD);
390 MPI_Send(&i, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
391 MPI_COMM_WORLD);
392 MPI_Recv(&myStatus, 1, MPI_INT, which_node, TAKE_THIS_TAG_INT,
393 MPI_COMM_WORLD, &istatus);
394
395 if(!myStatus) nodeZeroError();
396 }
397 }
398 myStatus = -1;
399 for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
400 MPI_Send( &myStatus, 1, MPI_INT, j,
401 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
402 }
403
404 } else {
405
406 MPI_Bcast(time, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
407 MPI_Bcast(boxMat, 9, MPI_DOUBLE, 0, MPI_COMM_WORLD);
408
409 done = 0;
410 while (!done) {
411
412 MPI_Recv(&myStatus, 1, MPI_INT, 0,
413 TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
414
415 if(!myStatus) anonymousNodeDie();
416
417 if(myStatus < 0) break;
418
419 MPI_Recv(read_buffer, BUFFERSIZE, MPI_CHAR, 0,
420 TAKE_THIS_TAG_CHAR, MPI_COMM_WORLD, &istatus);
421 MPI_Recv(&which_atom, 1, MPI_INT, 0,
422 TAKE_THIS_TAG_INT, MPI_COMM_WORLD, &istatus);
423
424 myStatus = 1;
425 parseErr = parseDumpLine( read_buffer, which_atom );
426 if( parseErr != NULL ){
427 strcpy( painCave.errMsg, parseErr );
428 myStatus = 0;;
429 simError();
430 }
431
432 MPI_Send( &myStatus, 1, MPI_INT, 0,
433 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
434
435 }
436 }
437
438 // last thing last, enable fatalities.
439 painCave.isEventLoop = 0;
440
441 simnfo->setTime( time );
442
443 for(i=0;i<3;i++)
444 for(j=0;j<3;j++) theBoxMat3[i][j] = boxMat[3*j+i];
445
446 simnfo->setBoxM( theBoxMat3 );
447
448
449 #endif
450 }
451
452 char* DumpReader::parseDumpLine(char* readLine, int globalIndex){
453
454 char *foo; // the pointer to the current string token
455
456 double rx, ry, rz; // position place holders
457 double vx, vy, vz; // velocity placeholders
458 double q[4]; // the quaternions
459 double jx, jy, jz; // angular velocity placeholders;
460 double qSqr, qLength; // needed to normalize the quaternion vector.
461
462 Atom **atoms = simnfo->atoms;
463 DirectionalAtom* dAtom;
464
465 int j, n_atoms, atomIndex;
466
467 #ifdef IS_MPI
468 n_atoms = mpiSim->getTotAtoms();
469 atomIndex=-1;
470 for (j=0; j < mpiSim->getMyNlocal(); j++) {
471 if (atoms[j]->getGlobalIndex() == globalIndex) atomIndex = j;
472 }
473 if (atomIndex == -1) {
474 sprintf( painCave.errMsg,
475 "Initialize from file error. Atom at index %d "
476 "in file %s does not exist on processor %d .\n",
477 globalIndex, inName, mpiSim->getMyNode() );
478 return strdup( painCave.errMsg );
479 }
480 #else
481 n_atoms = simnfo->n_atoms;
482 atomIndex = globalIndex;
483 #endif // is_mpi
484
485 // set the string tokenizer
486
487 foo = strtok(readLine, " ,;\t");
488
489 // check the atom name to the current atom
490
491 if( strcmp( foo, atoms[atomIndex]->getType() ) ){
492 sprintf( painCave.errMsg,
493 "Initialize from file error. Atom %s at index %d "
494 "in file %s does not"
495 " match the BASS atom %s.\n",
496 foo, atomIndex, inName, atoms[atomIndex]->getType() );
497 return strdup( painCave.errMsg );
498 }
499
500 // get the positions
501
502 foo = strtok(NULL, " ,;\t");
503 if(foo == NULL){
504 sprintf( painCave.errMsg,
505 "error in reading postition x from %s\n"
506 "natoms = %d, index = %d\n",
507 inName, n_atoms, atomIndex );
508 return strdup( painCave.errMsg );
509 }
510 rx = atof( foo );
511
512 foo = strtok(NULL, " ,;\t");
513 if(foo == NULL){
514 sprintf( painCave.errMsg,
515 "error in reading postition y from %s\n"
516 "natoms = %d, index = %d\n",
517 inName, n_atoms, atomIndex );
518 return strdup( painCave.errMsg );
519 }
520 ry = atof( foo );
521
522 foo = strtok(NULL, " ,;\t");
523 if(foo == NULL){
524 sprintf( painCave.errMsg,
525 "error in reading postition z from %s\n"
526 "natoms = %d, index = %d\n",
527 inName, n_atoms, atomIndex );
528 return strdup( painCave.errMsg );
529 }
530 rz = atof( foo );
531
532
533 // get the velocities
534
535 foo = strtok(NULL, " ,;\t");
536 if(foo == NULL){
537 sprintf( painCave.errMsg,
538 "error in reading velocity x from %s\n"
539 "natoms = %d, index = %d\n",
540 inName, n_atoms, atomIndex );
541 return strdup( painCave.errMsg );
542 }
543 vx = atof( foo );
544
545 foo = strtok(NULL, " ,;\t");
546 if(foo == NULL){
547 sprintf( painCave.errMsg,
548 "error in reading velocity y from %s\n"
549 "natoms = %d, index = %d\n",
550 inName, n_atoms, atomIndex );
551 return strdup( painCave.errMsg );
552 }
553 vy = atof( foo );
554
555 foo = strtok(NULL, " ,;\t");
556 if(foo == NULL){
557 sprintf( painCave.errMsg,
558 "error in reading velocity z from %s\n"
559 "natoms = %d, index = %d\n",
560 inName, n_atoms, atomIndex );
561 return strdup( painCave.errMsg );
562 }
563 vz = atof( foo );
564
565
566 // get the quaternions
567
568 if( atoms[atomIndex]->isDirectional() ){
569
570 foo = strtok(NULL, " ,;\t");
571 if(foo == NULL){
572 sprintf(painCave.errMsg,
573 "error in reading quaternion 0 from %s\n"
574 "natoms = %d, index = %d\n",
575 inName, n_atoms, atomIndex );
576 return strdup( painCave.errMsg );
577 }
578 q[0] = atof( foo );
579
580 foo = strtok(NULL, " ,;\t");
581 if(foo == NULL){
582 sprintf( painCave.errMsg,
583 "error in reading quaternion 1 from %s\n"
584 "natoms = %d, index = %d\n",
585 inName, n_atoms, atomIndex );
586 return strdup( painCave.errMsg );
587 }
588 q[1] = atof( foo );
589
590 foo = strtok(NULL, " ,;\t");
591 if(foo == NULL){
592 sprintf( painCave.errMsg,
593 "error in reading quaternion 2 from %s\n"
594 "natoms = %d, index = %d\n",
595 inName, n_atoms, atomIndex );
596 return strdup( painCave.errMsg );
597 }
598 q[2] = atof( foo );
599
600 foo = strtok(NULL, " ,;\t");
601 if(foo == NULL){
602 sprintf( painCave.errMsg,
603 "error in reading quaternion 3 from %s\n"
604 "natoms = %d, index = %d\n",
605 inName, n_atoms, atomIndex );
606 return strdup( painCave.errMsg );
607 }
608 q[3] = atof( foo );
609
610 // get the angular velocities
611
612 foo = strtok(NULL, " ,;\t");
613 if(foo == NULL){
614 sprintf( painCave.errMsg,
615 "error in reading angular momentum jx from %s\n"
616 "natoms = %d, index = %d\n",
617 inName, n_atoms, atomIndex );
618 return strdup( painCave.errMsg );
619 }
620 jx = atof( foo );
621
622 foo = strtok(NULL, " ,;\t");
623 if(foo == NULL){
624 sprintf( painCave.errMsg,
625 "error in reading angular momentum jy from %s\n"
626 "natoms = %d, index = %d\n",
627 inName, n_atoms, atomIndex );
628 return strdup( painCave.errMsg );
629 }
630 jy = atof(foo );
631
632 foo = strtok(NULL, " ,;\t");
633 if(foo == NULL){
634 sprintf( painCave.errMsg,
635 "error in reading angular momentum jz from %s\n"
636 "natoms = %d, index = %d\n",
637 inName, n_atoms, atomIndex );
638 return strdup( painCave.errMsg );
639 }
640 jz = atof( foo );
641
642 dAtom = ( DirectionalAtom* )atoms[atomIndex];
643
644 // check that the quaternion vector is normalized
645
646 qSqr = (q[0] * q[0]) + (q[1] * q[1]) + (q[2] * q[2]) + (q[3] * q[3]);
647
648 qLength = sqrt( qSqr );
649 q[0] = q[0] / qLength;
650 q[1] = q[1] / qLength;
651 q[2] = q[2] / qLength;
652 q[3] = q[3] / qLength;
653
654 dAtom->setQ( q );
655
656 // add the angular velocities
657
658 dAtom->setJx( jx );
659 dAtom->setJy( jy );
660 dAtom->setJz( jz );
661 }
662
663 // add the positions and velocities to the atom
664
665 atoms[atomIndex]->setX( rx );
666 atoms[atomIndex]->setY( ry );
667 atoms[atomIndex]->setZ( rz );
668
669 atoms[atomIndex]->set_vx( vx );
670 atoms[atomIndex]->set_vy( vy );
671 atoms[atomIndex]->set_vz( vz );
672
673 return NULL;
674 }
675
676
677 char* DumpReader::parseCommentLine(char* readLine, double time,
678 double boxMat[9]){
679
680 char *foo; // the pointer to the current string token
681 int j;
682
683 // set the string tokenizer
684
685 foo = strtok(readLine, " ,;\t");
686 if(foo == NULL){
687 sprintf( painCave.errMsg,
688 "error in reading time from %s\n",
689 inName );
690 return strdup( painCave.errMsg );
691 }
692 time = atof( foo );
693
694 // get the Hx vector
695
696 foo = strtok(NULL, " ,;\t");
697 if(foo == NULL){
698 sprintf( painCave.errMsg,
699 "error in reading Hx[0] from %s\n",
700 inName );
701 return strdup( painCave.errMsg );
702 }
703 boxMat[0] = atof( foo );
704
705 foo = strtok(NULL, " ,;\t");
706 if(foo == NULL){
707 sprintf( painCave.errMsg,
708 "error in reading Hx[1] from %s\n",
709 inName );
710 return strdup( painCave.errMsg );
711 }
712 boxMat[1] = atof( foo );
713
714 foo = strtok(NULL, " ,;\t");
715 if(foo == NULL){
716 sprintf( painCave.errMsg,
717 "error in reading Hx[2] from %s\n",
718 inName );
719 return strdup( painCave.errMsg );
720 }
721 boxMat[2] = atof( foo );
722
723 // get the Hy vector
724
725 foo = strtok(NULL, " ,;\t");
726 if(foo == NULL){
727 sprintf( painCave.errMsg,
728 "error in reading Hy[0] from %s\n",
729 inName );
730 return strdup( painCave.errMsg );
731 }
732 boxMat[3] = atof( foo );
733
734 foo = strtok(NULL, " ,;\t");
735 if(foo == NULL){
736 sprintf( painCave.errMsg,
737 "error in reading Hy[1] from %s\n",
738 inName );
739 return strdup( painCave.errMsg );
740 }
741 boxMat[4] = atof( foo );
742
743 foo = strtok(NULL, " ,;\t");
744 if(foo == NULL){
745 sprintf( painCave.errMsg,
746 "error in reading Hy[2] from %s\n",
747 inName );
748 return strdup( painCave.errMsg );
749 }
750 boxMat[5] = atof( foo );
751
752 // get the Hz vector
753
754 foo = strtok(NULL, " ,;\t");
755 if(foo == NULL){
756 sprintf( painCave.errMsg,
757 "error in reading Hz[0] from %s\n",
758 inName );
759 return strdup( painCave.errMsg );
760 }
761 boxMat[6] = atof( foo );
762
763 foo = strtok(NULL, " ,;\t");
764 if(foo == NULL){
765 sprintf( painCave.errMsg,
766 "error in reading Hz[1] from %s\n",
767 inName );
768 return strdup( painCave.errMsg );
769 }
770 boxMat[7] = atof( foo );
771
772 foo = strtok(NULL, " ,;\t");
773 if(foo == NULL){
774 sprintf( painCave.errMsg,
775 "error in reading Hz[2] from %s\n",
776 inName );
777 return strdup( painCave.errMsg );
778 }
779 boxMat[8] = atof( foo );
780
781 return NULL;
782 }
783
784
785 #ifdef IS_MPI
786
787 // a couple of functions to let us escape the read loop
788
789 void initFile::nodeZeroError( void ){
790 int j, myStatus;
791
792 myStatus = 0;
793 for (j = 0; j < mpiSim->getNumberProcessors(); j++) {
794 MPI_Send( &myStatus, 1, MPI_INT, j,
795 TAKE_THIS_TAG_INT, MPI_COMM_WORLD);
796 }
797
798
799 MPI_Finalize();
800 exit (0);
801
802 }
803
804 void initFile::anonymousNodeDie( void ){
805
806 MPI_Finalize();
807 exit (0);
808 }
809
810 #endif //is_mpi