xref: /phasta/phastaIO/phastaIO.cc (revision f42e0444da478b8804248904e7be4eb4355d183d)
159599516SKenneth E. Jansen /*
259599516SKenneth E. Jansen  *
359599516SKenneth E. Jansen  * This is the SyncIO library that uses MPI-IO collective functions to
459599516SKenneth E. Jansen  * implement a flexible I/O checkpoint solution for a large number of
559599516SKenneth E. Jansen  * processors.
659599516SKenneth E. Jansen  *
759599516SKenneth E. Jansen  * Previous developer: Ning Liu         (liun2@cs.rpi.edu)
859599516SKenneth E. Jansen  *                     Jing Fu          (fuj@cs.rpi.edu)
959599516SKenneth E. Jansen  * Current developers: Michel Rasquin   (Michel.Rasquin@colorado.edu),
1059599516SKenneth E. Jansen  *                     Ben Matthews     (benjamin.a.matthews@colorado.edu)
1159599516SKenneth E. Jansen  *
1259599516SKenneth E. Jansen  */
1359599516SKenneth E. Jansen 
1459599516SKenneth E. Jansen #include <map>
1559599516SKenneth E. Jansen #include <vector>
1659599516SKenneth E. Jansen #include <string>
17a746c838SCameron Smith #include <cstdarg>
1859599516SKenneth E. Jansen #include <string.h>
1959599516SKenneth E. Jansen #include <ctype.h>
2059599516SKenneth E. Jansen #include <stdlib.h>
2159599516SKenneth E. Jansen #include <stdio.h>
2259599516SKenneth E. Jansen #include <math.h>
23d3337298SCameron Smith #include <sstream>
2459599516SKenneth E. Jansen #include "phastaIO.h"
2559599516SKenneth E. Jansen #include "mpi.h"
2659599516SKenneth E. Jansen #include "phiotmrc.h"
27*f42e0444SCameron Smith #include "phiostats.h"
2897a07b0aSCameron Smith #include <assert.h>
2959599516SKenneth E. Jansen 
3059599516SKenneth E. Jansen #define VERSION_INFO_HEADER_SIZE 8192
3159599516SKenneth E. Jansen #define DB_HEADER_SIZE 1024
3259599516SKenneth E. Jansen #define ONE_MEGABYTE 1048576
3359599516SKenneth E. Jansen #define TWO_MEGABYTE 2097152
3459599516SKenneth E. Jansen #define ENDIAN_TEST_NUMBER 12180 // Troy's Zip Code!!
3559599516SKenneth E. Jansen #define MAX_PHASTA_FILES 64
3659599516SKenneth E. Jansen #define MAX_PHASTA_FILE_NAME_LENGTH 1024
3759599516SKenneth E. Jansen #define MAX_FIELDS_NUMBER ((VERSION_INFO_HEADER_SIZE/MAX_FIELDS_NAME_LENGTH)-4) // The meta data include - MPI_IO_Tag, nFields, nFields*names of the fields, nppf
3859599516SKenneth E. Jansen // -3 for MPI_IO_Tag, nFields and nppf, -4 for extra security (former nFiles)
3959599516SKenneth E. Jansen #define MAX_FIELDS_NAME_LENGTH 128
4059599516SKenneth E. Jansen #define DefaultMHSize (4*ONE_MEGABYTE)
4159599516SKenneth E. Jansen //#define DefaultMHSize (8350) //For test
4259599516SKenneth E. Jansen #define LATEST_WRITE_VERSION 1
4359599516SKenneth E. Jansen #define inv1024sq 953.674316406e-9 // = 1/1024/1024
4459599516SKenneth E. Jansen int MasterHeaderSize = -1;
4559599516SKenneth E. Jansen 
4659599516SKenneth E. Jansen bool PRINT_PERF = false; // default print no perf results
4759599516SKenneth E. Jansen int irank = -1; // global rank, should never be manually manipulated
4859599516SKenneth E. Jansen int mysize = -1;
4959599516SKenneth E. Jansen 
5059599516SKenneth E. Jansen // Static variables are bad but used here to store the subcommunicator and associated variables
5159599516SKenneth E. Jansen // Prevent MPI_Comm_split to be called more than once, especially on BGQ with the V1R2M1 driver (leak detected in MPI_Comm_split - IBM working on it)
5259599516SKenneth E. Jansen static int s_assign_local_comm = 0;
5359599516SKenneth E. Jansen static MPI_Comm s_local_comm;
5459599516SKenneth E. Jansen static int s_local_size = -1;
5559599516SKenneth E. Jansen static int s_local_rank = -1;
5659599516SKenneth E. Jansen 
5759599516SKenneth E. Jansen // the following defines are for debug printf
5859599516SKenneth E. Jansen #define PHASTAIO_DEBUG 0 //default to not print any debugging info
5959599516SKenneth E. Jansen 
60fe03b87fSCameron Smith void phprintf(const char* fmt, ...) {
61be3da47bSCameron Smith   (void)fmt;
6259599516SKenneth E. Jansen #if PHASTAIO_DEBUG
63fe03b87fSCameron Smith   char format[1024];
64fe03b87fSCameron Smith   snprintf(format, sizeof(format), "phastaIO debug: %s", fmt);
65fe03b87fSCameron Smith   va_list ap;
66fe03b87fSCameron Smith   va_start(ap,fmt);
67a746c838SCameron Smith   vprintf(format,ap);
68fe03b87fSCameron Smith   va_end(ap);
6959599516SKenneth E. Jansen #endif
70fe03b87fSCameron Smith }
71fe03b87fSCameron Smith 
72fe03b87fSCameron Smith void phprintf_0(const char* fmt, ...) {
73be3da47bSCameron Smith   (void)fmt;
74fe03b87fSCameron Smith #if PHASTAIO_DEBUG
75fe03b87fSCameron Smith   int rank = 0;
76fe03b87fSCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
77fe03b87fSCameron Smith   if(rank == 0){
78fe03b87fSCameron Smith     char format[1024];
79fe03b87fSCameron Smith     snprintf(format, sizeof(format), "phastaIO debug: irank=0 %s", fmt);
80fe03b87fSCameron Smith     va_list ap;
81a746c838SCameron Smith     va_start(ap,fmt);
82fe03b87fSCameron Smith     vprintf(format, ap);
83fe03b87fSCameron Smith     va_end(ap);
84fe03b87fSCameron Smith   }
85fe03b87fSCameron Smith #endif
86fe03b87fSCameron Smith }
8759599516SKenneth E. Jansen 
8859599516SKenneth E. Jansen enum PhastaIO_Errors
8959599516SKenneth E. Jansen {
9059599516SKenneth E. Jansen   MAX_PHASTA_FILES_EXCEEDED = -1,
9159599516SKenneth E. Jansen   UNABLE_TO_OPEN_FILE = -2,
9259599516SKenneth E. Jansen   NOT_A_MPI_FILE = -3,
9359599516SKenneth E. Jansen   GPID_EXCEEDED = -4,
9410d56689SCameron Smith   DATA_TYPE_ILLEGAL = -5
9559599516SKenneth E. Jansen };
9659599516SKenneth E. Jansen 
9759599516SKenneth E. Jansen using namespace std;
9859599516SKenneth E. Jansen 
9959599516SKenneth E. Jansen namespace{
10059599516SKenneth E. Jansen 
101961a4ff6SCameron Smith   map<int, std::string> LastHeaderKey;
10259599516SKenneth E. Jansen   vector< FILE* > fileArray;
10359599516SKenneth E. Jansen   vector< bool > byte_order;
10459599516SKenneth E. Jansen   vector< int > header_type;
10559599516SKenneth E. Jansen   int DataSize=0;
10659599516SKenneth E. Jansen   bool LastHeaderNotFound = false;
10759599516SKenneth E. Jansen   bool Wrong_Endian = false ;
10859599516SKenneth E. Jansen   bool Strict_Error = false ;
10959599516SKenneth E. Jansen   bool binary_format = true;
11059599516SKenneth E. Jansen 
11159599516SKenneth E. Jansen   /***********************************************************************/
11259599516SKenneth E. Jansen   /***************** NEW PHASTA IO CODE STARTS HERE **********************/
11359599516SKenneth E. Jansen   /***********************************************************************/
11459599516SKenneth E. Jansen 
11559599516SKenneth E. Jansen   typedef struct
11659599516SKenneth E. Jansen   {
11759599516SKenneth E. Jansen     char filename[MAX_PHASTA_FILE_NAME_LENGTH];   /* defafults to 1024 */
1182dd307a1SCameron Smith     unsigned long my_offset;
1192dd307a1SCameron Smith     unsigned long next_start_address;
1202dd307a1SCameron Smith     unsigned long **my_offset_table;
1212dd307a1SCameron Smith     unsigned long **my_read_table;
12259599516SKenneth E. Jansen 
12359599516SKenneth E. Jansen     double * double_chunk;
12459599516SKenneth E. Jansen     double * read_double_chunk;
12559599516SKenneth E. Jansen 
12659599516SKenneth E. Jansen     int field_count;
12759599516SKenneth E. Jansen     int part_count;
12859599516SKenneth E. Jansen     int read_field_count;
12959599516SKenneth E. Jansen     int read_part_count;
13059599516SKenneth E. Jansen     int GPid;
13159599516SKenneth E. Jansen     int start_id;
13259599516SKenneth E. Jansen 
13359599516SKenneth E. Jansen     int mhsize;
13459599516SKenneth E. Jansen 
13559599516SKenneth E. Jansen     int myrank;
13659599516SKenneth E. Jansen     int numprocs;
13759599516SKenneth E. Jansen     int local_myrank;
13859599516SKenneth E. Jansen     int local_numprocs;
13959599516SKenneth E. Jansen 
14059599516SKenneth E. Jansen     int nppp;
14159599516SKenneth E. Jansen     int nPPF;
14259599516SKenneth E. Jansen     int nFiles;
14359599516SKenneth E. Jansen     int nFields;
14459599516SKenneth E. Jansen 
14559599516SKenneth E. Jansen     int * int_chunk;
14659599516SKenneth E. Jansen     int * read_int_chunk;
14759599516SKenneth E. Jansen 
14859599516SKenneth E. Jansen     int Wrong_Endian; /* default to false */
14959599516SKenneth E. Jansen     char * master_header;
15059599516SKenneth E. Jansen     MPI_File file_handle;
15159599516SKenneth E. Jansen     MPI_Comm local_comm;
15259599516SKenneth E. Jansen   } phastaio_file_t;
15359599516SKenneth E. Jansen 
15459599516SKenneth E. Jansen   typedef struct
15559599516SKenneth E. Jansen   {
15659599516SKenneth E. Jansen     int nppf, nfields;
15759599516SKenneth E. Jansen     char * masterHeader;
15859599516SKenneth E. Jansen   }serial_file;
15959599516SKenneth E. Jansen 
16059599516SKenneth E. Jansen   serial_file *SerialFile;
16159599516SKenneth E. Jansen   phastaio_file_t *PhastaIOActiveFiles[MAX_PHASTA_FILES];
16259599516SKenneth E. Jansen   int PhastaIONextActiveIndex = 0; /* indicates next index to allocate */
16359599516SKenneth E. Jansen 
16459599516SKenneth E. Jansen   // the caller has the responsibility to delete the returned string
16559599516SKenneth E. Jansen   // TODO: StringStipper("nbc value? ") returns NULL?
166103be424SCameron Smith   char* StringStripper( const char  istring[] ) {
16759599516SKenneth E. Jansen     int length = strlen( istring );
16859599516SKenneth E. Jansen     char* dest = (char *)malloc( length + 1 );
16959599516SKenneth E. Jansen     strcpy( dest, istring );
17059599516SKenneth E. Jansen     dest[ length ] = '\0';
17159599516SKenneth E. Jansen     if ( char* p = strpbrk( dest, " ") )
17259599516SKenneth E. Jansen       *p = '\0';
17359599516SKenneth E. Jansen     return dest;
17459599516SKenneth E. Jansen   }
17559599516SKenneth E. Jansen 
176103be424SCameron Smith   inline int cscompare( const char teststring[], const char targetstring[] ) {
17759599516SKenneth E. Jansen     char* s1 = const_cast<char*>(teststring);
17859599516SKenneth E. Jansen     char* s2 = const_cast<char*>(targetstring);
17959599516SKenneth E. Jansen 
18059599516SKenneth E. Jansen     while( *s1 == ' ') s1++;
18159599516SKenneth E. Jansen     while( *s2 == ' ') s2++;
18259599516SKenneth E. Jansen     while( ( *s1 )
18359599516SKenneth E. Jansen         && ( *s2 )
18459599516SKenneth E. Jansen         && ( *s2 != '?')
18559599516SKenneth E. Jansen         && ( tolower( *s1 )==tolower( *s2 ) ) ) {
18659599516SKenneth E. Jansen       s1++;
18759599516SKenneth E. Jansen       s2++;
18859599516SKenneth E. Jansen       while( *s1 == ' ') s1++;
18959599516SKenneth E. Jansen       while( *s2 == ' ') s2++;
19059599516SKenneth E. Jansen     }
19159599516SKenneth E. Jansen     if ( !( *s1 ) || ( *s1 == '?') ) return 1;
19259599516SKenneth E. Jansen     else return 0;
19359599516SKenneth E. Jansen   }
19459599516SKenneth E. Jansen 
195103be424SCameron Smith   inline void isBinary( const char iotype[] ) {
19659599516SKenneth E. Jansen     char* fname = StringStripper( iotype );
19759599516SKenneth E. Jansen     if ( cscompare( fname, "binary" ) ) binary_format = true;
19859599516SKenneth E. Jansen     else binary_format = false;
19959599516SKenneth E. Jansen     free (fname);
20059599516SKenneth E. Jansen 
20159599516SKenneth E. Jansen   }
20259599516SKenneth E. Jansen 
203103be424SCameron Smith   inline size_t typeSize( const char typestring[] ) {
20459599516SKenneth E. Jansen     char* ts1 = StringStripper( typestring );
20559599516SKenneth E. Jansen     if ( cscompare( "integer", ts1 ) ) {
20659599516SKenneth E. Jansen       free (ts1);
20759599516SKenneth E. Jansen       return sizeof(int);
20859599516SKenneth E. Jansen     } else if ( cscompare( "double", ts1 ) ) {
20959599516SKenneth E. Jansen       free (ts1);
21059599516SKenneth E. Jansen       return sizeof( double );
21159599516SKenneth E. Jansen     } else {
21259599516SKenneth E. Jansen       free (ts1);
21359599516SKenneth E. Jansen       fprintf(stderr,"unknown type : %s\n",ts1);
21459599516SKenneth E. Jansen       return 0;
21559599516SKenneth E. Jansen     }
21659599516SKenneth E. Jansen   }
21759599516SKenneth E. Jansen 
218103be424SCameron Smith   int readHeader(
219103be424SCameron Smith       FILE*       fileObject,
22059599516SKenneth E. Jansen         const char  phrase[],
22159599516SKenneth E. Jansen         int*        params,
22259599516SKenneth E. Jansen         int         expect ) {
22359599516SKenneth E. Jansen     char* text_header;
22459599516SKenneth E. Jansen     char* token;
22597a07b0aSCameron Smith     char Line[1024] = "\0";
22659599516SKenneth E. Jansen     char junk;
22759599516SKenneth E. Jansen     bool FOUND = false ;
22859599516SKenneth E. Jansen     int real_length;
22959599516SKenneth E. Jansen     int skip_size, integer_value;
23059599516SKenneth E. Jansen     int rewind_count=0;
23159599516SKenneth E. Jansen 
23259599516SKenneth E. Jansen     if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
23359599516SKenneth E. Jansen       rewind( fileObject );
23459599516SKenneth E. Jansen       clearerr( fileObject );
23559599516SKenneth E. Jansen       rewind_count++;
23659599516SKenneth E. Jansen       fgets( Line, 1024, fileObject );
23759599516SKenneth E. Jansen     }
23859599516SKenneth E. Jansen 
23959599516SKenneth E. Jansen     while( !FOUND  && ( rewind_count < 2 ) )  {
24059599516SKenneth E. Jansen       if ( ( Line[0] != '\n' ) && ( real_length = strcspn( Line, "#" )) ) {
24159599516SKenneth E. Jansen         text_header = (char *)malloc( real_length + 1 );
24259599516SKenneth E. Jansen         strncpy( text_header, Line, real_length );
24359599516SKenneth E. Jansen         text_header[ real_length ] =static_cast<char>(NULL);
24459599516SKenneth E. Jansen         token = strtok ( text_header, ":" );
24597a07b0aSCameron Smith         assert(token);
24697a07b0aSCameron Smith         if( cscompare( phrase , token ) ) {
24759599516SKenneth E. Jansen           FOUND = true ;
24859599516SKenneth E. Jansen           token = strtok( NULL, " ,;<>" );
24997a07b0aSCameron Smith           assert(token);
25059599516SKenneth E. Jansen           skip_size = atoi( token );
25159599516SKenneth E. Jansen           int i;
25259599516SKenneth E. Jansen           for( i=0; i < expect && ( token = strtok( NULL," ,;<>") ); i++) {
25397a07b0aSCameron Smith             assert(token);
25459599516SKenneth E. Jansen             params[i] = atoi( token );
25559599516SKenneth E. Jansen           }
25659599516SKenneth E. Jansen           if ( i < expect ) {
257103be424SCameron Smith             fprintf(stderr,
258103be424SCameron Smith                 "Aloha Expected # of ints not found for: %s\n",phrase );
25959599516SKenneth E. Jansen           }
26059599516SKenneth E. Jansen         } else if ( cscompare(token,"byteorder magic number") ) {
26159599516SKenneth E. Jansen           if ( binary_format ) {
26259599516SKenneth E. Jansen             fread((void*)&integer_value,sizeof(int),1,fileObject);
26359599516SKenneth E. Jansen             fread( &junk, sizeof(char), 1 , fileObject );
26459599516SKenneth E. Jansen             if ( 362436 != integer_value ) Wrong_Endian = true;
26559599516SKenneth E. Jansen           } else{
26659599516SKenneth E. Jansen             fscanf(fileObject, "%d\n", &integer_value );
26759599516SKenneth E. Jansen           }
26859599516SKenneth E. Jansen         } else {
26959599516SKenneth E. Jansen           /* some other header, so just skip over */
27059599516SKenneth E. Jansen           token = strtok( NULL, " ,;<>" );
27197a07b0aSCameron Smith           assert(token);
27259599516SKenneth E. Jansen           skip_size = atoi( token );
27359599516SKenneth E. Jansen           if ( binary_format)
27459599516SKenneth E. Jansen             fseek( fileObject, skip_size, SEEK_CUR );
27559599516SKenneth E. Jansen           else
27659599516SKenneth E. Jansen             for( int gama=0; gama < skip_size; gama++ )
27759599516SKenneth E. Jansen               fgets( Line, 1024, fileObject );
27859599516SKenneth E. Jansen         }
27959599516SKenneth E. Jansen         free (text_header);
28059599516SKenneth E. Jansen       } // end of if before while loop
28159599516SKenneth E. Jansen 
28259599516SKenneth E. Jansen       if ( !FOUND )
28359599516SKenneth E. Jansen         if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
28459599516SKenneth E. Jansen           rewind( fileObject );
28559599516SKenneth E. Jansen           clearerr( fileObject );
28659599516SKenneth E. Jansen           rewind_count++;
28759599516SKenneth E. Jansen           fgets( Line, 1024, fileObject );
28859599516SKenneth E. Jansen         }
28959599516SKenneth E. Jansen     }
29059599516SKenneth E. Jansen 
29159599516SKenneth E. Jansen     if ( !FOUND ) {
29259599516SKenneth E. Jansen       //fprintf(stderr, "Error: Could not find: %s\n", phrase);
29359599516SKenneth E. Jansen       if(irank == 0) printf("WARNING: Could not find: %s\n", phrase);
29459599516SKenneth E. Jansen       return 1;
29559599516SKenneth E. Jansen     }
29659599516SKenneth E. Jansen     return 0;
29759599516SKenneth E. Jansen   }
29859599516SKenneth E. Jansen 
29959599516SKenneth E. Jansen } // end unnamed namespace
30059599516SKenneth E. Jansen 
30159599516SKenneth E. Jansen 
30259599516SKenneth E. Jansen // begin of publicly visible functions
30359599516SKenneth E. Jansen 
30459599516SKenneth E. Jansen /**
30559599516SKenneth E. Jansen  * This function takes a long long pointer and assign (start) phiotmrc value to it
30659599516SKenneth E. Jansen  */
30759599516SKenneth E. Jansen void startTimer(double* start) {
30859599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
30959599516SKenneth E. Jansen   MPI_Barrier(MPI_COMM_WORLD);
31059599516SKenneth E. Jansen   *start =  phiotmrc();
31159599516SKenneth E. Jansen }
31259599516SKenneth E. Jansen 
31359599516SKenneth E. Jansen /**
31459599516SKenneth E. Jansen  * This function takes a long long pointer and assign (end) phiotmrc value to it
31559599516SKenneth E. Jansen  */
31659599516SKenneth E. Jansen void endTimer(double* end) {
31759599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
31859599516SKenneth E. Jansen   *end = phiotmrc();
31959599516SKenneth E. Jansen   MPI_Barrier(MPI_COMM_WORLD);
32059599516SKenneth E. Jansen }
32159599516SKenneth E. Jansen 
32259599516SKenneth E. Jansen /**
32359599516SKenneth E. Jansen  * choose to print some performance results (or not) according to
32459599516SKenneth E. Jansen  * the PRINT_PERF macro
32559599516SKenneth E. Jansen  */
32659599516SKenneth E. Jansen void printPerf(
32759599516SKenneth E. Jansen     const char* func_name,
32859599516SKenneth E. Jansen     double start,
32959599516SKenneth E. Jansen     double end,
3302dd307a1SCameron Smith     unsigned long datasize,
33159599516SKenneth E. Jansen     int printdatainfo,
33259599516SKenneth E. Jansen     const char* extra_msg) {
33359599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
3342dd307a1SCameron Smith   unsigned long data_size = datasize;
33559599516SKenneth E. Jansen   double time = end - start;
3362dd307a1SCameron Smith   unsigned long isizemin,isizemax,isizetot;
33759599516SKenneth E. Jansen   double sizemin,sizemax,sizeavg,sizetot,rate;
33859599516SKenneth E. Jansen   double tmin, tmax, tavg, ttot;
33959599516SKenneth E. Jansen 
34059599516SKenneth E. Jansen   MPI_Allreduce(&time, &tmin,1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
34159599516SKenneth E. Jansen   MPI_Allreduce(&time, &tmax,1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
34259599516SKenneth E. Jansen   MPI_Allreduce(&time, &ttot,1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
34359599516SKenneth E. Jansen   tavg = ttot/mysize;
34459599516SKenneth E. Jansen 
34559599516SKenneth E. Jansen   if(irank == 0) {
34659599516SKenneth E. Jansen     if ( PhastaIONextActiveIndex == 0 ) printf("** 1PFPP ");
34759599516SKenneth E. Jansen     else  printf("** syncIO ");
34859599516SKenneth E. Jansen     printf("%s(): Tmax = %f sec, Tmin = %f sec, Tavg = %f sec", func_name, tmax, tmin, tavg);
34959599516SKenneth E. Jansen   }
35059599516SKenneth E. Jansen 
35159599516SKenneth E. Jansen   if(printdatainfo == 1) { // if printdatainfo ==1, compute I/O rate and block size
35259599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizemin,1,MPI_LONG_LONG_INT,MPI_MIN,MPI_COMM_WORLD);
35359599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizemax,1,MPI_LONG_LONG_INT,MPI_MAX,MPI_COMM_WORLD);
35459599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizetot,1,MPI_LONG_LONG_INT,MPI_SUM,MPI_COMM_WORLD);
35559599516SKenneth E. Jansen 
35659599516SKenneth E. Jansen     sizemin=(double)(isizemin*inv1024sq);
35759599516SKenneth E. Jansen     sizemax=(double)(isizemax*inv1024sq);
35859599516SKenneth E. Jansen     sizetot=(double)(isizetot*inv1024sq);
35959599516SKenneth E. Jansen     sizeavg=(double)(1.0*sizetot/mysize);
36059599516SKenneth E. Jansen     rate=(double)(1.0*sizetot/tmax);
36159599516SKenneth E. Jansen 
36259599516SKenneth E. Jansen     if( irank == 0) {
363103be424SCameron Smith       printf(", Rate = %f MB/s [%s] \n \t\t\t"
364103be424SCameron Smith              " block size: Min= %f MB; Max= %f MB; Avg= %f MB; Tot= %f MB\n",
365103be424SCameron Smith              rate, extra_msg, sizemin,sizemax,sizeavg,sizetot);
36659599516SKenneth E. Jansen     }
36759599516SKenneth E. Jansen   }
36859599516SKenneth E. Jansen   else {
36959599516SKenneth E. Jansen     if(irank == 0) {
37059599516SKenneth E. Jansen       printf(" \n");
37159599516SKenneth E. Jansen       //printf(" (%s) \n", extra_msg);
37259599516SKenneth E. Jansen     }
37359599516SKenneth E. Jansen   }
37459599516SKenneth E. Jansen }
37559599516SKenneth E. Jansen 
37659599516SKenneth E. Jansen /**
37759599516SKenneth E. Jansen  * This function is normally called at the beginning of a read operation, before
37859599516SKenneth E. Jansen  * init function.
37959599516SKenneth E. Jansen  * This function (uses rank 0) reads out nfields, nppf, master header size,
38059599516SKenneth E. Jansen  * endianess and allocates for masterHeader string.
38159599516SKenneth E. Jansen  * These values are essential for following read operations. Rank 0 will bcast
38259599516SKenneth E. Jansen  * these values to other ranks in the commm world
38359599516SKenneth E. Jansen  *
38459599516SKenneth E. Jansen  * If the file set is of old POSIX format, it would throw error and exit
38559599516SKenneth E. Jansen  */
38659599516SKenneth E. Jansen void queryphmpiio(const char filename[],int *nfields, int *nppf)
38759599516SKenneth E. Jansen {
38859599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &irank);
38959599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &mysize);
39059599516SKenneth E. Jansen 
39159599516SKenneth E. Jansen   if(irank == 0) {
39259599516SKenneth E. Jansen     FILE * fileHandle;
39359599516SKenneth E. Jansen     char* fname = StringStripper( filename );
39459599516SKenneth E. Jansen 
395*f42e0444SCameron Smith     double t0 = phiotmrc();
39659599516SKenneth E. Jansen     fileHandle = fopen (fname,"rb");
397*f42e0444SCameron Smith     phastaio_addOpenTime(phiotmrc()-t0);
39859599516SKenneth E. Jansen     if (fileHandle == NULL ) {
39959599516SKenneth E. Jansen       printf("\nError: File %s doesn't exist! Please check!\n",fname);
40059599516SKenneth E. Jansen     }
40159599516SKenneth E. Jansen     else {
40259599516SKenneth E. Jansen       SerialFile =(serial_file *)calloc( 1,  sizeof( serial_file) );
40359599516SKenneth E. Jansen       int meta_size_limit = VERSION_INFO_HEADER_SIZE;
40459599516SKenneth E. Jansen       SerialFile->masterHeader = (char *)malloc( meta_size_limit );
40559599516SKenneth E. Jansen       fread(SerialFile->masterHeader, 1, meta_size_limit, fileHandle);
40659599516SKenneth E. Jansen 
40759599516SKenneth E. Jansen       char read_out_tag[MAX_FIELDS_NAME_LENGTH];
40859599516SKenneth E. Jansen       char version[MAX_FIELDS_NAME_LENGTH/4];
40959599516SKenneth E. Jansen       int mhsize;
41059599516SKenneth E. Jansen       char * token;
41159599516SKenneth E. Jansen       int magic_number;
41259599516SKenneth E. Jansen 
41359599516SKenneth E. Jansen       memcpy( read_out_tag,
41459599516SKenneth E. Jansen           SerialFile->masterHeader,
41559599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
41659599516SKenneth E. Jansen 
41759599516SKenneth E. Jansen       if ( cscompare ("MPI_IO_Tag",read_out_tag) ) {
41859599516SKenneth E. Jansen         // Test endianess ...
41959599516SKenneth E. Jansen         memcpy (&magic_number,
42059599516SKenneth E. Jansen             SerialFile->masterHeader + sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
42159599516SKenneth E. Jansen             sizeof(int) );                                        // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
42259599516SKenneth E. Jansen 
42359599516SKenneth E. Jansen         if ( magic_number != ENDIAN_TEST_NUMBER ) {
42459599516SKenneth E. Jansen           printf("Endian is different!\n");
42559599516SKenneth E. Jansen           // Will do swap later
42659599516SKenneth E. Jansen         }
42759599516SKenneth E. Jansen 
42859599516SKenneth E. Jansen         // test version, old version, default masterheader size is 4M
42959599516SKenneth E. Jansen         // newer version masterheader size is read from first line
43059599516SKenneth E. Jansen         memcpy(version,
43159599516SKenneth E. Jansen             SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/2,
43259599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4 - 1); //TODO: why -1?
43359599516SKenneth E. Jansen 
43459599516SKenneth E. Jansen         if( cscompare ("version",version) ) {
43559599516SKenneth E. Jansen           // if there is "version" tag in the file, then it is newer format
43659599516SKenneth E. Jansen           // read master header size from here, otherwise use default
43759599516SKenneth E. Jansen           // Note: if version is "1", we know mhsize is at 3/4 place...
43859599516SKenneth E. Jansen 
43959599516SKenneth E. Jansen           token = strtok(version, ":");
44059599516SKenneth E. Jansen           token = strtok(NULL, " ,;<>" );
44159599516SKenneth E. Jansen           int iversion = atoi(token);
44259599516SKenneth E. Jansen 
44359599516SKenneth E. Jansen           if( iversion == 1) {
44459599516SKenneth E. Jansen             memcpy( &mhsize,
44559599516SKenneth E. Jansen                 SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/4*3 + sizeof("mhsize : ")-1,
44659599516SKenneth E. Jansen                 sizeof(int));
44759599516SKenneth E. Jansen             if ( magic_number != ENDIAN_TEST_NUMBER )
44859599516SKenneth E. Jansen               SwapArrayByteOrder(&mhsize, sizeof(int), 1);
44959599516SKenneth E. Jansen 
45059599516SKenneth E. Jansen             if( mhsize > DefaultMHSize ) {
45159599516SKenneth E. Jansen               //if actual headersize is larger than default, let's re-read
45259599516SKenneth E. Jansen               free(SerialFile->masterHeader);
45359599516SKenneth E. Jansen               SerialFile->masterHeader = (char *)malloc(mhsize);
45459599516SKenneth E. Jansen               fseek(fileHandle, 0, SEEK_SET); // reset the file stream position
45559599516SKenneth E. Jansen               fread(SerialFile->masterHeader,1,mhsize,fileHandle);
45659599516SKenneth E. Jansen             }
45759599516SKenneth E. Jansen           }
45859599516SKenneth E. Jansen           //TODO: check if this is a valid int??
45959599516SKenneth E. Jansen           MasterHeaderSize = mhsize;
46059599516SKenneth E. Jansen         }
46159599516SKenneth E. Jansen         else { // else it's version 0's format w/o version tag, implicating MHSize=4M
46259599516SKenneth E. Jansen           MasterHeaderSize = DefaultMHSize;
46359599516SKenneth E. Jansen         }
46459599516SKenneth E. Jansen 
46559599516SKenneth E. Jansen         memcpy( read_out_tag,
46659599516SKenneth E. Jansen             SerialFile->masterHeader+MAX_FIELDS_NAME_LENGTH+1,
46759599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH ); //TODO: why +1
46859599516SKenneth E. Jansen 
46959599516SKenneth E. Jansen         // Read in # fields ...
47059599516SKenneth E. Jansen         token = strtok( read_out_tag, ":" );
47159599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
47259599516SKenneth E. Jansen         *nfields = atoi( token );
47359599516SKenneth E. Jansen         if ( *nfields > MAX_FIELDS_NUMBER) {
47459599516SKenneth E. Jansen           printf("Error queryphmpiio: nfields is larger than MAX_FIELDS_NUMBER!\n");
47559599516SKenneth E. Jansen         }
47659599516SKenneth E. Jansen         SerialFile->nfields=*nfields; //TODO: sanity check of this int?
47759599516SKenneth E. Jansen 
47859599516SKenneth E. Jansen         memcpy( read_out_tag,
47959599516SKenneth E. Jansen             SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH * 2
48059599516SKenneth E. Jansen             + *nfields * MAX_FIELDS_NAME_LENGTH,
48159599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
48259599516SKenneth E. Jansen 
48359599516SKenneth E. Jansen         token = strtok( read_out_tag, ":" );
48459599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
48559599516SKenneth E. Jansen         *nppf = atoi( token );
48659599516SKenneth E. Jansen         SerialFile->nppf=*nppf; //TODO: sanity check of int
48759599516SKenneth E. Jansen       } // end of if("MPI_IO_TAG")
48859599516SKenneth E. Jansen       else {
48959599516SKenneth E. Jansen         printf("Error queryphmpiio: The file you opened is not of syncIO new format, please check! read_out_tag = %s\n",read_out_tag);
49059599516SKenneth E. Jansen         exit(1);
49159599516SKenneth E. Jansen       }
492*f42e0444SCameron Smith       double t0 = phiotmrc();
49359599516SKenneth E. Jansen       fclose(fileHandle);
494*f42e0444SCameron Smith       phastaio_addCloseTime(phiotmrc()-t0);
49559599516SKenneth E. Jansen       free(SerialFile->masterHeader);
49659599516SKenneth E. Jansen       free(SerialFile);
49759599516SKenneth E. Jansen     } //end of else
49859599516SKenneth E. Jansen     free(fname);
49959599516SKenneth E. Jansen   }
50059599516SKenneth E. Jansen 
50159599516SKenneth E. Jansen   // Bcast value to every one
50259599516SKenneth E. Jansen   MPI_Bcast( nfields, 1, MPI_INT, 0, MPI_COMM_WORLD );
50359599516SKenneth E. Jansen   MPI_Bcast( nppf, 1, MPI_INT, 0, MPI_COMM_WORLD );
50459599516SKenneth E. Jansen   MPI_Bcast( &MasterHeaderSize, 1, MPI_INT, 0, MPI_COMM_WORLD );
505a746c838SCameron Smith   phprintf("Info queryphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
50659599516SKenneth E. Jansen }
50759599516SKenneth E. Jansen 
50859599516SKenneth E. Jansen /**
50959599516SKenneth E. Jansen  * This function computes the right master header size (round to size of 2^n).
51059599516SKenneth E. Jansen  * This is only needed for file format version 1 in "write" mode.
51159599516SKenneth E. Jansen  */
51259599516SKenneth E. Jansen int computeMHSize(int nfields, int nppf, int version) {
51310d56689SCameron Smith   int mhsize=0;
51459599516SKenneth E. Jansen   if(version == 1) {
51559599516SKenneth E. Jansen     //int meta_info_size = (2+nfields+1) * MAX_FIELDS_NAME_LENGTH; // 2 is MPI_IO_TAG and nFields, the others 1 is nppf
51659599516SKenneth E. Jansen     int meta_info_size = VERSION_INFO_HEADER_SIZE;
5172dd307a1SCameron Smith     int actual_size =  nfields * nppf * sizeof(unsigned long) + meta_info_size;
51859599516SKenneth E. Jansen     //printf("actual_size = %d, offset table size = %d\n", actual_size,  nfields * nppf * sizeof(long long));
51959599516SKenneth E. Jansen     if (actual_size > DefaultMHSize) {
52059599516SKenneth E. Jansen       mhsize = (int) ceil( (double) actual_size/DefaultMHSize); // it's rounded to ceiling of this value
52159599516SKenneth E. Jansen       mhsize *= DefaultMHSize;
52259599516SKenneth E. Jansen     }
52359599516SKenneth E. Jansen     else {
52459599516SKenneth E. Jansen       mhsize = DefaultMHSize;
52559599516SKenneth E. Jansen     }
526b744cbcaSCameron Smith   } else {
527b744cbcaSCameron Smith     int rank = 0;
528b744cbcaSCameron Smith     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
529b744cbcaSCameron Smith     if(!rank) {
530b744cbcaSCameron Smith       fprintf(stderr,
531b744cbcaSCameron Smith           "ERROR invalid version passed to %s... exiting\n", __func__);
532b744cbcaSCameron Smith       exit(EXIT_FAILURE);
533b744cbcaSCameron Smith     }
53459599516SKenneth E. Jansen   }
53559599516SKenneth E. Jansen   return mhsize;
53659599516SKenneth E. Jansen }
53759599516SKenneth E. Jansen 
53859599516SKenneth E. Jansen /**
53959599516SKenneth E. Jansen  * Computes correct color of a rank according to number of files.
54059599516SKenneth E. Jansen  */
54159599516SKenneth E. Jansen extern "C" int computeColor( int myrank, int numprocs, int nfiles) {
54259599516SKenneth E. Jansen   int color =
54359599516SKenneth E. Jansen     (int)(myrank / (numprocs / nfiles));
54459599516SKenneth E. Jansen   return color;
54559599516SKenneth E. Jansen }
54659599516SKenneth E. Jansen 
54759599516SKenneth E. Jansen 
54859599516SKenneth E. Jansen /**
54959599516SKenneth E. Jansen  * Check the file descriptor.
55059599516SKenneth E. Jansen  */
551103be424SCameron Smith void checkFileDescriptor(const char fctname[], int*  fileDescriptor ) {
55259599516SKenneth E. Jansen   if ( *fileDescriptor < 0 ) {
55359599516SKenneth E. Jansen     printf("Error: File descriptor = %d in %s\n",*fileDescriptor,fctname);
55459599516SKenneth E. Jansen     exit(1);
55559599516SKenneth E. Jansen   }
55659599516SKenneth E. Jansen }
55759599516SKenneth E. Jansen 
55859599516SKenneth E. Jansen /**
55959599516SKenneth E. Jansen  * Initialize the file struct members and allocate space for file struct
56059599516SKenneth E. Jansen  * buffers.
56159599516SKenneth E. Jansen  *
56259599516SKenneth E. Jansen  * Note: this function is only called when we are using new format. Old POSIX
56359599516SKenneth E. Jansen  * format should skip this routine and call openfile() directly instead.
56459599516SKenneth E. Jansen  */
56559599516SKenneth E. Jansen int initphmpiio( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[])
56659599516SKenneth E. Jansen {
56759599516SKenneth E. Jansen   // we init irank again in case query not called (e.g. syncIO write case)
56859599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &irank);
56959599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &mysize);
57059599516SKenneth E. Jansen 
571a746c838SCameron Smith   phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d, nfields %d, nppf %d, nfiles %d\n", irank, MasterHeaderSize, *nfields, *nppf, *nfiles);
57259599516SKenneth E. Jansen 
57359599516SKenneth E. Jansen   double timer_start, timer_end;
57459599516SKenneth E. Jansen   startTimer(&timer_start);
57559599516SKenneth E. Jansen 
57659599516SKenneth E. Jansen   char* imode = StringStripper( mode );
57759599516SKenneth E. Jansen 
57859599516SKenneth E. Jansen   // Note: if it's read, we presume query was called prior to init and
57959599516SKenneth E. Jansen   // MasterHeaderSize is already set to correct value from parsing header
58059599516SKenneth E. Jansen   // otherwise it's write then it needs some computation to be set
58159599516SKenneth E. Jansen   if ( cscompare( "read", imode ) ) {
58259599516SKenneth E. Jansen     // do nothing
58359599516SKenneth E. Jansen   }
58459599516SKenneth E. Jansen   else if( cscompare( "write", imode ) ) {
58559599516SKenneth E. Jansen     MasterHeaderSize =  computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
58659599516SKenneth E. Jansen   }
58759599516SKenneth E. Jansen   else {
58859599516SKenneth E. Jansen     printf("Error initphmpiio: can't recognize the mode %s", imode);
58959599516SKenneth E. Jansen     exit(1);
59059599516SKenneth E. Jansen   }
59159599516SKenneth E. Jansen   free ( imode );
59259599516SKenneth E. Jansen 
593a746c838SCameron Smith   phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
59459599516SKenneth E. Jansen 
59559599516SKenneth E. Jansen   int i, j;
59659599516SKenneth E. Jansen 
59759599516SKenneth E. Jansen   if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
59859599516SKenneth E. Jansen     printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
59959599516SKenneth E. Jansen     endTimer(&timer_end);
60059599516SKenneth E. Jansen     printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
60159599516SKenneth E. Jansen     return MAX_PHASTA_FILES_EXCEEDED;
60259599516SKenneth E. Jansen   }
60359599516SKenneth E. Jansen   //		else if( PhastaIONextActiveIndex == 0 )  //Hang in debug mode on Intrepid
60459599516SKenneth E. Jansen   //		{
60559599516SKenneth E. Jansen   //			for( i = 0; i < MAX_PHASTA_FILES; i++ );
60659599516SKenneth E. Jansen   //			{
60759599516SKenneth E. Jansen   //				PhastaIOActiveFiles[i] = NULL;
60859599516SKenneth E. Jansen   //			}
60959599516SKenneth E. Jansen   //		}
61059599516SKenneth E. Jansen 
61159599516SKenneth E. Jansen 
61259599516SKenneth E. Jansen   PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1,  sizeof( phastaio_file_t) );
61359599516SKenneth E. Jansen 
61459599516SKenneth E. Jansen   i = PhastaIONextActiveIndex;
61559599516SKenneth E. Jansen   PhastaIONextActiveIndex++;
61659599516SKenneth E. Jansen 
61759599516SKenneth E. Jansen   //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
61859599516SKenneth E. Jansen 
61959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize;  // what does this mean??? TODO
62059599516SKenneth E. Jansen 
62159599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->Wrong_Endian = false;
62259599516SKenneth E. Jansen 
62359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFields = *nfields;
62459599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nPPF = *nppf;
62559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFiles = *nfiles;
62659599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->myrank));
62759599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->numprocs));
62859599516SKenneth E. Jansen 
62959599516SKenneth E. Jansen 
63059599516SKenneth E. Jansen   if( *nfiles > 1 ) { // split the ranks according to each mpiio file
63159599516SKenneth E. Jansen 
63259599516SKenneth E. Jansen     if ( s_assign_local_comm == 0) { // call mpi_comm_split for the first (and only) time
63359599516SKenneth E. Jansen 
63459599516SKenneth E. Jansen       if (PhastaIOActiveFiles[i]->myrank == 0) printf("Building subcommunicator\n");
63559599516SKenneth E. Jansen 
63659599516SKenneth E. Jansen       int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
63759599516SKenneth E. Jansen       MPI_Comm_split(MPI_COMM_WORLD,
63859599516SKenneth E. Jansen           color,
63959599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->myrank,
64059599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_comm));
64159599516SKenneth E. Jansen       MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
64259599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_numprocs));
64359599516SKenneth E. Jansen       MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
64459599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_myrank));
64559599516SKenneth E. Jansen 
64659599516SKenneth E. Jansen       // back up now these variables so that we do not need to call comm_split again
64759599516SKenneth E. Jansen       s_local_comm = PhastaIOActiveFiles[i]->local_comm;
64859599516SKenneth E. Jansen       s_local_size = PhastaIOActiveFiles[i]->local_numprocs;
64959599516SKenneth E. Jansen       s_local_rank = PhastaIOActiveFiles[i]->local_myrank;
65059599516SKenneth E. Jansen       s_assign_local_comm = 1;
65159599516SKenneth E. Jansen     }
65259599516SKenneth E. Jansen     else { // recycle the subcommunicator
65359599516SKenneth E. Jansen       if (PhastaIOActiveFiles[i]->myrank == 0) printf("Recycling subcommunicator\n");
65459599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_comm = s_local_comm;
65559599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_numprocs = s_local_size;
65659599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_myrank = s_local_rank;
65759599516SKenneth E. Jansen     }
65859599516SKenneth E. Jansen   }
65959599516SKenneth E. Jansen   else { // *nfiles == 1 here - no need to call mpi_comm_split here
66059599516SKenneth E. Jansen 
66159599516SKenneth E. Jansen     if (PhastaIOActiveFiles[i]->myrank == 0) printf("Bypassing subcommunicator\n");
66259599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_comm = MPI_COMM_WORLD;
66359599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_numprocs = PhastaIOActiveFiles[i]->numprocs;
66459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_myrank = PhastaIOActiveFiles[i]->myrank;
66559599516SKenneth E. Jansen 
66659599516SKenneth E. Jansen   }
66759599516SKenneth E. Jansen 
66859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nppp =
66959599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
67059599516SKenneth E. Jansen 
67159599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
67259599516SKenneth E. Jansen     (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
67359599516SKenneth E. Jansen     (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
67459599516SKenneth E. Jansen 
67559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_offset_table =
6762dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
67759599516SKenneth E. Jansen 
67859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_read_table =
6792dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
68059599516SKenneth E. Jansen 
68159599516SKenneth E. Jansen   for (j=0; j<*nfields; j++)
68259599516SKenneth E. Jansen   {
68359599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[j] =
6842dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
68559599516SKenneth E. Jansen 
68659599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_read_table[j] =
6872dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
68859599516SKenneth E. Jansen   }
68959599516SKenneth E. Jansen   *filehandle = i;
69059599516SKenneth E. Jansen 
69159599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
69259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
69359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
69459599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
69559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
69659599516SKenneth E. Jansen 
69759599516SKenneth E. Jansen   // Time monitoring
69859599516SKenneth E. Jansen   endTimer(&timer_end);
69959599516SKenneth E. Jansen   printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
70059599516SKenneth E. Jansen 
70159599516SKenneth E. Jansen   phprintf_0("Info initphmpiio: quiting function");
70259599516SKenneth E. Jansen 
70359599516SKenneth E. Jansen   return i;
70459599516SKenneth E. Jansen }
70559599516SKenneth E. Jansen 
70659599516SKenneth E. Jansen /**
70759599516SKenneth E. Jansen  * Destruct the file struct and free buffers allocated in init function.
70859599516SKenneth E. Jansen  */
70959599516SKenneth E. Jansen void finalizephmpiio( int *fileDescriptor )
71059599516SKenneth E. Jansen {
71159599516SKenneth E. Jansen   double timer_start, timer_end;
71259599516SKenneth E. Jansen   startTimer(&timer_start);
71359599516SKenneth E. Jansen 
71459599516SKenneth E. Jansen   int i, j;
71559599516SKenneth E. Jansen   i = *fileDescriptor;
71659599516SKenneth E. Jansen   //PhastaIONextActiveIndex--;
71759599516SKenneth E. Jansen 
71859599516SKenneth E. Jansen   /* //free the offset table for this phasta file */
71959599516SKenneth E. Jansen   //for(j=0; j<MAX_FIELDS_NUMBER; j++) //Danger: undefined behavior for my_*_table.[j] not allocated or not initialized to NULL
72059599516SKenneth E. Jansen   for(j=0; j<PhastaIOActiveFiles[i]->nFields; j++)
72159599516SKenneth E. Jansen   {
72259599516SKenneth E. Jansen     free( PhastaIOActiveFiles[i]->my_offset_table[j]);
72359599516SKenneth E. Jansen     free( PhastaIOActiveFiles[i]->my_read_table[j]);
72459599516SKenneth E. Jansen   }
72559599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->my_offset_table );
72659599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->my_read_table );
72759599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->master_header );
72859599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->double_chunk );
72959599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->int_chunk );
73059599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->read_double_chunk );
73159599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->read_int_chunk );
73259599516SKenneth E. Jansen 
733946a6cddSCameron Smith   if( PhastaIOActiveFiles[i]->nFiles > 1 && s_assign_local_comm ) { // the comm was split
734946a6cddSCameron Smith     if (PhastaIOActiveFiles[i]->myrank == 0) printf("Freeing subcommunicator\n");
735946a6cddSCameron Smith     s_assign_local_comm = 0;
736946a6cddSCameron Smith     MPI_Comm_free(&(PhastaIOActiveFiles[i]->local_comm));
737946a6cddSCameron Smith   }
738946a6cddSCameron Smith 
73959599516SKenneth E. Jansen   free( PhastaIOActiveFiles[i]);
74059599516SKenneth E. Jansen 
74159599516SKenneth E. Jansen   endTimer(&timer_end);
74259599516SKenneth E. Jansen   printPerf("finalizempiio", timer_start, timer_end, 0, 0, "");
74359599516SKenneth E. Jansen 
74459599516SKenneth E. Jansen   PhastaIONextActiveIndex--;
74559599516SKenneth E. Jansen }
74659599516SKenneth E. Jansen 
74759599516SKenneth E. Jansen 
74859599516SKenneth E. Jansen /**
74959599516SKenneth E. Jansen  * Special init for M2N in order to create a subcommunicator for the reduced solution (requires PRINT_PERF to be false for now)
75059599516SKenneth E. Jansen  * Initialize the file struct members and allocate space for file struct buffers.
75159599516SKenneth E. Jansen  *
75259599516SKenneth E. Jansen  * Note: this function is only called when we are using new format. Old POSIX
75359599516SKenneth E. Jansen  * format should skip this routine and call openfile() directly instead.
75459599516SKenneth E. Jansen  */
75559599516SKenneth E. Jansen int initphmpiiosub( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[],MPI_Comm my_local_comm)
75659599516SKenneth E. Jansen {
75759599516SKenneth E. Jansen   // we init irank again in case query not called (e.g. syncIO write case)
75859599516SKenneth E. Jansen 
75959599516SKenneth E. Jansen   MPI_Comm_rank(my_local_comm, &irank);
76059599516SKenneth E. Jansen   MPI_Comm_size(my_local_comm, &mysize);
76159599516SKenneth E. Jansen 
762a746c838SCameron Smith   phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
76359599516SKenneth E. Jansen 
76459599516SKenneth E. Jansen   double timer_start, timer_end;
76559599516SKenneth E. Jansen   startTimer(&timer_start);
76659599516SKenneth E. Jansen 
76759599516SKenneth E. Jansen   char* imode = StringStripper( mode );
76859599516SKenneth E. Jansen 
76959599516SKenneth E. Jansen   // Note: if it's read, we presume query was called prior to init and
77059599516SKenneth E. Jansen   // MasterHeaderSize is already set to correct value from parsing header
77159599516SKenneth E. Jansen   // otherwise it's write then it needs some computation to be set
77259599516SKenneth E. Jansen   if ( cscompare( "read", imode ) ) {
77359599516SKenneth E. Jansen     // do nothing
77459599516SKenneth E. Jansen   }
77559599516SKenneth E. Jansen   else if( cscompare( "write", imode ) ) {
77659599516SKenneth E. Jansen     MasterHeaderSize =  computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
77759599516SKenneth E. Jansen   }
77859599516SKenneth E. Jansen   else {
77959599516SKenneth E. Jansen     printf("Error initphmpiio: can't recognize the mode %s", imode);
78059599516SKenneth E. Jansen     exit(1);
78159599516SKenneth E. Jansen   }
78259599516SKenneth E. Jansen   free ( imode );
78359599516SKenneth E. Jansen 
784a746c838SCameron Smith   phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
78559599516SKenneth E. Jansen 
78659599516SKenneth E. Jansen   int i, j;
78759599516SKenneth E. Jansen 
78859599516SKenneth E. Jansen   if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
78959599516SKenneth E. Jansen     printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
79059599516SKenneth E. Jansen     endTimer(&timer_end);
79159599516SKenneth E. Jansen     printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
79259599516SKenneth E. Jansen     return MAX_PHASTA_FILES_EXCEEDED;
79359599516SKenneth E. Jansen   }
79459599516SKenneth E. Jansen   //		else if( PhastaIONextActiveIndex == 0 )  //Hang in debug mode on Intrepid
79559599516SKenneth E. Jansen   //		{
79659599516SKenneth E. Jansen   //			for( i = 0; i < MAX_PHASTA_FILES; i++ );
79759599516SKenneth E. Jansen   //			{
79859599516SKenneth E. Jansen   //				PhastaIOActiveFiles[i] = NULL;
79959599516SKenneth E. Jansen   //			}
80059599516SKenneth E. Jansen   //		}
80159599516SKenneth E. Jansen 
80259599516SKenneth E. Jansen 
80359599516SKenneth E. Jansen   PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1,  sizeof( phastaio_file_t) );
80459599516SKenneth E. Jansen 
80559599516SKenneth E. Jansen   i = PhastaIONextActiveIndex;
80659599516SKenneth E. Jansen   PhastaIONextActiveIndex++;
80759599516SKenneth E. Jansen 
80859599516SKenneth E. Jansen   //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
80959599516SKenneth E. Jansen 
81059599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize;  // what does this mean??? TODO
81159599516SKenneth E. Jansen 
81259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->Wrong_Endian = false;
81359599516SKenneth E. Jansen 
81459599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFields = *nfields;
81559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nPPF = *nppf;
81659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFiles = *nfiles;
81759599516SKenneth E. Jansen   MPI_Comm_rank(my_local_comm, &(PhastaIOActiveFiles[i]->myrank));
81859599516SKenneth E. Jansen   MPI_Comm_size(my_local_comm, &(PhastaIOActiveFiles[i]->numprocs));
81959599516SKenneth E. Jansen 
82059599516SKenneth E. Jansen   int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
82159599516SKenneth E. Jansen   MPI_Comm_split(my_local_comm,
82259599516SKenneth E. Jansen       color,
82359599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->myrank,
82459599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_comm));
82559599516SKenneth E. Jansen   MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
82659599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_numprocs));
82759599516SKenneth E. Jansen   MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
82859599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_myrank));
82959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nppp =
83059599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
83159599516SKenneth E. Jansen 
83259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
83359599516SKenneth E. Jansen     (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
83459599516SKenneth E. Jansen     (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
83559599516SKenneth E. Jansen 
83659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_offset_table =
8372dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
83859599516SKenneth E. Jansen 
83959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_read_table =
8402dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
84159599516SKenneth E. Jansen 
84259599516SKenneth E. Jansen   for (j=0; j<*nfields; j++)
84359599516SKenneth E. Jansen   {
84459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[j] =
8452dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
84659599516SKenneth E. Jansen 
84759599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_read_table[j] =
8482dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
84959599516SKenneth E. Jansen   }
85059599516SKenneth E. Jansen   *filehandle = i;
85159599516SKenneth E. Jansen 
85259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
85359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
85459599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
85559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
85659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
85759599516SKenneth E. Jansen 
85859599516SKenneth E. Jansen   // Time monitoring
85959599516SKenneth E. Jansen   endTimer(&timer_end);
86059599516SKenneth E. Jansen   printPerf("initphmpiiosub", timer_start, timer_end, 0, 0, "");
86159599516SKenneth E. Jansen 
86259599516SKenneth E. Jansen   phprintf_0("Info initphmpiiosub: quiting function");
86359599516SKenneth E. Jansen 
86459599516SKenneth E. Jansen   return i;
86559599516SKenneth E. Jansen }
86659599516SKenneth E. Jansen 
86759599516SKenneth E. Jansen 
86859599516SKenneth E. Jansen 
86959599516SKenneth E. Jansen /** open file for both POSIX and MPI-IO syncIO format.
87059599516SKenneth E. Jansen  *
87159599516SKenneth E. Jansen  * If it's old POSIX format, simply call posix fopen().
87259599516SKenneth E. Jansen  *
87359599516SKenneth E. Jansen  * If it's MPI-IO foramt:
87459599516SKenneth E. Jansen  * in "read" mode, it builds the header table that points to the offset of
87559599516SKenneth E. Jansen  * fields for parts;
87659599516SKenneth E. Jansen  * in "write" mode, it opens the file with MPI-IO open routine.
87759599516SKenneth E. Jansen  */
878103be424SCameron Smith void openfile(const char filename[], const char mode[], int*  fileDescriptor )
87959599516SKenneth E. Jansen {
88059599516SKenneth E. Jansen   phprintf_0("Info: entering openfile");
88159599516SKenneth E. Jansen 
88259599516SKenneth E. Jansen   double timer_start, timer_end;
88359599516SKenneth E. Jansen   startTimer(&timer_start);
88459599516SKenneth E. Jansen 
88559599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 )
88659599516SKenneth E. Jansen   {
88759599516SKenneth E. Jansen     FILE* file=NULL ;
88859599516SKenneth E. Jansen     *fileDescriptor = 0;
88959599516SKenneth E. Jansen     char* fname = StringStripper( filename );
89059599516SKenneth E. Jansen     char* imode = StringStripper( mode );
89159599516SKenneth E. Jansen 
892*f42e0444SCameron Smith     double t0 = phiotmrc();
89359599516SKenneth E. Jansen     if ( cscompare( "read", imode ) ) file = fopen(fname, "rb" );
89459599516SKenneth E. Jansen     else if( cscompare( "write", imode ) ) file = fopen(fname, "wb" );
89559599516SKenneth E. Jansen     else if( cscompare( "append", imode ) ) file = fopen(fname, "ab" );
896*f42e0444SCameron Smith     phastaio_addOpenTime(phiotmrc()-t0);
89759599516SKenneth E. Jansen 
89859599516SKenneth E. Jansen     if ( !file ){
89967701978SCameron Smith       fprintf(stderr,"Error openfile: unable to open file %s\n",fname);
90059599516SKenneth E. Jansen     } else {
90159599516SKenneth E. Jansen       fileArray.push_back( file );
90259599516SKenneth E. Jansen       byte_order.push_back( false );
90359599516SKenneth E. Jansen       header_type.push_back( sizeof(int) );
90459599516SKenneth E. Jansen       *fileDescriptor = fileArray.size();
90559599516SKenneth E. Jansen     }
90659599516SKenneth E. Jansen     free (fname);
90759599516SKenneth E. Jansen     free (imode);
90859599516SKenneth E. Jansen   }
90959599516SKenneth E. Jansen   else // else it would be parallel I/O, opposed to posix io
91059599516SKenneth E. Jansen   {
91159599516SKenneth E. Jansen     char* fname = StringStripper( filename );
91259599516SKenneth E. Jansen     char* imode = StringStripper( mode );
91359599516SKenneth E. Jansen     int rc;
91459599516SKenneth E. Jansen     int i = *fileDescriptor;
91559599516SKenneth E. Jansen     checkFileDescriptor("openfile",&i);
91659599516SKenneth E. Jansen     char* token;
91759599516SKenneth E. Jansen 
91859599516SKenneth E. Jansen     if ( cscompare( "read", imode ) )
91959599516SKenneth E. Jansen     {
92059599516SKenneth E. Jansen       //	      if (PhastaIOActiveFiles[i]->myrank == 0)
92159599516SKenneth E. Jansen       //                printf("\n **********\nRead open ... ... regular version\n");
92259599516SKenneth E. Jansen 
923*f42e0444SCameron Smith       double t0 = phiotmrc();
92459599516SKenneth E. Jansen       rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
92559599516SKenneth E. Jansen           fname,
92659599516SKenneth E. Jansen           MPI_MODE_RDONLY,
92759599516SKenneth E. Jansen           MPI_INFO_NULL,
92859599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->file_handle) );
929*f42e0444SCameron Smith       phastaio_addOpenTime(phiotmrc()-t0);
93059599516SKenneth E. Jansen 
93159599516SKenneth E. Jansen       if(rc)
93259599516SKenneth E. Jansen       {
93359599516SKenneth E. Jansen         *fileDescriptor = UNABLE_TO_OPEN_FILE;
9344d60bba2SCameron Smith         int error_string_length;
9354d60bba2SCameron Smith         char error_string[4096];
9364d60bba2SCameron Smith         MPI_Error_string(rc, error_string, &error_string_length);
9374d60bba2SCameron Smith         fprintf(stderr, "Error openfile: Unable to open file %s! MPI reports \"%s\"\n",fname,error_string);
93859599516SKenneth E. Jansen         endTimer(&timer_end);
93959599516SKenneth E. Jansen         printPerf("openfile", timer_start, timer_end, 0, 0, "");
94059599516SKenneth E. Jansen         return;
94159599516SKenneth E. Jansen       }
94259599516SKenneth E. Jansen 
94359599516SKenneth E. Jansen       MPI_Status read_tag_status;
94459599516SKenneth E. Jansen       char read_out_tag[MAX_FIELDS_NAME_LENGTH];
94559599516SKenneth E. Jansen       int j;
94659599516SKenneth E. Jansen       int magic_number;
94759599516SKenneth E. Jansen 
94859599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
94959599516SKenneth E. Jansen         MPI_File_read_at( PhastaIOActiveFiles[i]->file_handle,
95059599516SKenneth E. Jansen             0,
95159599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header,
95259599516SKenneth E. Jansen             MasterHeaderSize,
95359599516SKenneth E. Jansen             MPI_CHAR,
95459599516SKenneth E. Jansen             &read_tag_status );
95559599516SKenneth E. Jansen       }
95659599516SKenneth E. Jansen 
95759599516SKenneth E. Jansen       MPI_Bcast( PhastaIOActiveFiles[i]->master_header,
95859599516SKenneth E. Jansen           MasterHeaderSize,
95959599516SKenneth E. Jansen           MPI_CHAR,
96059599516SKenneth E. Jansen           0,
96159599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->local_comm );
96259599516SKenneth E. Jansen 
96359599516SKenneth E. Jansen       memcpy( read_out_tag,
96459599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->master_header,
96559599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
96659599516SKenneth E. Jansen 
96759599516SKenneth E. Jansen       if ( cscompare ("MPI_IO_Tag",read_out_tag) )
96859599516SKenneth E. Jansen       {
96959599516SKenneth E. Jansen         // Test endianess ...
97059599516SKenneth E. Jansen         memcpy ( &magic_number,
97159599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
97259599516SKenneth E. Jansen             sizeof(int) );                                                   // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
97359599516SKenneth E. Jansen 
97459599516SKenneth E. Jansen         if ( magic_number != ENDIAN_TEST_NUMBER )
97559599516SKenneth E. Jansen         {
97659599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->Wrong_Endian = true;
97759599516SKenneth E. Jansen         }
97859599516SKenneth E. Jansen 
97959599516SKenneth E. Jansen         memcpy( read_out_tag,
98059599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH+1, // TODO: WHY +1???
98159599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH );
98259599516SKenneth E. Jansen 
98359599516SKenneth E. Jansen         // Read in # fields ...
98459599516SKenneth E. Jansen         token = strtok ( read_out_tag, ":" );
98559599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
98659599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->nFields = atoi( token );
98759599516SKenneth E. Jansen 
9882dd307a1SCameron Smith         unsigned long **header_table;
9892dd307a1SCameron Smith         header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
99059599516SKenneth E. Jansen 
99159599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
99259599516SKenneth E. Jansen         {
9932dd307a1SCameron Smith           header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
99459599516SKenneth E. Jansen         }
99559599516SKenneth E. Jansen 
99659599516SKenneth E. Jansen         // Read in the offset table ...
99759599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
99859599516SKenneth E. Jansen         {
99959599516SKenneth E. Jansen           if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
100059599516SKenneth E. Jansen             memcpy( header_table[j],
100159599516SKenneth E. Jansen                 PhastaIOActiveFiles[i]->master_header +
100259599516SKenneth E. Jansen                 VERSION_INFO_HEADER_SIZE +
10032dd307a1SCameron Smith                 j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
10042dd307a1SCameron Smith                 PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
100559599516SKenneth E. Jansen           }
100659599516SKenneth E. Jansen 
100759599516SKenneth E. Jansen           MPI_Scatter( header_table[j],
100859599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->nppp,
100959599516SKenneth E. Jansen               MPI_LONG_LONG_INT,
101059599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->my_read_table[j],
101159599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->nppp,
101259599516SKenneth E. Jansen               MPI_LONG_LONG_INT,
101359599516SKenneth E. Jansen               0,
101459599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->local_comm );
101559599516SKenneth E. Jansen 
101659599516SKenneth E. Jansen           // Swap byte order if endianess is different ...
101759599516SKenneth E. Jansen           if ( PhastaIOActiveFiles[i]->Wrong_Endian ) {
101859599516SKenneth E. Jansen             SwapArrayByteOrder( PhastaIOActiveFiles[i]->my_read_table[j],
10194187599aSCameron Smith                 sizeof(unsigned long),
102059599516SKenneth E. Jansen                 PhastaIOActiveFiles[i]->nppp );
102159599516SKenneth E. Jansen           }
102259599516SKenneth E. Jansen         }
102359599516SKenneth E. Jansen 
102459599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
102559599516SKenneth E. Jansen           free ( header_table[j] );
102659599516SKenneth E. Jansen         }
102759599516SKenneth E. Jansen         free (header_table);
102859599516SKenneth E. Jansen 
102959599516SKenneth E. Jansen       } // end of if MPI_IO_TAG
103059599516SKenneth E. Jansen       else //else not valid MPI file
103159599516SKenneth E. Jansen       {
103259599516SKenneth E. Jansen         *fileDescriptor = NOT_A_MPI_FILE;
103359599516SKenneth E. Jansen         printf("Error openfile: The file %s you opened is not in syncIO new format, please check again! File descriptor = %d, MasterHeaderSize = %d, read_out_tag = %s\n",fname,*fileDescriptor,MasterHeaderSize,read_out_tag);
103459599516SKenneth E. Jansen         //Printing MasterHeaderSize is useful to test a compiler bug on Intrepid BGP
103559599516SKenneth E. Jansen         endTimer(&timer_end);
103659599516SKenneth E. Jansen         printPerf("openfile", timer_start, timer_end, 0, 0, "");
103759599516SKenneth E. Jansen         return;
103859599516SKenneth E. Jansen       }
103959599516SKenneth E. Jansen     } // end of if "read"
104059599516SKenneth E. Jansen     else if( cscompare( "write", imode ) )
104159599516SKenneth E. Jansen     {
1042*f42e0444SCameron Smith       double t0 = phiotmrc();
104359599516SKenneth E. Jansen       rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
104459599516SKenneth E. Jansen           fname,
104559599516SKenneth E. Jansen           MPI_MODE_WRONLY | MPI_MODE_CREATE,
104659599516SKenneth E. Jansen           MPI_INFO_NULL,
104759599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->file_handle) );
1048*f42e0444SCameron Smith       phastaio_addOpenTime(phiotmrc()-t0);
10494d60bba2SCameron Smith       if(rc != MPI_SUCCESS)
105059599516SKenneth E. Jansen       {
105159599516SKenneth E. Jansen         *fileDescriptor = UNABLE_TO_OPEN_FILE;
10524d60bba2SCameron Smith         int error_string_length;
10534d60bba2SCameron Smith         char error_string[4096];
10544d60bba2SCameron Smith         MPI_Error_string(rc, error_string, &error_string_length);
10554d60bba2SCameron Smith         fprintf(stderr, "Error openfile: Unable to open file %s! MPI reports \"%s\"\n",fname,error_string);
105659599516SKenneth E. Jansen         return;
105759599516SKenneth E. Jansen       }
105859599516SKenneth E. Jansen     } // end of if "write"
105959599516SKenneth E. Jansen     free (fname);
106059599516SKenneth E. Jansen     free (imode);
106159599516SKenneth E. Jansen   } // end of if FileIndex != 0
106259599516SKenneth E. Jansen 
106359599516SKenneth E. Jansen   endTimer(&timer_end);
106459599516SKenneth E. Jansen   printPerf("openfile", timer_start, timer_end, 0, 0, "");
106559599516SKenneth E. Jansen }
106659599516SKenneth E. Jansen 
106759599516SKenneth E. Jansen /** close file for both POSIX and MPI-IO syncIO format.
106859599516SKenneth E. Jansen  *
106959599516SKenneth E. Jansen  * If it's old POSIX format, simply call posix fclose().
107059599516SKenneth E. Jansen  *
107159599516SKenneth E. Jansen  * If it's MPI-IO foramt:
107259599516SKenneth E. Jansen  * in "read" mode, it simply close file with MPI-IO close routine.
107359599516SKenneth E. Jansen  * in "write" mode, rank 0 in each group will re-assemble the master header and
107459599516SKenneth E. Jansen  * offset table and write to the beginning of file, then close the file.
107559599516SKenneth E. Jansen  */
1076103be424SCameron Smith void closefile( int* fileDescriptor, const char mode[] )
107759599516SKenneth E. Jansen {
107859599516SKenneth E. Jansen   double timer_start, timer_end;
107959599516SKenneth E. Jansen   startTimer(&timer_start);
108059599516SKenneth E. Jansen 
108159599516SKenneth E. Jansen   int i = *fileDescriptor;
108259599516SKenneth E. Jansen   checkFileDescriptor("closefile",&i);
108359599516SKenneth E. Jansen 
108459599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
108559599516SKenneth E. Jansen     char* imode = StringStripper( mode );
108659599516SKenneth E. Jansen 
108759599516SKenneth E. Jansen     if( cscompare( "write", imode )
108859599516SKenneth E. Jansen         || cscompare( "append", imode ) ) {
108959599516SKenneth E. Jansen       fflush( fileArray[ *fileDescriptor - 1 ] );
109059599516SKenneth E. Jansen     }
109159599516SKenneth E. Jansen 
1092*f42e0444SCameron Smith     double t0 = phiotmrc();
109359599516SKenneth E. Jansen     fclose( fileArray[ *fileDescriptor - 1 ] );
1094*f42e0444SCameron Smith     phastaio_addCloseTime(phiotmrc()-t0);
109559599516SKenneth E. Jansen     free (imode);
109659599516SKenneth E. Jansen   }
109759599516SKenneth E. Jansen   else {
109859599516SKenneth E. Jansen     char* imode = StringStripper( mode );
109959599516SKenneth E. Jansen 
110059599516SKenneth E. Jansen     //write master header here:
110159599516SKenneth E. Jansen     if ( cscompare( "write", imode ) ) {
110259599516SKenneth E. Jansen       //	      if ( PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields < 2*ONE_MEGABYTE/8 )  //SHOULD BE CHECKED
110359599516SKenneth E. Jansen       //		MasterHeaderSize = 4*ONE_MEGABYTE;
110459599516SKenneth E. Jansen       //	      else
110559599516SKenneth E. Jansen       //		MasterHeaderSize = 4*ONE_MEGABYTE + PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields * 8 - 2*ONE_MEGABYTE;
110659599516SKenneth E. Jansen 
110759599516SKenneth E. Jansen       MasterHeaderSize = computeMHSize( PhastaIOActiveFiles[i]->nFields, PhastaIOActiveFiles[i]->nPPF, LATEST_WRITE_VERSION);
110859599516SKenneth E. Jansen       phprintf_0("Info closefile: myrank = %d, MasterHeaderSize = %d\n", PhastaIOActiveFiles[i]->myrank, MasterHeaderSize);
110959599516SKenneth E. Jansen 
111059599516SKenneth E. Jansen       MPI_Status write_header_status;
111159599516SKenneth E. Jansen       char mpi_tag[MAX_FIELDS_NAME_LENGTH];
111259599516SKenneth E. Jansen       char version[MAX_FIELDS_NAME_LENGTH/4];
111359599516SKenneth E. Jansen       char mhsize[MAX_FIELDS_NAME_LENGTH/4];
111459599516SKenneth E. Jansen       int magic_number = ENDIAN_TEST_NUMBER;
111559599516SKenneth E. Jansen 
111659599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 )
111759599516SKenneth E. Jansen       {
111859599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
111959599516SKenneth E. Jansen         sprintf(mpi_tag, "MPI_IO_Tag : ");
112059599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header,
112159599516SKenneth E. Jansen             mpi_tag,
112259599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
112359599516SKenneth E. Jansen 
112459599516SKenneth E. Jansen         bzero((void*)version,MAX_FIELDS_NAME_LENGTH/4);
112559599516SKenneth E. Jansen         // this version is "1", print version in ASCII
112659599516SKenneth E. Jansen         sprintf(version, "version : %d",1);
112759599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/2,
112859599516SKenneth E. Jansen             version,
112959599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4);
113059599516SKenneth E. Jansen 
113159599516SKenneth E. Jansen         // master header size is computed using the formula above
113259599516SKenneth E. Jansen         bzero((void*)mhsize,MAX_FIELDS_NAME_LENGTH/4);
113359599516SKenneth E. Jansen         sprintf(mhsize, "mhsize : ");
113459599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/4*3,
113559599516SKenneth E. Jansen             mhsize,
113659599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4);
113759599516SKenneth E. Jansen 
113859599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
113959599516SKenneth E. Jansen         sprintf(mpi_tag,
114059599516SKenneth E. Jansen             "\nnFields : %d\n",
114159599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nFields);
114259599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH,
114359599516SKenneth E. Jansen             mpi_tag,
114459599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
114559599516SKenneth E. Jansen 
114659599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
114759599516SKenneth E. Jansen         sprintf(mpi_tag, "\nnPPF : %d\n", PhastaIOActiveFiles[i]->nPPF);
114859599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+
114959599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nFields *
115059599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH +
115159599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH * 2,
115259599516SKenneth E. Jansen             mpi_tag,
115359599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
115459599516SKenneth E. Jansen 
115559599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
115659599516SKenneth E. Jansen             &magic_number,
115759599516SKenneth E. Jansen             sizeof(int));
115859599516SKenneth E. Jansen 
115959599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("mhsize : ") -1 + MAX_FIELDS_NAME_LENGTH/4*3,
116059599516SKenneth E. Jansen             &MasterHeaderSize,
116159599516SKenneth E. Jansen             sizeof(int));
116259599516SKenneth E. Jansen       }
116359599516SKenneth E. Jansen 
116459599516SKenneth E. Jansen       int j = 0;
11652dd307a1SCameron Smith       unsigned long **header_table;
11662dd307a1SCameron Smith       header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
116759599516SKenneth E. Jansen 
116859599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
11692dd307a1SCameron Smith         header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
117059599516SKenneth E. Jansen       }
117159599516SKenneth E. Jansen 
117259599516SKenneth E. Jansen       //if( irank == 0 ) printf("gonna mpi_gather, myrank = %d\n", irank);
117359599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
117459599516SKenneth E. Jansen         MPI_Gather( PhastaIOActiveFiles[i]->my_offset_table[j],
117559599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nppp,
117659599516SKenneth E. Jansen             MPI_LONG_LONG_INT,
117759599516SKenneth E. Jansen             header_table[j],
117859599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nppp,
117959599516SKenneth E. Jansen             MPI_LONG_LONG_INT,
118059599516SKenneth E. Jansen             0,
118159599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->local_comm );
118259599516SKenneth E. Jansen       }
118359599516SKenneth E. Jansen 
118459599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
118559599516SKenneth E. Jansen 
118659599516SKenneth E. Jansen         //if( irank == 0 ) printf("gonna memcpy for every procs, myrank = %d\n", irank);
118759599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
118859599516SKenneth E. Jansen           memcpy ( PhastaIOActiveFiles[i]->master_header +
118959599516SKenneth E. Jansen               VERSION_INFO_HEADER_SIZE +
11902dd307a1SCameron Smith               j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
119159599516SKenneth E. Jansen               header_table[j],
11922dd307a1SCameron Smith               PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
119359599516SKenneth E. Jansen         }
119459599516SKenneth E. Jansen 
119559599516SKenneth E. Jansen         //if( irank == 0 ) printf("gonna file_write_at(), myrank = %d\n", irank);
119659599516SKenneth E. Jansen         MPI_File_write_at( PhastaIOActiveFiles[i]->file_handle,
119759599516SKenneth E. Jansen             0,
119859599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header,
119959599516SKenneth E. Jansen             MasterHeaderSize,
120059599516SKenneth E. Jansen             MPI_CHAR,
120159599516SKenneth E. Jansen             &write_header_status );
120259599516SKenneth E. Jansen       }
120359599516SKenneth E. Jansen 
120459599516SKenneth E. Jansen       ////free(PhastaIOActiveFiles[i]->master_header);
120559599516SKenneth E. Jansen 
120659599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
120759599516SKenneth E. Jansen         free ( header_table[j] );
120859599516SKenneth E. Jansen       }
120959599516SKenneth E. Jansen       free (header_table);
121059599516SKenneth E. Jansen     }
121159599516SKenneth E. Jansen 
121259599516SKenneth E. Jansen     //if( irank == 0 ) printf("gonna file_close(), myrank = %d\n", irank);
1213*f42e0444SCameron Smith     double t0 = phiotmrc();
121459599516SKenneth E. Jansen     MPI_File_close( &( PhastaIOActiveFiles[i]->file_handle ) );
1215*f42e0444SCameron Smith     phastaio_addCloseTime(phiotmrc()-t0);
121659599516SKenneth E. Jansen     free ( imode );
121759599516SKenneth E. Jansen   }
121859599516SKenneth E. Jansen 
121959599516SKenneth E. Jansen   endTimer(&timer_end);
122059599516SKenneth E. Jansen   printPerf("closefile_", timer_start, timer_end, 0, 0, "");
122159599516SKenneth E. Jansen }
122259599516SKenneth E. Jansen 
12233872e963SCameron Smith int readHeader( FILE* f, const char phrase[],
12243872e963SCameron Smith     int* params, int numParams, const char* iotype) {
12253872e963SCameron Smith   isBinary(iotype);
12263872e963SCameron Smith   return readHeader(f,phrase,params,numParams);
12273872e963SCameron Smith }
12283872e963SCameron Smith 
1229103be424SCameron Smith void readheader(
1230103be424SCameron Smith     int* fileDescriptor,
123159599516SKenneth E. Jansen     const  char keyphrase[],
123259599516SKenneth E. Jansen     void* valueArray,
123359599516SKenneth E. Jansen     int*  nItems,
123459599516SKenneth E. Jansen     const char  datatype[],
123559599516SKenneth E. Jansen     const char  iotype[] )
123659599516SKenneth E. Jansen {
123759599516SKenneth E. Jansen   double timer_start, timer_end;
1238d3337298SCameron Smith 
123959599516SKenneth E. Jansen   startTimer(&timer_start);
124059599516SKenneth E. Jansen 
124159599516SKenneth E. Jansen   int i = *fileDescriptor;
124259599516SKenneth E. Jansen   checkFileDescriptor("readheader",&i);
124359599516SKenneth E. Jansen 
124459599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
124559599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
124659599516SKenneth E. Jansen     FILE* fileObject;
124759599516SKenneth E. Jansen     int* valueListInt;
124859599516SKenneth E. Jansen 
124959599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
125059599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
125159599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
125259599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
125359599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
125459599516SKenneth E. Jansen       endTimer(&timer_end);
125559599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
125659599516SKenneth E. Jansen       return;
125759599516SKenneth E. Jansen     }
125859599516SKenneth E. Jansen 
1259961a4ff6SCameron Smith     LastHeaderKey[filePtr] = std::string(keyphrase);
126059599516SKenneth E. Jansen     LastHeaderNotFound = false;
126159599516SKenneth E. Jansen 
126259599516SKenneth E. Jansen     fileObject = fileArray[ filePtr ] ;
126359599516SKenneth E. Jansen     Wrong_Endian = byte_order[ filePtr ];
126459599516SKenneth E. Jansen 
126559599516SKenneth E. Jansen     isBinary( iotype );
126659599516SKenneth E. Jansen     typeSize( datatype );   //redundant call, just avoid a compiler warning.
126759599516SKenneth E. Jansen 
126859599516SKenneth E. Jansen     // right now we are making the assumption that we will only write integers
126959599516SKenneth E. Jansen     // on the header line.
127059599516SKenneth E. Jansen 
127159599516SKenneth E. Jansen     valueListInt = static_cast< int* >( valueArray );
127259599516SKenneth E. Jansen     int ierr = readHeader( fileObject ,
127359599516SKenneth E. Jansen         keyphrase,
127459599516SKenneth E. Jansen         valueListInt,
127559599516SKenneth E. Jansen         *nItems ) ;
127659599516SKenneth E. Jansen 
127759599516SKenneth E. Jansen     byte_order[ filePtr ] = Wrong_Endian ;
127859599516SKenneth E. Jansen 
127959599516SKenneth E. Jansen     if ( ierr ) LastHeaderNotFound = true;
128059599516SKenneth E. Jansen 
128159599516SKenneth E. Jansen     //return ; // don't return, go to the end to print perf
128259599516SKenneth E. Jansen   }
128359599516SKenneth E. Jansen   else {
128459599516SKenneth E. Jansen     int* valueListInt;
128559599516SKenneth E. Jansen     valueListInt = static_cast <int*>(valueArray);
1286400e9fc0SCameron Smith     char* token = NULL;
128759599516SKenneth E. Jansen     bool FOUND = false ;
128859599516SKenneth E. Jansen     isBinary( iotype );
128959599516SKenneth E. Jansen 
129059599516SKenneth E. Jansen     MPI_Status read_offset_status;
129159599516SKenneth E. Jansen     char read_out_tag[MAX_FIELDS_NAME_LENGTH];
1292400e9fc0SCameron Smith     memset(read_out_tag, '\0', MAX_FIELDS_NAME_LENGTH);
129359599516SKenneth E. Jansen     char readouttag[MAX_FIELDS_NUMBER][MAX_FIELDS_NAME_LENGTH];
129459599516SKenneth E. Jansen     int j;
129559599516SKenneth E. Jansen 
129659599516SKenneth E. Jansen     int string_length = strlen( keyphrase );
129759599516SKenneth E. Jansen     char* buffer = (char*) malloc ( string_length+1 );
129859599516SKenneth E. Jansen 
129959599516SKenneth E. Jansen     strcpy ( buffer, keyphrase );
130059599516SKenneth E. Jansen     buffer[ string_length ] = '\0';
130159599516SKenneth E. Jansen 
130259599516SKenneth E. Jansen     char* st2 = strtok ( buffer, "@" );
130359599516SKenneth E. Jansen     st2 = strtok (NULL, "@");
130459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = atoi(st2);
130559599516SKenneth E. Jansen     if ( char* p = strpbrk(buffer, "@") )
130659599516SKenneth E. Jansen       *p = '\0';
130759599516SKenneth E. Jansen 
130859599516SKenneth E. Jansen     // Check if the user has input the right GPid
130959599516SKenneth E. Jansen     if ( ( PhastaIOActiveFiles[i]->GPid <=
131059599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->myrank *
131159599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->nppp )||
131259599516SKenneth E. Jansen         ( PhastaIOActiveFiles[i]->GPid >
131359599516SKenneth E. Jansen           ( PhastaIOActiveFiles[i]->myrank + 1 ) *
131459599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->nppp ) )
131559599516SKenneth E. Jansen     {
131659599516SKenneth E. Jansen       *fileDescriptor = NOT_A_MPI_FILE;
131759599516SKenneth E. Jansen       printf("Error readheader: The file is not in syncIO new format, please check! myrank = %d, GPid = %d, nppp = %d, keyphrase = %s\n", PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->GPid, PhastaIOActiveFiles[i]->nppp, keyphrase);
131859599516SKenneth E. Jansen       // It is possible atoi could not generate a clear integer from st2 because of additional garbage character in keyphrase
131959599516SKenneth E. Jansen       endTimer(&timer_end);
132059599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
132159599516SKenneth E. Jansen       return;
132259599516SKenneth E. Jansen     }
132359599516SKenneth E. Jansen 
132459599516SKenneth E. Jansen     // Find the field we want ...
132559599516SKenneth E. Jansen     //for ( j = 0; j<MAX_FIELDS_NUMBER; j++ )
132659599516SKenneth E. Jansen     for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
132759599516SKenneth E. Jansen     {
132859599516SKenneth E. Jansen       memcpy( readouttag[j],
132959599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->master_header + j*MAX_FIELDS_NAME_LENGTH+MAX_FIELDS_NAME_LENGTH*2+1,
133059599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
133159599516SKenneth E. Jansen     }
133259599516SKenneth E. Jansen 
133359599516SKenneth E. Jansen     for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
133459599516SKenneth E. Jansen     {
133559599516SKenneth E. Jansen       token = strtok ( readouttag[j], ":" );
133659599516SKenneth E. Jansen 
133759599516SKenneth E. Jansen       //if ( cscompare( buffer, token ) )
133859599516SKenneth E. Jansen       if ( cscompare( token , buffer ) && cscompare( buffer, token ) )
133959599516SKenneth E. Jansen         // This double comparison is required for the field "number of nodes" and all the other fields that start with "number of nodes" (i.g. number of nodes in the mesh").
134059599516SKenneth E. Jansen         // Would be safer to rename "number of nodes" by "number of nodes in the part" so that the name are completely unique. But much more work to do that (Nspre, phParAdapt, etc).
134159599516SKenneth E. Jansen         // Since the field name are unique in SyncIO (as it includes part ID), this should be safe and there should be no issue with the "?" trailing character.
134259599516SKenneth E. Jansen       {
134359599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->read_field_count = j;
134459599516SKenneth E. Jansen         FOUND = true;
134559599516SKenneth E. Jansen         //printf("buffer: %s | token: %s | j: %d\n",buffer,token,j);
134659599516SKenneth E. Jansen         break;
134759599516SKenneth E. Jansen       }
134859599516SKenneth E. Jansen     }
134959599516SKenneth E. Jansen     free(buffer);
135059599516SKenneth E. Jansen 
135159599516SKenneth E. Jansen     if (!FOUND)
135259599516SKenneth E. Jansen     {
135359599516SKenneth E. Jansen       //if(irank==0) printf("Warning readheader: Not found %s \n",keyphrase); //PhastaIOActiveFiles[i]->myrank is certainly initialized here.
135459599516SKenneth E. Jansen       if(PhastaIOActiveFiles[i]->myrank == 0) printf("WARNING readheader: Not found %s\n",keyphrase);
135559599516SKenneth E. Jansen       endTimer(&timer_end);
135659599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
135759599516SKenneth E. Jansen       return;
135859599516SKenneth E. Jansen     }
135959599516SKenneth E. Jansen 
136059599516SKenneth E. Jansen     // Find the part we want ...
136159599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->read_part_count = PhastaIOActiveFiles[i]->GPid -
136259599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->myrank * PhastaIOActiveFiles[i]->nppp - 1;
136359599516SKenneth E. Jansen 
136459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset =
136559599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->my_read_table[PhastaIOActiveFiles[i]->read_field_count][PhastaIOActiveFiles[i]->read_part_count];
136659599516SKenneth E. Jansen 
136759599516SKenneth E. Jansen     // 	  printf("****Rank %d offset is %d\n",PhastaIOActiveFiles[i]->myrank,PhastaIOActiveFiles[i]->my_offset);
136859599516SKenneth E. Jansen 
136959599516SKenneth E. Jansen     // Read each datablock header here ...
137059599516SKenneth E. Jansen 
137159599516SKenneth E. Jansen     MPI_File_read_at_all( PhastaIOActiveFiles[i]->file_handle,
137259599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->my_offset+1,
137359599516SKenneth E. Jansen         read_out_tag,
137459599516SKenneth E. Jansen         MAX_FIELDS_NAME_LENGTH-1,
137559599516SKenneth E. Jansen         MPI_CHAR,
137659599516SKenneth E. Jansen         &read_offset_status );
137759599516SKenneth E. Jansen     token = strtok ( read_out_tag, ":" );
137859599516SKenneth E. Jansen 
137959599516SKenneth E. Jansen     // 	  printf("&&&&Rank %d read_out_tag is %s\n",PhastaIOActiveFiles[i]->myrank,read_out_tag);
138059599516SKenneth E. Jansen 
138159599516SKenneth E. Jansen     if( cscompare( keyphrase , token ) ) //No need to compare also token with keyphrase like above. We should already have the right one. Otherwise there is a problem.
138259599516SKenneth E. Jansen     {
138359599516SKenneth E. Jansen       FOUND = true ;
138459599516SKenneth E. Jansen       token = strtok( NULL, " ,;<>" );
138559599516SKenneth E. Jansen       for( j=0; j < *nItems && ( token = strtok( NULL," ,;<>") ); j++ )
138659599516SKenneth E. Jansen         valueListInt[j] = atoi( token );
138759599516SKenneth E. Jansen 
138859599516SKenneth E. Jansen       if ( j < *nItems )
138959599516SKenneth E. Jansen       {
139059599516SKenneth E. Jansen         fprintf( stderr, "Expected # of ints not found for: %s\n", keyphrase );
139159599516SKenneth E. Jansen       }
139259599516SKenneth E. Jansen     }
139359599516SKenneth E. Jansen     else {
139459599516SKenneth E. Jansen       //if(irank==0)
139559599516SKenneth E. Jansen       if(PhastaIOActiveFiles[i]->myrank == 0)
139659599516SKenneth E. Jansen         // If we enter this if, there is a problem with the name of some fields
139759599516SKenneth E. Jansen       {
139859599516SKenneth E. Jansen         printf("Error readheader: Unexpected mismatch between keyphrase = %s and token = %s\n",keyphrase,token);
139959599516SKenneth E. Jansen       }
140059599516SKenneth E. Jansen     }
140159599516SKenneth E. Jansen   }
140259599516SKenneth E. Jansen 
140359599516SKenneth E. Jansen   endTimer(&timer_end);
140459599516SKenneth E. Jansen   printPerf("readheader", timer_start, timer_end, 0, 0, "");
140559599516SKenneth E. Jansen 
140659599516SKenneth E. Jansen }
140759599516SKenneth E. Jansen 
14083956dcfeSCameron Smith void readDataBlock(
14093956dcfeSCameron Smith     FILE* fileObject,
14103956dcfeSCameron Smith     void* valueArray,
14113956dcfeSCameron Smith     int nItems,
14123956dcfeSCameron Smith     const char  datatype[],
14133956dcfeSCameron Smith     const char  iotype[] )
14143956dcfeSCameron Smith {
14153956dcfeSCameron Smith   isBinary(iotype);
14163956dcfeSCameron Smith   size_t type_size = typeSize( datatype );
1417*f42e0444SCameron Smith   double t0 = phiotmrc();
14183956dcfeSCameron Smith   if ( binary_format ) {
14193956dcfeSCameron Smith     char junk = '\0';
14203956dcfeSCameron Smith     fread( valueArray, type_size, nItems, fileObject );
14213956dcfeSCameron Smith     fread( &junk, sizeof(char), 1 , fileObject );
14223956dcfeSCameron Smith     if ( Wrong_Endian ) SwapArrayByteOrder( valueArray, type_size, nItems );
14233956dcfeSCameron Smith   } else {
14243956dcfeSCameron Smith     char* ts1 = StringStripper( datatype );
14253956dcfeSCameron Smith     if ( cscompare( "integer", ts1 ) ) {
14263956dcfeSCameron Smith       for( int n=0; n < nItems ; n++ )
14273956dcfeSCameron Smith         fscanf(fileObject, "%d\n",(int*)((int*)valueArray+n) );
14283956dcfeSCameron Smith     } else if ( cscompare( "double", ts1 ) ) {
14293956dcfeSCameron Smith       for( int n=0; n < nItems ; n++ )
14303956dcfeSCameron Smith         fscanf(fileObject, "%lf\n",(double*)((double*)valueArray+n) );
14313956dcfeSCameron Smith     }
14323956dcfeSCameron Smith     free (ts1);
14333956dcfeSCameron Smith   }
1434*f42e0444SCameron Smith   phastaio_addReadTime(phiotmrc()-t0);
1435*f42e0444SCameron Smith   phastaio_addReadBytes(nItems*type_size);
14363956dcfeSCameron Smith }
14373956dcfeSCameron Smith 
1438103be424SCameron Smith void readdatablock(
1439103be424SCameron Smith     int*  fileDescriptor,
144059599516SKenneth E. Jansen     const char keyphrase[],
144159599516SKenneth E. Jansen     void* valueArray,
144259599516SKenneth E. Jansen     int*  nItems,
144359599516SKenneth E. Jansen     const char  datatype[],
144459599516SKenneth E. Jansen     const char  iotype[] )
144559599516SKenneth E. Jansen {
144659599516SKenneth E. Jansen   //if(irank == 0) printf("entering readdatablock()\n");
14472dd307a1SCameron Smith   unsigned long data_size = 0;
144859599516SKenneth E. Jansen   double timer_start, timer_end;
144959599516SKenneth E. Jansen   startTimer(&timer_start);
145059599516SKenneth E. Jansen 
145159599516SKenneth E. Jansen   int i = *fileDescriptor;
145259599516SKenneth E. Jansen   checkFileDescriptor("readdatablock",&i);
145359599516SKenneth E. Jansen 
145459599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
145559599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
145659599516SKenneth E. Jansen     FILE* fileObject;
145759599516SKenneth E. Jansen 
145859599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
145959599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
146059599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before\n") ;
146159599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
146259599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
146359599516SKenneth E. Jansen       endTimer(&timer_end);
146459599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
146559599516SKenneth E. Jansen       return;
146659599516SKenneth E. Jansen     }
146759599516SKenneth E. Jansen 
146859599516SKenneth E. Jansen     // error check..
146959599516SKenneth E. Jansen     // since we require that a consistant header always preceed the data block
147059599516SKenneth E. Jansen     // let us check to see that it is actually the case.
147159599516SKenneth E. Jansen 
1472961a4ff6SCameron Smith     if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
147359599516SKenneth E. Jansen       fprintf(stderr, "Header not consistant with data block\n");
1474961a4ff6SCameron Smith       fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
147559599516SKenneth E. Jansen       fprintf(stderr, "DataBlock: %s\n ", keyphrase );
147659599516SKenneth E. Jansen       fprintf(stderr, "Please recheck read sequence \n");
147759599516SKenneth E. Jansen       if( Strict_Error ) {
147859599516SKenneth E. Jansen         fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
147959599516SKenneth E. Jansen         endTimer(&timer_end);
148059599516SKenneth E. Jansen         printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
148159599516SKenneth E. Jansen         return;
148259599516SKenneth E. Jansen       }
148359599516SKenneth E. Jansen     }
148459599516SKenneth E. Jansen 
148559599516SKenneth E. Jansen     if ( LastHeaderNotFound ) {
148659599516SKenneth E. Jansen       endTimer(&timer_end);
148759599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
148859599516SKenneth E. Jansen       return;
148959599516SKenneth E. Jansen     }
149059599516SKenneth E. Jansen     fileObject = fileArray[ filePtr ];
149159599516SKenneth E. Jansen     Wrong_Endian = byte_order[ filePtr ];
1492bae968b9SCameron Smith     LastHeaderKey.erase(filePtr);
14933956dcfeSCameron Smith     readDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
149459599516SKenneth E. Jansen 
149559599516SKenneth E. Jansen     //return;
149659599516SKenneth E. Jansen   }
149759599516SKenneth E. Jansen   else {
149859599516SKenneth E. Jansen     // 	  printf("read data block\n");
149959599516SKenneth E. Jansen     MPI_Status read_data_status;
150059599516SKenneth E. Jansen     size_t type_size = typeSize( datatype );
150159599516SKenneth E. Jansen     int nUnits = *nItems;
150259599516SKenneth E. Jansen     isBinary( iotype );
150359599516SKenneth E. Jansen 
150459599516SKenneth E. Jansen     // read datablock then
150559599516SKenneth E. Jansen     //MR CHANGE
150659599516SKenneth E. Jansen     //             if ( cscompare ( datatype, "double"))
150759599516SKenneth E. Jansen     char* ts2 = StringStripper( datatype );
150859599516SKenneth E. Jansen     if ( cscompare ( "double" , ts2))
150959599516SKenneth E. Jansen       //MR CHANGE END
151059599516SKenneth E. Jansen     {
151159599516SKenneth E. Jansen 
1512*f42e0444SCameron Smith       double t0 = phiotmrc();
151359599516SKenneth E. Jansen       MPI_File_read_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
151459599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
151559599516SKenneth E. Jansen           valueArray,
151659599516SKenneth E. Jansen           nUnits,
151759599516SKenneth E. Jansen           MPI_DOUBLE );
151859599516SKenneth E. Jansen       MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
151959599516SKenneth E. Jansen           valueArray,
152059599516SKenneth E. Jansen           &read_data_status );
152159599516SKenneth E. Jansen       data_size=8*nUnits;
1522*f42e0444SCameron Smith       phastaio_addReadTime(phiotmrc()-t0);
1523*f42e0444SCameron Smith       phastaio_addReadBytes(nUnits*sizeof(double));
152459599516SKenneth E. Jansen     }
152559599516SKenneth E. Jansen     //MR CHANGE
152659599516SKenneth E. Jansen     //             else if ( cscompare ( datatype, "integer"))
152759599516SKenneth E. Jansen     else if ( cscompare ( "integer" , ts2))
152859599516SKenneth E. Jansen       //MR CHANGE END
152959599516SKenneth E. Jansen     {
1530*f42e0444SCameron Smith       double t0 = phiotmrc();
153159599516SKenneth E. Jansen       MPI_File_read_at_all_begin(PhastaIOActiveFiles[i]->file_handle,
153259599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
153359599516SKenneth E. Jansen           valueArray,
153459599516SKenneth E. Jansen           nUnits,
153559599516SKenneth E. Jansen           MPI_INT );
153659599516SKenneth E. Jansen       MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
153759599516SKenneth E. Jansen           valueArray,
153859599516SKenneth E. Jansen           &read_data_status );
153959599516SKenneth E. Jansen       data_size=4*nUnits;
1540*f42e0444SCameron Smith       phastaio_addReadTime(phiotmrc()-t0);
1541*f42e0444SCameron Smith       phastaio_addReadBytes(nUnits*sizeof(int));
154259599516SKenneth E. Jansen     }
154359599516SKenneth E. Jansen     else
154459599516SKenneth E. Jansen     {
154559599516SKenneth E. Jansen       *fileDescriptor = DATA_TYPE_ILLEGAL;
154659599516SKenneth E. Jansen       printf("readdatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
154759599516SKenneth E. Jansen       endTimer(&timer_end);
154859599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
154959599516SKenneth E. Jansen       return;
155059599516SKenneth E. Jansen     }
155159599516SKenneth E. Jansen     free(ts2);
155259599516SKenneth E. Jansen 
155359599516SKenneth E. Jansen 
155459599516SKenneth E. Jansen     // 	  printf("%d Read finishe\n",PhastaIOActiveFiles[i]->myrank);
155559599516SKenneth E. Jansen 
155659599516SKenneth E. Jansen     // Swap data byte order if endianess is different ...
155759599516SKenneth E. Jansen     if ( PhastaIOActiveFiles[i]->Wrong_Endian )
155859599516SKenneth E. Jansen     {
155959599516SKenneth E. Jansen       SwapArrayByteOrder( valueArray, type_size, nUnits );
156059599516SKenneth E. Jansen     }
156159599516SKenneth E. Jansen   }
156259599516SKenneth E. Jansen 
156359599516SKenneth E. Jansen   endTimer(&timer_end);
156459599516SKenneth E. Jansen   char extra_msg[1024];
156559599516SKenneth E. Jansen   memset(extra_msg, '\0', 1024);
156659599516SKenneth E. Jansen   char* key = StringStripper(keyphrase);
156759599516SKenneth E. Jansen   sprintf(extra_msg, " field is %s ", key);
156859599516SKenneth E. Jansen   printPerf("readdatablock", timer_start, timer_end, data_size, 1, extra_msg);
156959599516SKenneth E. Jansen   free(key);
157059599516SKenneth E. Jansen 
157159599516SKenneth E. Jansen }
157259599516SKenneth E. Jansen 
1573103be424SCameron Smith void writeHeader(
1574103be424SCameron Smith     FILE* f,
1575ea868eb1SCameron Smith     const char keyphrase[],
1576ea868eb1SCameron Smith     const void* valueArray,
1577ea868eb1SCameron Smith     const int nItems,
1578ea868eb1SCameron Smith     const int ndataItems,
1579ea868eb1SCameron Smith     const char datatype[],
1580ea868eb1SCameron Smith     const char iotype[])
1581ea868eb1SCameron Smith {
1582ea868eb1SCameron Smith   isBinary( iotype );
1583ea868eb1SCameron Smith 
1584ea868eb1SCameron Smith   const int _newline =
1585ea868eb1SCameron Smith     ( ndataItems > 0 ) ? sizeof( char ) : 0;
1586ea868eb1SCameron Smith   int size_of_nextblock =
1587ea868eb1SCameron Smith     ( binary_format ) ? typeSize(datatype) * ndataItems + _newline : ndataItems;
1588ea868eb1SCameron Smith 
1589ea868eb1SCameron Smith   fprintf( f, "%s : < %d > ", keyphrase, size_of_nextblock );
1590ea868eb1SCameron Smith   for( int i = 0; i < nItems; i++ )
1591ea868eb1SCameron Smith     fprintf( f, "%d ", *((int*)((int*)valueArray+i)));
1592ea868eb1SCameron Smith   fprintf( f, "\n");
1593ea868eb1SCameron Smith }
1594ea868eb1SCameron Smith 
1595103be424SCameron Smith void writeheader(
1596103be424SCameron Smith     const int* fileDescriptor,
159759599516SKenneth E. Jansen     const char keyphrase[],
159859599516SKenneth E. Jansen     const void* valueArray,
159959599516SKenneth E. Jansen     const int* nItems,
160059599516SKenneth E. Jansen     const int* ndataItems,
160159599516SKenneth E. Jansen     const char datatype[],
160259599516SKenneth E. Jansen     const char iotype[])
160359599516SKenneth E. Jansen {
160459599516SKenneth E. Jansen 
160559599516SKenneth E. Jansen   //if(irank == 0) printf("entering writeheader()\n");
160659599516SKenneth E. Jansen 
160759599516SKenneth E. Jansen   double timer_start, timer_end;
160859599516SKenneth E. Jansen   startTimer(&timer_start);
160959599516SKenneth E. Jansen 
161059599516SKenneth E. Jansen   int i = *fileDescriptor;
161159599516SKenneth E. Jansen   checkFileDescriptor("writeheader",&i);
161259599516SKenneth E. Jansen 
161359599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
161459599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
161559599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
161659599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
161759599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
161859599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
161959599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
162059599516SKenneth E. Jansen       endTimer(&timer_end);
162159599516SKenneth E. Jansen       printPerf("writeheader", timer_start, timer_end, 0, 0, "");
162259599516SKenneth E. Jansen       return;
162359599516SKenneth E. Jansen     }
162459599516SKenneth E. Jansen 
1625961a4ff6SCameron Smith     LastHeaderKey[filePtr] = std::string(keyphrase);
162659599516SKenneth E. Jansen     DataSize = *ndataItems;
1627ea868eb1SCameron Smith     FILE* fileObject = fileArray[ filePtr ] ;
1628ea868eb1SCameron Smith     header_type[ filePtr ] = typeSize( datatype );
1629ea868eb1SCameron Smith     writeHeader(fileObject,keyphrase,valueArray,*nItems,
1630ea868eb1SCameron Smith         *ndataItems,datatype,iotype);
163159599516SKenneth E. Jansen   }
163259599516SKenneth E. Jansen   else { // else it's parallel I/O
163359599516SKenneth E. Jansen     DataSize = *ndataItems;
163459599516SKenneth E. Jansen     size_t type_size = typeSize( datatype );
163559599516SKenneth E. Jansen     isBinary( iotype );
163659599516SKenneth E. Jansen     char mpi_tag[MAX_FIELDS_NAME_LENGTH];
163759599516SKenneth E. Jansen 
163859599516SKenneth E. Jansen     int string_length = strlen( keyphrase );
163959599516SKenneth E. Jansen     char* buffer = (char*) malloc ( string_length+1 );
164059599516SKenneth E. Jansen 
164159599516SKenneth E. Jansen     strcpy ( buffer, keyphrase);
164259599516SKenneth E. Jansen     buffer[ string_length ] = '\0';
164359599516SKenneth E. Jansen 
164459599516SKenneth E. Jansen     char* st2 = strtok ( buffer, "@" );
164559599516SKenneth E. Jansen     st2 = strtok (NULL, "@");
164659599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = atoi(st2);
164759599516SKenneth E. Jansen 
164859599516SKenneth E. Jansen     if ( char* p = strpbrk(buffer, "@") )
164959599516SKenneth E. Jansen       *p = '\0';
165059599516SKenneth E. Jansen 
165159599516SKenneth E. Jansen     bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
165259599516SKenneth E. Jansen     sprintf(mpi_tag, "\n%s : %d\n", buffer, PhastaIOActiveFiles[i]->field_count);
16532dd307a1SCameron Smith     unsigned long offset_value;
165459599516SKenneth E. Jansen 
165559599516SKenneth E. Jansen     int temp = *ndataItems;
16562dd307a1SCameron Smith     unsigned long number_of_items = (unsigned long)temp;
165759599516SKenneth E. Jansen     MPI_Barrier(PhastaIOActiveFiles[i]->local_comm);
165859599516SKenneth E. Jansen 
165959599516SKenneth E. Jansen     MPI_Scan( &number_of_items,
166059599516SKenneth E. Jansen         &offset_value,
166159599516SKenneth E. Jansen         1,
166259599516SKenneth E. Jansen         MPI_LONG_LONG_INT,
166359599516SKenneth E. Jansen         MPI_SUM,
166459599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_comm );
166559599516SKenneth E. Jansen 
166659599516SKenneth E. Jansen     offset_value = (offset_value - number_of_items) * type_size;
166759599516SKenneth E. Jansen 
166859599516SKenneth E. Jansen     offset_value += PhastaIOActiveFiles[i]->local_myrank *
166959599516SKenneth E. Jansen       DB_HEADER_SIZE +
167059599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->next_start_address;
167159599516SKenneth E. Jansen     // This offset is the starting address of each datablock header...
167259599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset = offset_value;
167359599516SKenneth E. Jansen 
167459599516SKenneth E. Jansen     // Write in my offset table ...
167559599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[PhastaIOActiveFiles[i]->field_count][PhastaIOActiveFiles[i]->part_count] =
167659599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->my_offset;
167759599516SKenneth E. Jansen 
167859599516SKenneth E. Jansen     // Update the next-start-address ...
167959599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->next_start_address = offset_value +
168059599516SKenneth E. Jansen       number_of_items * type_size +
168159599516SKenneth E. Jansen       DB_HEADER_SIZE;
168259599516SKenneth E. Jansen     MPI_Bcast( &(PhastaIOActiveFiles[i]->next_start_address),
168359599516SKenneth E. Jansen         1,
168459599516SKenneth E. Jansen         MPI_LONG_LONG_INT,
168559599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_numprocs-1,
168659599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_comm );
168759599516SKenneth E. Jansen 
168859599516SKenneth E. Jansen     // Prepare datablock header ...
168959599516SKenneth E. Jansen     int _newline = (*ndataItems>0)?sizeof(char):0;
169059599516SKenneth E. Jansen     unsigned int size_of_nextblock = type_size * (*ndataItems) + _newline;
169159599516SKenneth E. Jansen 
169259599516SKenneth E. Jansen     //char datablock_header[255];
169359599516SKenneth E. Jansen     //bzero((void*)datablock_header,255);
169459599516SKenneth E. Jansen     char datablock_header[DB_HEADER_SIZE];
169559599516SKenneth E. Jansen     bzero((void*)datablock_header,DB_HEADER_SIZE);
169659599516SKenneth E. Jansen 
169759599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = PhastaIOActiveFiles[i]->nppp*PhastaIOActiveFiles[i]->myrank+PhastaIOActiveFiles[i]->part_count;
169859599516SKenneth E. Jansen     sprintf( datablock_header,
169959599516SKenneth E. Jansen         "\n%s : < %u >",
170059599516SKenneth E. Jansen         keyphrase,
170159599516SKenneth E. Jansen         size_of_nextblock );
170259599516SKenneth E. Jansen 
170359599516SKenneth E. Jansen     for ( int j = 0; j < *nItems; j++ )
170459599516SKenneth E. Jansen     {
170559599516SKenneth E. Jansen       sprintf( datablock_header,
170659599516SKenneth E. Jansen           "%s %d ",
170759599516SKenneth E. Jansen           datablock_header,
170859599516SKenneth E. Jansen           *((int*)((int*)valueArray+j)));
170959599516SKenneth E. Jansen     }
171059599516SKenneth E. Jansen     sprintf( datablock_header,
171159599516SKenneth E. Jansen         "%s\n ",
171259599516SKenneth E. Jansen         datablock_header );
171359599516SKenneth E. Jansen 
171459599516SKenneth E. Jansen     // Write datablock header ...
171559599516SKenneth E. Jansen     //MR CHANGE
171659599516SKenneth E. Jansen     // 	if ( cscompare(datatype,"double") )
171759599516SKenneth E. Jansen     char* ts1 = StringStripper( datatype );
171859599516SKenneth E. Jansen     if ( cscompare("double",ts1) )
171959599516SKenneth E. Jansen       //MR CHANGE END
172059599516SKenneth E. Jansen     {
172159599516SKenneth E. Jansen       free ( PhastaIOActiveFiles[i]->double_chunk );
172259599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->double_chunk = ( double * )malloc( (sizeof( double )*number_of_items+ DB_HEADER_SIZE));
172359599516SKenneth E. Jansen 
172459599516SKenneth E. Jansen       double * aa = ( double * )datablock_header;
172559599516SKenneth E. Jansen       memcpy(PhastaIOActiveFiles[i]->double_chunk, aa, DB_HEADER_SIZE);
172659599516SKenneth E. Jansen     }
172759599516SKenneth E. Jansen     //MR CHANGE
172859599516SKenneth E. Jansen     // 	if  ( cscompare(datatype,"integer") )
172959599516SKenneth E. Jansen     else if ( cscompare("integer",ts1) )
173059599516SKenneth E. Jansen       //MR CHANGE END
173159599516SKenneth E. Jansen     {
173259599516SKenneth E. Jansen       free ( PhastaIOActiveFiles[i]->int_chunk );
173359599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->int_chunk = ( int * )malloc( (sizeof( int )*number_of_items+ DB_HEADER_SIZE));
173459599516SKenneth E. Jansen 
173559599516SKenneth E. Jansen       int * aa = ( int * )datablock_header;
173659599516SKenneth E. Jansen       memcpy(PhastaIOActiveFiles[i]->int_chunk, aa, DB_HEADER_SIZE);
173759599516SKenneth E. Jansen     }
173859599516SKenneth E. Jansen     else {
173959599516SKenneth E. Jansen       //             *fileDescriptor = DATA_TYPE_ILLEGAL;
174059599516SKenneth E. Jansen       printf("writeheader - DATA_TYPE_ILLEGAL - %s\n",datatype);
174159599516SKenneth E. Jansen       endTimer(&timer_end);
174259599516SKenneth E. Jansen       printPerf("writeheader", timer_start, timer_end, 0, 0, "");
174359599516SKenneth E. Jansen       return;
174459599516SKenneth E. Jansen     }
174559599516SKenneth E. Jansen     free(ts1);
174659599516SKenneth E. Jansen 
174759599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->part_count++;
174859599516SKenneth E. Jansen     if ( PhastaIOActiveFiles[i]->part_count == PhastaIOActiveFiles[i]->nppp ) {
174959599516SKenneth E. Jansen       //A new field will be written
175059599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
175159599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header +
175259599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->field_count *
175359599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH +
175459599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH * 2,
175559599516SKenneth E. Jansen             mpi_tag,
175659599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH-1);
175759599516SKenneth E. Jansen       }
175859599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->field_count++;
175959599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->part_count=0;
176059599516SKenneth E. Jansen     }
176159599516SKenneth E. Jansen     free(buffer);
176259599516SKenneth E. Jansen   }
176359599516SKenneth E. Jansen 
176459599516SKenneth E. Jansen   endTimer(&timer_end);
176559599516SKenneth E. Jansen   printPerf("writeheader", timer_start, timer_end, 0, 0, "");
176659599516SKenneth E. Jansen }
176759599516SKenneth E. Jansen 
1768103be424SCameron Smith void writeDataBlock(
1769103be424SCameron Smith     FILE* f,
1770ea868eb1SCameron Smith     const void* valueArray,
1771ea868eb1SCameron Smith     const int   nItems,
1772ea868eb1SCameron Smith     const char  datatype[],
1773ea868eb1SCameron Smith     const char  iotype[]  )
1774ea868eb1SCameron Smith {
1775ea868eb1SCameron Smith   isBinary( iotype );
1776ea868eb1SCameron Smith   size_t type_size = typeSize( datatype );
1777*f42e0444SCameron Smith   double t0 = phiotmrc();
1778ea868eb1SCameron Smith   if ( binary_format ) {
1779ea868eb1SCameron Smith     fwrite( valueArray, type_size, nItems, f );
1780ea868eb1SCameron Smith     fprintf( f,"\n");
1781ea868eb1SCameron Smith   } else {
1782ea868eb1SCameron Smith     char* ts1 = StringStripper( datatype );
1783ea868eb1SCameron Smith     if ( cscompare( "integer", ts1 ) ) {
1784be3da47bSCameron Smith       const int* vals = (int*) valueArray;
1785ea868eb1SCameron Smith       for( int n=0; n < nItems ; n++ )
1786be3da47bSCameron Smith         fprintf(f,"%d\n",vals[n]);
1787ea868eb1SCameron Smith     } else if ( cscompare( "double", ts1 ) ) {
1788be3da47bSCameron Smith       const double* vals = (double*) valueArray;
1789ea868eb1SCameron Smith       for( int n=0; n < nItems ; n++ )
1790be3da47bSCameron Smith         fprintf(f,"%f\n",vals[n]);
1791ea868eb1SCameron Smith     }
1792ea868eb1SCameron Smith     free (ts1);
1793ea868eb1SCameron Smith   }
1794*f42e0444SCameron Smith   phastaio_addWriteTime(phiotmrc()-t0);
1795*f42e0444SCameron Smith   phastaio_addWriteBytes(nItems*type_size);
1796ea868eb1SCameron Smith }
1797ea868eb1SCameron Smith 
1798103be424SCameron Smith void writedatablock(
1799103be424SCameron Smith     const int* fileDescriptor,
180059599516SKenneth E. Jansen     const char keyphrase[],
180159599516SKenneth E. Jansen     const void* valueArray,
180259599516SKenneth E. Jansen     const int* nItems,
180359599516SKenneth E. Jansen     const char datatype[],
180459599516SKenneth E. Jansen     const char iotype[] )
180559599516SKenneth E. Jansen {
180659599516SKenneth E. Jansen   //if(irank == 0) printf("entering writedatablock()\n");
180759599516SKenneth E. Jansen 
18082dd307a1SCameron Smith   unsigned long data_size = 0;
180959599516SKenneth E. Jansen   double timer_start, timer_end;
181059599516SKenneth E. Jansen   startTimer(&timer_start);
181159599516SKenneth E. Jansen 
181259599516SKenneth E. Jansen   int i = *fileDescriptor;
181359599516SKenneth E. Jansen   checkFileDescriptor("writedatablock",&i);
181459599516SKenneth E. Jansen 
181559599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
181659599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
181759599516SKenneth E. Jansen 
181859599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
181959599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
182059599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
182159599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
182259599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
182359599516SKenneth E. Jansen       endTimer(&timer_end);
182459599516SKenneth E. Jansen       printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
182559599516SKenneth E. Jansen       return;
182659599516SKenneth E. Jansen     }
182759599516SKenneth E. Jansen     // since we require that a consistant header always preceed the data block
182859599516SKenneth E. Jansen     // let us check to see that it is actually the case.
182959599516SKenneth E. Jansen 
1830961a4ff6SCameron Smith     if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
183159599516SKenneth E. Jansen       fprintf(stderr, "Header not consistant with data block\n");
1832961a4ff6SCameron Smith       fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
183359599516SKenneth E. Jansen       fprintf(stderr, "DataBlock: %s\n ", keyphrase );
183459599516SKenneth E. Jansen       fprintf(stderr, "Please recheck write sequence \n");
183559599516SKenneth E. Jansen       if( Strict_Error ) {
183659599516SKenneth E. Jansen         fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
183759599516SKenneth E. Jansen         endTimer(&timer_end);
183859599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
183959599516SKenneth E. Jansen         return;
184059599516SKenneth E. Jansen       }
184159599516SKenneth E. Jansen     }
184259599516SKenneth E. Jansen 
184359599516SKenneth E. Jansen     FILE* fileObject =  fileArray[ filePtr ] ;
184459599516SKenneth E. Jansen     size_t type_size=typeSize( datatype );
184559599516SKenneth E. Jansen     isBinary( iotype );
184659599516SKenneth E. Jansen 
1847bae968b9SCameron Smith     LastHeaderKey.erase(filePtr);
1848bae968b9SCameron Smith 
184959599516SKenneth E. Jansen     if ( header_type[filePtr] != (int)type_size ) {
185059599516SKenneth E. Jansen       fprintf(stderr,"header and datablock differ on typeof data in the block for\n");
185159599516SKenneth E. Jansen       fprintf(stderr,"keyphrase : %s\n", keyphrase);
185259599516SKenneth E. Jansen       if( Strict_Error ) {
185359599516SKenneth E. Jansen         fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
185459599516SKenneth E. Jansen         endTimer(&timer_end);
185559599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
185659599516SKenneth E. Jansen         return;
185759599516SKenneth E. Jansen       }
185859599516SKenneth E. Jansen     }
185959599516SKenneth E. Jansen 
186059599516SKenneth E. Jansen     int nUnits = *nItems;
186159599516SKenneth E. Jansen 
186259599516SKenneth E. Jansen     if ( nUnits != DataSize ) {
186359599516SKenneth E. Jansen       fprintf(stderr,"header and datablock differ on number of data items for\n");
186459599516SKenneth E. Jansen       fprintf(stderr,"keyphrase : %s\n", keyphrase);
186559599516SKenneth E. Jansen       if( Strict_Error ) {
186659599516SKenneth E. Jansen         fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
186759599516SKenneth E. Jansen         endTimer(&timer_end);
186859599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
186959599516SKenneth E. Jansen         return;
187059599516SKenneth E. Jansen       }
187159599516SKenneth E. Jansen     }
1872ea868eb1SCameron Smith     writeDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
187359599516SKenneth E. Jansen   }
187459599516SKenneth E. Jansen   else {  // syncIO case
187559599516SKenneth E. Jansen     MPI_Status write_data_status;
187659599516SKenneth E. Jansen     isBinary( iotype );
187759599516SKenneth E. Jansen     int nUnits = *nItems;
187859599516SKenneth E. Jansen 
187959599516SKenneth E. Jansen     //MR CHANGE
188059599516SKenneth E. Jansen     // 	if ( cscompare(datatype,"double") )
188159599516SKenneth E. Jansen     char* ts1 = StringStripper( datatype );
188259599516SKenneth E. Jansen     if ( cscompare("double",ts1) )
188359599516SKenneth E. Jansen       //MR CHANGE END
188459599516SKenneth E. Jansen     {
188559599516SKenneth E. Jansen       memcpy((PhastaIOActiveFiles[i]->double_chunk+DB_HEADER_SIZE/sizeof(double)), valueArray, nUnits*sizeof(double));
1886*f42e0444SCameron Smith       double t0 = phiotmrc();
188759599516SKenneth E. Jansen       MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
188859599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset,
188959599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->double_chunk,
189059599516SKenneth E. Jansen           //BLOCK_SIZE/sizeof(double),
189159599516SKenneth E. Jansen           nUnits+DB_HEADER_SIZE/sizeof(double),
189259599516SKenneth E. Jansen           MPI_DOUBLE );
189359599516SKenneth E. Jansen       MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
189459599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->double_chunk,
189559599516SKenneth E. Jansen           &write_data_status );
189659599516SKenneth E. Jansen       data_size=8*nUnits;
1897*f42e0444SCameron Smith       phastaio_addWriteTime(phiotmrc()-t0);
1898*f42e0444SCameron Smith       phastaio_addWriteBytes((nUnits*sizeof(double))+DB_HEADER_SIZE);
189959599516SKenneth E. Jansen     }
190059599516SKenneth E. Jansen     //MR CHANGE
190159599516SKenneth E. Jansen     // 	else if ( cscompare ( datatype, "integer"))
190259599516SKenneth E. Jansen     else if ( cscompare("integer",ts1) )
190359599516SKenneth E. Jansen       //MR CHANGE END
190459599516SKenneth E. Jansen     {
190559599516SKenneth E. Jansen       memcpy((PhastaIOActiveFiles[i]->int_chunk+DB_HEADER_SIZE/sizeof(int)), valueArray, nUnits*sizeof(int));
1906*f42e0444SCameron Smith       double t0 = phiotmrc();
190759599516SKenneth E. Jansen       MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
190859599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset,
190959599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->int_chunk,
191059599516SKenneth E. Jansen           nUnits+DB_HEADER_SIZE/sizeof(int),
191159599516SKenneth E. Jansen           MPI_INT );
191259599516SKenneth E. Jansen       MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
191359599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->int_chunk,
191459599516SKenneth E. Jansen           &write_data_status );
191559599516SKenneth E. Jansen       data_size=4*nUnits;
1916*f42e0444SCameron Smith       phastaio_addWriteTime(phiotmrc()-t0);
1917*f42e0444SCameron Smith       phastaio_addWriteBytes((nUnits*sizeof(int))+DB_HEADER_SIZE);
191859599516SKenneth E. Jansen     }
191959599516SKenneth E. Jansen     else {
192059599516SKenneth E. Jansen       printf("Error: writedatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
192159599516SKenneth E. Jansen       endTimer(&timer_end);
192259599516SKenneth E. Jansen       printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
192359599516SKenneth E. Jansen       return;
192459599516SKenneth E. Jansen     }
192559599516SKenneth E. Jansen     free(ts1);
192659599516SKenneth E. Jansen   }
192759599516SKenneth E. Jansen 
192859599516SKenneth E. Jansen   endTimer(&timer_end);
192959599516SKenneth E. Jansen   char extra_msg[1024];
193059599516SKenneth E. Jansen   memset(extra_msg, '\0', 1024);
193159599516SKenneth E. Jansen   char* key = StringStripper(keyphrase);
193259599516SKenneth E. Jansen   sprintf(extra_msg, " field is %s ", key);
193359599516SKenneth E. Jansen   printPerf("writedatablock", timer_start, timer_end, data_size, 1, extra_msg);
193459599516SKenneth E. Jansen   free(key);
193559599516SKenneth E. Jansen 
193659599516SKenneth E. Jansen }
193759599516SKenneth E. Jansen 
1938103be424SCameron Smith void SwapArrayByteOrder( void* array, int nbytes, int nItems )
193959599516SKenneth E. Jansen {
194059599516SKenneth E. Jansen   /* This swaps the byte order for the array of nItems each
194159599516SKenneth E. Jansen      of size nbytes , This will be called only locally  */
194259599516SKenneth E. Jansen   int i,j;
194359599516SKenneth E. Jansen   unsigned char* ucDst = (unsigned char*)array;
194459599516SKenneth E. Jansen   for(i=0; i < nItems; i++) {
194559599516SKenneth E. Jansen     for(j=0; j < (nbytes/2); j++)
194659599516SKenneth E. Jansen       std::swap( ucDst[j] , ucDst[(nbytes - 1) - j] );
194759599516SKenneth E. Jansen     ucDst += nbytes;
194859599516SKenneth E. Jansen   }
194959599516SKenneth E. Jansen }
195059599516SKenneth E. Jansen 
1951103be424SCameron Smith void writestring( int* fileDescriptor, const char inString[] )
195259599516SKenneth E. Jansen {
195359599516SKenneth E. Jansen   int filePtr = *fileDescriptor - 1;
195459599516SKenneth E. Jansen   FILE* fileObject = fileArray[filePtr] ;
195559599516SKenneth E. Jansen   fprintf(fileObject,"%s",inString );
195659599516SKenneth E. Jansen   return;
195759599516SKenneth E. Jansen }
195859599516SKenneth E. Jansen 
1959103be424SCameron Smith void Gather_Headers( int* fileDescriptor, vector< string >& headers )
196059599516SKenneth E. Jansen {
196159599516SKenneth E. Jansen   FILE* fileObject;
196259599516SKenneth E. Jansen   char Line[1024];
196359599516SKenneth E. Jansen 
196459599516SKenneth E. Jansen   fileObject = fileArray[ (*fileDescriptor)-1 ];
196559599516SKenneth E. Jansen 
196659599516SKenneth E. Jansen   while( !feof(fileObject) ) {
196759599516SKenneth E. Jansen     fgets( Line, 1024, fileObject);
196859599516SKenneth E. Jansen     if ( Line[0] == '#' ) {
196959599516SKenneth E. Jansen       headers.push_back( Line );
197059599516SKenneth E. Jansen     } else {
197159599516SKenneth E. Jansen       break;
197259599516SKenneth E. Jansen     }
197359599516SKenneth E. Jansen   }
197459599516SKenneth E. Jansen   rewind( fileObject );
197559599516SKenneth E. Jansen   clearerr( fileObject );
197659599516SKenneth E. Jansen }
197759599516SKenneth E. Jansen 
1978103be424SCameron Smith void isWrong( void ) {
1979103be424SCameron Smith   (Wrong_Endian) ? fprintf(stdout,"YES\n") : fprintf(stdout,"NO\n");
1980103be424SCameron Smith }
198159599516SKenneth E. Jansen 
1982103be424SCameron Smith void togglestrictmode( void ) { Strict_Error = !Strict_Error; }
1983103be424SCameron Smith 
1984103be424SCameron Smith int isLittleEndian( void )
198559599516SKenneth E. Jansen {
198659599516SKenneth E. Jansen   // this function returns a 1 if the current running architecture is
198759599516SKenneth E. Jansen   // LittleEndian Byte Ordered, else it returns a zero
198859599516SKenneth E. Jansen   union  {
198959599516SKenneth E. Jansen     long a;
199059599516SKenneth E. Jansen     char c[sizeof( long )];
199159599516SKenneth E. Jansen   } endianUnion;
199259599516SKenneth E. Jansen   endianUnion.a = 1 ;
199359599516SKenneth E. Jansen   if ( endianUnion.c[sizeof(long)-1] != 1 ) return 1 ;
199459599516SKenneth E. Jansen   else return 0;
199559599516SKenneth E. Jansen }
199659599516SKenneth E. Jansen 
199759599516SKenneth E. Jansen namespace PHASTA {
199859599516SKenneth E. Jansen   const char* const PhastaIO_traits<int>::type_string = "integer";
199959599516SKenneth E. Jansen   const char* const PhastaIO_traits<double>::type_string = "double";
200059599516SKenneth E. Jansen }
2001