xref: /phasta/phastaIO/phastaIO.cc (revision be3da47be76b4ea61b64417cc61d4240f594ed07)
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>
1759599516SKenneth E. Jansen #include <string.h>
1859599516SKenneth E. Jansen #include <ctype.h>
1959599516SKenneth E. Jansen #include <stdlib.h>
2059599516SKenneth E. Jansen #include <stdio.h>
2159599516SKenneth E. Jansen #include <math.h>
22d3337298SCameron Smith #include <sstream>
2359599516SKenneth E. Jansen #include "phastaIO.h"
2459599516SKenneth E. Jansen #include "mpi.h"
2559599516SKenneth E. Jansen #include "phiotmrc.h"
2697a07b0aSCameron Smith #include <assert.h>
2759599516SKenneth E. Jansen 
2859599516SKenneth E. Jansen #define VERSION_INFO_HEADER_SIZE 8192
2959599516SKenneth E. Jansen #define DB_HEADER_SIZE 1024
3059599516SKenneth E. Jansen #define ONE_MEGABYTE 1048576
3159599516SKenneth E. Jansen #define TWO_MEGABYTE 2097152
3259599516SKenneth E. Jansen #define ENDIAN_TEST_NUMBER 12180 // Troy's Zip Code!!
3359599516SKenneth E. Jansen #define MAX_PHASTA_FILES 64
3459599516SKenneth E. Jansen #define MAX_PHASTA_FILE_NAME_LENGTH 1024
3559599516SKenneth 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
3659599516SKenneth E. Jansen // -3 for MPI_IO_Tag, nFields and nppf, -4 for extra security (former nFiles)
3759599516SKenneth E. Jansen #define MAX_FIELDS_NAME_LENGTH 128
3859599516SKenneth E. Jansen #define DefaultMHSize (4*ONE_MEGABYTE)
3959599516SKenneth E. Jansen //#define DefaultMHSize (8350) //For test
4059599516SKenneth E. Jansen #define LATEST_WRITE_VERSION 1
4159599516SKenneth E. Jansen #define inv1024sq 953.674316406e-9 // = 1/1024/1024
4259599516SKenneth E. Jansen int MasterHeaderSize = -1;
4359599516SKenneth E. Jansen 
4459599516SKenneth E. Jansen bool PRINT_PERF = false; // default print no perf results
4559599516SKenneth E. Jansen int irank = -1; // global rank, should never be manually manipulated
4659599516SKenneth E. Jansen int mysize = -1;
4759599516SKenneth E. Jansen 
4859599516SKenneth E. Jansen // Static variables are bad but used here to store the subcommunicator and associated variables
4959599516SKenneth 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)
5059599516SKenneth E. Jansen static int s_assign_local_comm = 0;
5159599516SKenneth E. Jansen static MPI_Comm s_local_comm;
5259599516SKenneth E. Jansen static int s_local_size = -1;
5359599516SKenneth E. Jansen static int s_local_rank = -1;
5459599516SKenneth E. Jansen 
5559599516SKenneth E. Jansen // the following defines are for debug printf
5659599516SKenneth E. Jansen #define PHASTAIO_DEBUG 0 //default to not print any debugging info
5759599516SKenneth E. Jansen 
58fe03b87fSCameron Smith void phprintf(const char* fmt, ...) {
59*be3da47bSCameron Smith   (void)fmt;
6059599516SKenneth E. Jansen #if PHASTAIO_DEBUG
61fe03b87fSCameron Smith   char format[1024];
62fe03b87fSCameron Smith   snprintf(format, sizeof(format), "phastaIO debug: %s", fmt);
63fe03b87fSCameron Smith   va_list ap;
64fe03b87fSCameron Smith   va_start(ap,fmt);
65fe03b87fSCameron Smith   vprintf(format,ap)
66fe03b87fSCameron Smith   va_end(ap);
6759599516SKenneth E. Jansen #endif
68fe03b87fSCameron Smith }
69fe03b87fSCameron Smith 
70fe03b87fSCameron Smith void phprintf_0(const char* fmt, ...) {
71*be3da47bSCameron Smith   (void)fmt;
72fe03b87fSCameron Smith #if PHASTAIO_DEBUG
73fe03b87fSCameron Smith   int rank = 0;
74fe03b87fSCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
75fe03b87fSCameron Smith   if(rank == 0){
76fe03b87fSCameron Smith     char format[1024];
77fe03b87fSCameron Smith     snprintf(format, sizeof(format), "phastaIO debug: irank=0 %s", fmt);
78fe03b87fSCameron Smith     va_list ap;
79fe03b87fSCameron Smith     va_start(ap,s);
80fe03b87fSCameron Smith     vprintf(format, ap);
81fe03b87fSCameron Smith     va_end(ap);
82fe03b87fSCameron Smith   }
83fe03b87fSCameron Smith #endif
84fe03b87fSCameron Smith }
8559599516SKenneth E. Jansen 
8659599516SKenneth E. Jansen enum PhastaIO_Errors
8759599516SKenneth E. Jansen {
8859599516SKenneth E. Jansen   MAX_PHASTA_FILES_EXCEEDED = -1,
8959599516SKenneth E. Jansen   UNABLE_TO_OPEN_FILE = -2,
9059599516SKenneth E. Jansen   NOT_A_MPI_FILE = -3,
9159599516SKenneth E. Jansen   GPID_EXCEEDED = -4,
9210d56689SCameron Smith   DATA_TYPE_ILLEGAL = -5
9359599516SKenneth E. Jansen };
9459599516SKenneth E. Jansen 
9559599516SKenneth E. Jansen using namespace std;
9659599516SKenneth E. Jansen 
9759599516SKenneth E. Jansen namespace{
9859599516SKenneth E. Jansen 
99961a4ff6SCameron Smith   map<int, std::string> LastHeaderKey;
10059599516SKenneth E. Jansen   vector< FILE* > fileArray;
10159599516SKenneth E. Jansen   vector< bool > byte_order;
10259599516SKenneth E. Jansen   vector< int > header_type;
10359599516SKenneth E. Jansen   int DataSize=0;
10459599516SKenneth E. Jansen   bool LastHeaderNotFound = false;
10559599516SKenneth E. Jansen   bool Wrong_Endian = false ;
10659599516SKenneth E. Jansen   bool Strict_Error = false ;
10759599516SKenneth E. Jansen   bool binary_format = true;
10859599516SKenneth E. Jansen 
10959599516SKenneth E. Jansen   /***********************************************************************/
11059599516SKenneth E. Jansen   /***************** NEW PHASTA IO CODE STARTS HERE **********************/
11159599516SKenneth E. Jansen   /***********************************************************************/
11259599516SKenneth E. Jansen 
11359599516SKenneth E. Jansen   typedef struct
11459599516SKenneth E. Jansen   {
11559599516SKenneth E. Jansen     char filename[MAX_PHASTA_FILE_NAME_LENGTH];   /* defafults to 1024 */
1162dd307a1SCameron Smith     unsigned long my_offset;
1172dd307a1SCameron Smith     unsigned long next_start_address;
1182dd307a1SCameron Smith     unsigned long **my_offset_table;
1192dd307a1SCameron Smith     unsigned long **my_read_table;
12059599516SKenneth E. Jansen 
12159599516SKenneth E. Jansen     double * double_chunk;
12259599516SKenneth E. Jansen     double * read_double_chunk;
12359599516SKenneth E. Jansen 
12459599516SKenneth E. Jansen     int field_count;
12559599516SKenneth E. Jansen     int part_count;
12659599516SKenneth E. Jansen     int read_field_count;
12759599516SKenneth E. Jansen     int read_part_count;
12859599516SKenneth E. Jansen     int GPid;
12959599516SKenneth E. Jansen     int start_id;
13059599516SKenneth E. Jansen 
13159599516SKenneth E. Jansen     int mhsize;
13259599516SKenneth E. Jansen 
13359599516SKenneth E. Jansen     int myrank;
13459599516SKenneth E. Jansen     int numprocs;
13559599516SKenneth E. Jansen     int local_myrank;
13659599516SKenneth E. Jansen     int local_numprocs;
13759599516SKenneth E. Jansen 
13859599516SKenneth E. Jansen     int nppp;
13959599516SKenneth E. Jansen     int nPPF;
14059599516SKenneth E. Jansen     int nFiles;
14159599516SKenneth E. Jansen     int nFields;
14259599516SKenneth E. Jansen 
14359599516SKenneth E. Jansen     int * int_chunk;
14459599516SKenneth E. Jansen     int * read_int_chunk;
14559599516SKenneth E. Jansen 
14659599516SKenneth E. Jansen     int Wrong_Endian; /* default to false */
14759599516SKenneth E. Jansen     char * master_header;
14859599516SKenneth E. Jansen     MPI_File file_handle;
14959599516SKenneth E. Jansen     MPI_Comm local_comm;
15059599516SKenneth E. Jansen   } phastaio_file_t;
15159599516SKenneth E. Jansen 
15259599516SKenneth E. Jansen   typedef struct
15359599516SKenneth E. Jansen   {
15459599516SKenneth E. Jansen     int nppf, nfields;
15559599516SKenneth E. Jansen     char * masterHeader;
15659599516SKenneth E. Jansen   }serial_file;
15759599516SKenneth E. Jansen 
15859599516SKenneth E. Jansen   serial_file *SerialFile;
15959599516SKenneth E. Jansen   phastaio_file_t *PhastaIOActiveFiles[MAX_PHASTA_FILES];
16059599516SKenneth E. Jansen   int PhastaIONextActiveIndex = 0; /* indicates next index to allocate */
16159599516SKenneth E. Jansen 
16259599516SKenneth E. Jansen   // the caller has the responsibility to delete the returned string
16359599516SKenneth E. Jansen   // TODO: StringStipper("nbc value? ") returns NULL?
164103be424SCameron Smith   char* StringStripper( const char  istring[] ) {
16559599516SKenneth E. Jansen     int length = strlen( istring );
16659599516SKenneth E. Jansen     char* dest = (char *)malloc( length + 1 );
16759599516SKenneth E. Jansen     strcpy( dest, istring );
16859599516SKenneth E. Jansen     dest[ length ] = '\0';
16959599516SKenneth E. Jansen     if ( char* p = strpbrk( dest, " ") )
17059599516SKenneth E. Jansen       *p = '\0';
17159599516SKenneth E. Jansen     return dest;
17259599516SKenneth E. Jansen   }
17359599516SKenneth E. Jansen 
174103be424SCameron Smith   inline int cscompare( const char teststring[], const char targetstring[] ) {
17559599516SKenneth E. Jansen     char* s1 = const_cast<char*>(teststring);
17659599516SKenneth E. Jansen     char* s2 = const_cast<char*>(targetstring);
17759599516SKenneth E. Jansen 
17859599516SKenneth E. Jansen     while( *s1 == ' ') s1++;
17959599516SKenneth E. Jansen     while( *s2 == ' ') s2++;
18059599516SKenneth E. Jansen     while( ( *s1 )
18159599516SKenneth E. Jansen         && ( *s2 )
18259599516SKenneth E. Jansen         && ( *s2 != '?')
18359599516SKenneth E. Jansen         && ( tolower( *s1 )==tolower( *s2 ) ) ) {
18459599516SKenneth E. Jansen       s1++;
18559599516SKenneth E. Jansen       s2++;
18659599516SKenneth E. Jansen       while( *s1 == ' ') s1++;
18759599516SKenneth E. Jansen       while( *s2 == ' ') s2++;
18859599516SKenneth E. Jansen     }
18959599516SKenneth E. Jansen     if ( !( *s1 ) || ( *s1 == '?') ) return 1;
19059599516SKenneth E. Jansen     else return 0;
19159599516SKenneth E. Jansen   }
19259599516SKenneth E. Jansen 
193103be424SCameron Smith   inline void isBinary( const char iotype[] ) {
19459599516SKenneth E. Jansen     char* fname = StringStripper( iotype );
19559599516SKenneth E. Jansen     if ( cscompare( fname, "binary" ) ) binary_format = true;
19659599516SKenneth E. Jansen     else binary_format = false;
19759599516SKenneth E. Jansen     free (fname);
19859599516SKenneth E. Jansen 
19959599516SKenneth E. Jansen   }
20059599516SKenneth E. Jansen 
201103be424SCameron Smith   inline size_t typeSize( const char typestring[] ) {
20259599516SKenneth E. Jansen     char* ts1 = StringStripper( typestring );
20359599516SKenneth E. Jansen     if ( cscompare( "integer", ts1 ) ) {
20459599516SKenneth E. Jansen       free (ts1);
20559599516SKenneth E. Jansen       return sizeof(int);
20659599516SKenneth E. Jansen     } else if ( cscompare( "double", ts1 ) ) {
20759599516SKenneth E. Jansen       free (ts1);
20859599516SKenneth E. Jansen       return sizeof( double );
20959599516SKenneth E. Jansen     } else {
21059599516SKenneth E. Jansen       free (ts1);
21159599516SKenneth E. Jansen       fprintf(stderr,"unknown type : %s\n",ts1);
21259599516SKenneth E. Jansen       return 0;
21359599516SKenneth E. Jansen     }
21459599516SKenneth E. Jansen   }
21559599516SKenneth E. Jansen 
216103be424SCameron Smith   int readHeader(
217103be424SCameron Smith       FILE*       fileObject,
21859599516SKenneth E. Jansen         const char  phrase[],
21959599516SKenneth E. Jansen         int*        params,
22059599516SKenneth E. Jansen         int         expect ) {
22159599516SKenneth E. Jansen     char* text_header;
22259599516SKenneth E. Jansen     char* token;
22397a07b0aSCameron Smith     char Line[1024] = "\0";
22459599516SKenneth E. Jansen     char junk;
22559599516SKenneth E. Jansen     bool FOUND = false ;
22659599516SKenneth E. Jansen     int real_length;
22759599516SKenneth E. Jansen     int skip_size, integer_value;
22859599516SKenneth E. Jansen     int rewind_count=0;
22959599516SKenneth E. Jansen 
23059599516SKenneth E. Jansen     if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
23159599516SKenneth E. Jansen       rewind( fileObject );
23259599516SKenneth E. Jansen       clearerr( fileObject );
23359599516SKenneth E. Jansen       rewind_count++;
23459599516SKenneth E. Jansen       fgets( Line, 1024, fileObject );
23559599516SKenneth E. Jansen     }
23659599516SKenneth E. Jansen 
23759599516SKenneth E. Jansen     while( !FOUND  && ( rewind_count < 2 ) )  {
23859599516SKenneth E. Jansen       if ( ( Line[0] != '\n' ) && ( real_length = strcspn( Line, "#" )) ) {
23959599516SKenneth E. Jansen         text_header = (char *)malloc( real_length + 1 );
24059599516SKenneth E. Jansen         strncpy( text_header, Line, real_length );
24159599516SKenneth E. Jansen         text_header[ real_length ] =static_cast<char>(NULL);
24259599516SKenneth E. Jansen         token = strtok ( text_header, ":" );
24397a07b0aSCameron Smith         assert(token);
24497a07b0aSCameron Smith         if( cscompare( phrase , token ) ) {
24559599516SKenneth E. Jansen           FOUND = true ;
24659599516SKenneth E. Jansen           token = strtok( NULL, " ,;<>" );
24797a07b0aSCameron Smith           assert(token);
24859599516SKenneth E. Jansen           skip_size = atoi( token );
24959599516SKenneth E. Jansen           int i;
25059599516SKenneth E. Jansen           for( i=0; i < expect && ( token = strtok( NULL," ,;<>") ); i++) {
25197a07b0aSCameron Smith             assert(token);
25259599516SKenneth E. Jansen             params[i] = atoi( token );
25359599516SKenneth E. Jansen           }
25459599516SKenneth E. Jansen           if ( i < expect ) {
255103be424SCameron Smith             fprintf(stderr,
256103be424SCameron Smith                 "Aloha Expected # of ints not found for: %s\n",phrase );
25759599516SKenneth E. Jansen           }
25859599516SKenneth E. Jansen         } else if ( cscompare(token,"byteorder magic number") ) {
25959599516SKenneth E. Jansen           if ( binary_format ) {
26059599516SKenneth E. Jansen             fread((void*)&integer_value,sizeof(int),1,fileObject);
26159599516SKenneth E. Jansen             fread( &junk, sizeof(char), 1 , fileObject );
26259599516SKenneth E. Jansen             if ( 362436 != integer_value ) Wrong_Endian = true;
26359599516SKenneth E. Jansen           } else{
26459599516SKenneth E. Jansen             fscanf(fileObject, "%d\n", &integer_value );
26559599516SKenneth E. Jansen           }
26659599516SKenneth E. Jansen         } else {
26759599516SKenneth E. Jansen           /* some other header, so just skip over */
26859599516SKenneth E. Jansen           token = strtok( NULL, " ,;<>" );
26997a07b0aSCameron Smith           assert(token);
27059599516SKenneth E. Jansen           skip_size = atoi( token );
27159599516SKenneth E. Jansen           if ( binary_format)
27259599516SKenneth E. Jansen             fseek( fileObject, skip_size, SEEK_CUR );
27359599516SKenneth E. Jansen           else
27459599516SKenneth E. Jansen             for( int gama=0; gama < skip_size; gama++ )
27559599516SKenneth E. Jansen               fgets( Line, 1024, fileObject );
27659599516SKenneth E. Jansen         }
27759599516SKenneth E. Jansen         free (text_header);
27859599516SKenneth E. Jansen       } // end of if before while loop
27959599516SKenneth E. Jansen 
28059599516SKenneth E. Jansen       if ( !FOUND )
28159599516SKenneth E. Jansen         if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
28259599516SKenneth E. Jansen           rewind( fileObject );
28359599516SKenneth E. Jansen           clearerr( fileObject );
28459599516SKenneth E. Jansen           rewind_count++;
28559599516SKenneth E. Jansen           fgets( Line, 1024, fileObject );
28659599516SKenneth E. Jansen         }
28759599516SKenneth E. Jansen     }
28859599516SKenneth E. Jansen 
28959599516SKenneth E. Jansen     if ( !FOUND ) {
29059599516SKenneth E. Jansen       //fprintf(stderr, "Error: Could not find: %s\n", phrase);
29159599516SKenneth E. Jansen       if(irank == 0) printf("WARNING: Could not find: %s\n", phrase);
29259599516SKenneth E. Jansen       return 1;
29359599516SKenneth E. Jansen     }
29459599516SKenneth E. Jansen     return 0;
29559599516SKenneth E. Jansen   }
29659599516SKenneth E. Jansen 
29759599516SKenneth E. Jansen } // end unnamed namespace
29859599516SKenneth E. Jansen 
29959599516SKenneth E. Jansen 
30059599516SKenneth E. Jansen // begin of publicly visible functions
30159599516SKenneth E. Jansen 
30259599516SKenneth E. Jansen /**
30359599516SKenneth E. Jansen  * This function takes a long long pointer and assign (start) phiotmrc value to it
30459599516SKenneth E. Jansen  */
30559599516SKenneth E. Jansen void startTimer(double* start) {
30659599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
30759599516SKenneth E. Jansen   MPI_Barrier(MPI_COMM_WORLD);
30859599516SKenneth E. Jansen   *start =  phiotmrc();
30959599516SKenneth E. Jansen }
31059599516SKenneth E. Jansen 
31159599516SKenneth E. Jansen /**
31259599516SKenneth E. Jansen  * This function takes a long long pointer and assign (end) phiotmrc value to it
31359599516SKenneth E. Jansen  */
31459599516SKenneth E. Jansen void endTimer(double* end) {
31559599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
31659599516SKenneth E. Jansen   *end = phiotmrc();
31759599516SKenneth E. Jansen   MPI_Barrier(MPI_COMM_WORLD);
31859599516SKenneth E. Jansen }
31959599516SKenneth E. Jansen 
32059599516SKenneth E. Jansen /**
32159599516SKenneth E. Jansen  * choose to print some performance results (or not) according to
32259599516SKenneth E. Jansen  * the PRINT_PERF macro
32359599516SKenneth E. Jansen  */
32459599516SKenneth E. Jansen void printPerf(
32559599516SKenneth E. Jansen     const char* func_name,
32659599516SKenneth E. Jansen     double start,
32759599516SKenneth E. Jansen     double end,
3282dd307a1SCameron Smith     unsigned long datasize,
32959599516SKenneth E. Jansen     int printdatainfo,
33059599516SKenneth E. Jansen     const char* extra_msg) {
33159599516SKenneth E. Jansen   if( !PRINT_PERF ) return;
3322dd307a1SCameron Smith   unsigned long data_size = datasize;
33359599516SKenneth E. Jansen   double time = end - start;
3342dd307a1SCameron Smith   unsigned long isizemin,isizemax,isizetot;
33559599516SKenneth E. Jansen   double sizemin,sizemax,sizeavg,sizetot,rate;
33659599516SKenneth E. Jansen   double tmin, tmax, tavg, ttot;
33759599516SKenneth E. Jansen 
33859599516SKenneth E. Jansen   MPI_Allreduce(&time, &tmin,1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
33959599516SKenneth E. Jansen   MPI_Allreduce(&time, &tmax,1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
34059599516SKenneth E. Jansen   MPI_Allreduce(&time, &ttot,1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
34159599516SKenneth E. Jansen   tavg = ttot/mysize;
34259599516SKenneth E. Jansen 
34359599516SKenneth E. Jansen   if(irank == 0) {
34459599516SKenneth E. Jansen     if ( PhastaIONextActiveIndex == 0 ) printf("** 1PFPP ");
34559599516SKenneth E. Jansen     else  printf("** syncIO ");
34659599516SKenneth E. Jansen     printf("%s(): Tmax = %f sec, Tmin = %f sec, Tavg = %f sec", func_name, tmax, tmin, tavg);
34759599516SKenneth E. Jansen   }
34859599516SKenneth E. Jansen 
34959599516SKenneth E. Jansen   if(printdatainfo == 1) { // if printdatainfo ==1, compute I/O rate and block size
35059599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizemin,1,MPI_LONG_LONG_INT,MPI_MIN,MPI_COMM_WORLD);
35159599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizemax,1,MPI_LONG_LONG_INT,MPI_MAX,MPI_COMM_WORLD);
35259599516SKenneth E. Jansen     MPI_Allreduce(&data_size,&isizetot,1,MPI_LONG_LONG_INT,MPI_SUM,MPI_COMM_WORLD);
35359599516SKenneth E. Jansen 
35459599516SKenneth E. Jansen     sizemin=(double)(isizemin*inv1024sq);
35559599516SKenneth E. Jansen     sizemax=(double)(isizemax*inv1024sq);
35659599516SKenneth E. Jansen     sizetot=(double)(isizetot*inv1024sq);
35759599516SKenneth E. Jansen     sizeavg=(double)(1.0*sizetot/mysize);
35859599516SKenneth E. Jansen     rate=(double)(1.0*sizetot/tmax);
35959599516SKenneth E. Jansen 
36059599516SKenneth E. Jansen     if( irank == 0) {
361103be424SCameron Smith       printf(", Rate = %f MB/s [%s] \n \t\t\t"
362103be424SCameron Smith              " block size: Min= %f MB; Max= %f MB; Avg= %f MB; Tot= %f MB\n",
363103be424SCameron Smith              rate, extra_msg, sizemin,sizemax,sizeavg,sizetot);
36459599516SKenneth E. Jansen     }
36559599516SKenneth E. Jansen   }
36659599516SKenneth E. Jansen   else {
36759599516SKenneth E. Jansen     if(irank == 0) {
36859599516SKenneth E. Jansen       printf(" \n");
36959599516SKenneth E. Jansen       //printf(" (%s) \n", extra_msg);
37059599516SKenneth E. Jansen     }
37159599516SKenneth E. Jansen   }
37259599516SKenneth E. Jansen }
37359599516SKenneth E. Jansen 
37459599516SKenneth E. Jansen /**
37559599516SKenneth E. Jansen  * This function is normally called at the beginning of a read operation, before
37659599516SKenneth E. Jansen  * init function.
37759599516SKenneth E. Jansen  * This function (uses rank 0) reads out nfields, nppf, master header size,
37859599516SKenneth E. Jansen  * endianess and allocates for masterHeader string.
37959599516SKenneth E. Jansen  * These values are essential for following read operations. Rank 0 will bcast
38059599516SKenneth E. Jansen  * these values to other ranks in the commm world
38159599516SKenneth E. Jansen  *
38259599516SKenneth E. Jansen  * If the file set is of old POSIX format, it would throw error and exit
38359599516SKenneth E. Jansen  */
38459599516SKenneth E. Jansen void queryphmpiio(const char filename[],int *nfields, int *nppf)
38559599516SKenneth E. Jansen {
38659599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &irank);
38759599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &mysize);
38859599516SKenneth E. Jansen 
38959599516SKenneth E. Jansen   if(irank == 0) {
39059599516SKenneth E. Jansen     FILE * fileHandle;
39159599516SKenneth E. Jansen     char* fname = StringStripper( filename );
39259599516SKenneth E. Jansen 
39359599516SKenneth E. Jansen     fileHandle = fopen (fname,"rb");
39459599516SKenneth E. Jansen     if (fileHandle == NULL ) {
39559599516SKenneth E. Jansen       printf("\nError: File %s doesn't exist! Please check!\n",fname);
39659599516SKenneth E. Jansen     }
39759599516SKenneth E. Jansen     else {
39859599516SKenneth E. Jansen       SerialFile =(serial_file *)calloc( 1,  sizeof( serial_file) );
39959599516SKenneth E. Jansen       int meta_size_limit = VERSION_INFO_HEADER_SIZE;
40059599516SKenneth E. Jansen       SerialFile->masterHeader = (char *)malloc( meta_size_limit );
40159599516SKenneth E. Jansen       fread(SerialFile->masterHeader, 1, meta_size_limit, fileHandle);
40259599516SKenneth E. Jansen 
40359599516SKenneth E. Jansen       char read_out_tag[MAX_FIELDS_NAME_LENGTH];
40459599516SKenneth E. Jansen       char version[MAX_FIELDS_NAME_LENGTH/4];
40559599516SKenneth E. Jansen       int mhsize;
40659599516SKenneth E. Jansen       char * token;
40759599516SKenneth E. Jansen       int magic_number;
40859599516SKenneth E. Jansen 
40959599516SKenneth E. Jansen       memcpy( read_out_tag,
41059599516SKenneth E. Jansen           SerialFile->masterHeader,
41159599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
41259599516SKenneth E. Jansen 
41359599516SKenneth E. Jansen       if ( cscompare ("MPI_IO_Tag",read_out_tag) ) {
41459599516SKenneth E. Jansen         // Test endianess ...
41559599516SKenneth E. Jansen         memcpy (&magic_number,
41659599516SKenneth E. Jansen             SerialFile->masterHeader + sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
41759599516SKenneth E. Jansen             sizeof(int) );                                        // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
41859599516SKenneth E. Jansen 
41959599516SKenneth E. Jansen         if ( magic_number != ENDIAN_TEST_NUMBER ) {
42059599516SKenneth E. Jansen           printf("Endian is different!\n");
42159599516SKenneth E. Jansen           // Will do swap later
42259599516SKenneth E. Jansen         }
42359599516SKenneth E. Jansen 
42459599516SKenneth E. Jansen         // test version, old version, default masterheader size is 4M
42559599516SKenneth E. Jansen         // newer version masterheader size is read from first line
42659599516SKenneth E. Jansen         memcpy(version,
42759599516SKenneth E. Jansen             SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/2,
42859599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4 - 1); //TODO: why -1?
42959599516SKenneth E. Jansen 
43059599516SKenneth E. Jansen         if( cscompare ("version",version) ) {
43159599516SKenneth E. Jansen           // if there is "version" tag in the file, then it is newer format
43259599516SKenneth E. Jansen           // read master header size from here, otherwise use default
43359599516SKenneth E. Jansen           // Note: if version is "1", we know mhsize is at 3/4 place...
43459599516SKenneth E. Jansen 
43559599516SKenneth E. Jansen           token = strtok(version, ":");
43659599516SKenneth E. Jansen           token = strtok(NULL, " ,;<>" );
43759599516SKenneth E. Jansen           int iversion = atoi(token);
43859599516SKenneth E. Jansen 
43959599516SKenneth E. Jansen           if( iversion == 1) {
44059599516SKenneth E. Jansen             memcpy( &mhsize,
44159599516SKenneth E. Jansen                 SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/4*3 + sizeof("mhsize : ")-1,
44259599516SKenneth E. Jansen                 sizeof(int));
44359599516SKenneth E. Jansen             if ( magic_number != ENDIAN_TEST_NUMBER )
44459599516SKenneth E. Jansen               SwapArrayByteOrder(&mhsize, sizeof(int), 1);
44559599516SKenneth E. Jansen 
44659599516SKenneth E. Jansen             if( mhsize > DefaultMHSize ) {
44759599516SKenneth E. Jansen               //if actual headersize is larger than default, let's re-read
44859599516SKenneth E. Jansen               free(SerialFile->masterHeader);
44959599516SKenneth E. Jansen               SerialFile->masterHeader = (char *)malloc(mhsize);
45059599516SKenneth E. Jansen               fseek(fileHandle, 0, SEEK_SET); // reset the file stream position
45159599516SKenneth E. Jansen               fread(SerialFile->masterHeader,1,mhsize,fileHandle);
45259599516SKenneth E. Jansen             }
45359599516SKenneth E. Jansen           }
45459599516SKenneth E. Jansen           //TODO: check if this is a valid int??
45559599516SKenneth E. Jansen           MasterHeaderSize = mhsize;
45659599516SKenneth E. Jansen         }
45759599516SKenneth E. Jansen         else { // else it's version 0's format w/o version tag, implicating MHSize=4M
45859599516SKenneth E. Jansen           MasterHeaderSize = DefaultMHSize;
45959599516SKenneth E. Jansen         }
46059599516SKenneth E. Jansen 
46159599516SKenneth E. Jansen         memcpy( read_out_tag,
46259599516SKenneth E. Jansen             SerialFile->masterHeader+MAX_FIELDS_NAME_LENGTH+1,
46359599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH ); //TODO: why +1
46459599516SKenneth E. Jansen 
46559599516SKenneth E. Jansen         // Read in # fields ...
46659599516SKenneth E. Jansen         token = strtok( read_out_tag, ":" );
46759599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
46859599516SKenneth E. Jansen         *nfields = atoi( token );
46959599516SKenneth E. Jansen         if ( *nfields > MAX_FIELDS_NUMBER) {
47059599516SKenneth E. Jansen           printf("Error queryphmpiio: nfields is larger than MAX_FIELDS_NUMBER!\n");
47159599516SKenneth E. Jansen         }
47259599516SKenneth E. Jansen         SerialFile->nfields=*nfields; //TODO: sanity check of this int?
47359599516SKenneth E. Jansen 
47459599516SKenneth E. Jansen         memcpy( read_out_tag,
47559599516SKenneth E. Jansen             SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH * 2
47659599516SKenneth E. Jansen             + *nfields * MAX_FIELDS_NAME_LENGTH,
47759599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
47859599516SKenneth E. Jansen 
47959599516SKenneth E. Jansen         token = strtok( read_out_tag, ":" );
48059599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
48159599516SKenneth E. Jansen         *nppf = atoi( token );
48259599516SKenneth E. Jansen         SerialFile->nppf=*nppf; //TODO: sanity check of int
48359599516SKenneth E. Jansen       } // end of if("MPI_IO_TAG")
48459599516SKenneth E. Jansen       else {
48559599516SKenneth 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);
48659599516SKenneth E. Jansen         exit(1);
48759599516SKenneth E. Jansen       }
48859599516SKenneth E. Jansen       fclose(fileHandle);
48959599516SKenneth E. Jansen       free(SerialFile->masterHeader);
49059599516SKenneth E. Jansen       free(SerialFile);
49159599516SKenneth E. Jansen     } //end of else
49259599516SKenneth E. Jansen     free(fname);
49359599516SKenneth E. Jansen   }
49459599516SKenneth E. Jansen 
49559599516SKenneth E. Jansen   // Bcast value to every one
49659599516SKenneth E. Jansen   MPI_Bcast( nfields, 1, MPI_INT, 0, MPI_COMM_WORLD );
49759599516SKenneth E. Jansen   MPI_Bcast( nppf, 1, MPI_INT, 0, MPI_COMM_WORLD );
49859599516SKenneth E. Jansen   MPI_Bcast( &MasterHeaderSize, 1, MPI_INT, 0, MPI_COMM_WORLD );
49959599516SKenneth E. Jansen   phprintf("Info queryphmpiio: myrank = %d, MasterHeaderSize = %d", irank, MasterHeaderSize);
50059599516SKenneth E. Jansen }
50159599516SKenneth E. Jansen 
50259599516SKenneth E. Jansen /**
50359599516SKenneth E. Jansen  * This function computes the right master header size (round to size of 2^n).
50459599516SKenneth E. Jansen  * This is only needed for file format version 1 in "write" mode.
50559599516SKenneth E. Jansen  */
50659599516SKenneth E. Jansen int computeMHSize(int nfields, int nppf, int version) {
50710d56689SCameron Smith   int mhsize=0;
50859599516SKenneth E. Jansen   if(version == 1) {
50959599516SKenneth 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
51059599516SKenneth E. Jansen     int meta_info_size = VERSION_INFO_HEADER_SIZE;
5112dd307a1SCameron Smith     int actual_size =  nfields * nppf * sizeof(unsigned long) + meta_info_size;
51259599516SKenneth E. Jansen     //printf("actual_size = %d, offset table size = %d\n", actual_size,  nfields * nppf * sizeof(long long));
51359599516SKenneth E. Jansen     if (actual_size > DefaultMHSize) {
51459599516SKenneth E. Jansen       mhsize = (int) ceil( (double) actual_size/DefaultMHSize); // it's rounded to ceiling of this value
51559599516SKenneth E. Jansen       mhsize *= DefaultMHSize;
51659599516SKenneth E. Jansen     }
51759599516SKenneth E. Jansen     else {
51859599516SKenneth E. Jansen       mhsize = DefaultMHSize;
51959599516SKenneth E. Jansen     }
520b744cbcaSCameron Smith   } else {
521b744cbcaSCameron Smith     int rank = 0;
522b744cbcaSCameron Smith     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
523b744cbcaSCameron Smith     if(!rank) {
524b744cbcaSCameron Smith       fprintf(stderr,
525b744cbcaSCameron Smith           "ERROR invalid version passed to %s... exiting\n", __func__);
526b744cbcaSCameron Smith       exit(EXIT_FAILURE);
527b744cbcaSCameron Smith     }
52859599516SKenneth E. Jansen   }
52959599516SKenneth E. Jansen   return mhsize;
53059599516SKenneth E. Jansen }
53159599516SKenneth E. Jansen 
53259599516SKenneth E. Jansen /**
53359599516SKenneth E. Jansen  * Computes correct color of a rank according to number of files.
53459599516SKenneth E. Jansen  */
53559599516SKenneth E. Jansen extern "C" int computeColor( int myrank, int numprocs, int nfiles) {
53659599516SKenneth E. Jansen   int color =
53759599516SKenneth E. Jansen     (int)(myrank / (numprocs / nfiles));
53859599516SKenneth E. Jansen   return color;
53959599516SKenneth E. Jansen }
54059599516SKenneth E. Jansen 
54159599516SKenneth E. Jansen 
54259599516SKenneth E. Jansen /**
54359599516SKenneth E. Jansen  * Check the file descriptor.
54459599516SKenneth E. Jansen  */
545103be424SCameron Smith void checkFileDescriptor(const char fctname[], int*  fileDescriptor ) {
54659599516SKenneth E. Jansen   if ( *fileDescriptor < 0 ) {
54759599516SKenneth E. Jansen     printf("Error: File descriptor = %d in %s\n",*fileDescriptor,fctname);
54859599516SKenneth E. Jansen     exit(1);
54959599516SKenneth E. Jansen   }
55059599516SKenneth E. Jansen }
55159599516SKenneth E. Jansen 
55259599516SKenneth E. Jansen /**
55359599516SKenneth E. Jansen  * Initialize the file struct members and allocate space for file struct
55459599516SKenneth E. Jansen  * buffers.
55559599516SKenneth E. Jansen  *
55659599516SKenneth E. Jansen  * Note: this function is only called when we are using new format. Old POSIX
55759599516SKenneth E. Jansen  * format should skip this routine and call openfile() directly instead.
55859599516SKenneth E. Jansen  */
55959599516SKenneth E. Jansen int initphmpiio( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[])
56059599516SKenneth E. Jansen {
56159599516SKenneth E. Jansen   // we init irank again in case query not called (e.g. syncIO write case)
56259599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &irank);
56359599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &mysize);
56459599516SKenneth E. Jansen 
56559599516SKenneth E. Jansen   phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d", irank, MasterHeaderSize);
56659599516SKenneth E. Jansen 
56759599516SKenneth E. Jansen   double timer_start, timer_end;
56859599516SKenneth E. Jansen   startTimer(&timer_start);
56959599516SKenneth E. Jansen 
57059599516SKenneth E. Jansen   char* imode = StringStripper( mode );
57159599516SKenneth E. Jansen 
57259599516SKenneth E. Jansen   // Note: if it's read, we presume query was called prior to init and
57359599516SKenneth E. Jansen   // MasterHeaderSize is already set to correct value from parsing header
57459599516SKenneth E. Jansen   // otherwise it's write then it needs some computation to be set
57559599516SKenneth E. Jansen   if ( cscompare( "read", imode ) ) {
57659599516SKenneth E. Jansen     // do nothing
57759599516SKenneth E. Jansen   }
57859599516SKenneth E. Jansen   else if( cscompare( "write", imode ) ) {
57959599516SKenneth E. Jansen     MasterHeaderSize =  computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
58059599516SKenneth E. Jansen   }
58159599516SKenneth E. Jansen   else {
58259599516SKenneth E. Jansen     printf("Error initphmpiio: can't recognize the mode %s", imode);
58359599516SKenneth E. Jansen     exit(1);
58459599516SKenneth E. Jansen   }
58559599516SKenneth E. Jansen   free ( imode );
58659599516SKenneth E. Jansen 
58759599516SKenneth E. Jansen   phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d", irank, MasterHeaderSize);
58859599516SKenneth E. Jansen 
58959599516SKenneth E. Jansen   int i, j;
59059599516SKenneth E. Jansen 
59159599516SKenneth E. Jansen   if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
59259599516SKenneth E. Jansen     printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
59359599516SKenneth E. Jansen     endTimer(&timer_end);
59459599516SKenneth E. Jansen     printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
59559599516SKenneth E. Jansen     return MAX_PHASTA_FILES_EXCEEDED;
59659599516SKenneth E. Jansen   }
59759599516SKenneth E. Jansen   //		else if( PhastaIONextActiveIndex == 0 )  //Hang in debug mode on Intrepid
59859599516SKenneth E. Jansen   //		{
59959599516SKenneth E. Jansen   //			for( i = 0; i < MAX_PHASTA_FILES; i++ );
60059599516SKenneth E. Jansen   //			{
60159599516SKenneth E. Jansen   //				PhastaIOActiveFiles[i] = NULL;
60259599516SKenneth E. Jansen   //			}
60359599516SKenneth E. Jansen   //		}
60459599516SKenneth E. Jansen 
60559599516SKenneth E. Jansen 
60659599516SKenneth E. Jansen   PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1,  sizeof( phastaio_file_t) );
60759599516SKenneth E. Jansen 
60859599516SKenneth E. Jansen   i = PhastaIONextActiveIndex;
60959599516SKenneth E. Jansen   PhastaIONextActiveIndex++;
61059599516SKenneth E. Jansen 
61159599516SKenneth E. Jansen   //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
61259599516SKenneth E. Jansen 
61359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize;  // what does this mean??? TODO
61459599516SKenneth E. Jansen 
61559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->Wrong_Endian = false;
61659599516SKenneth E. Jansen 
61759599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFields = *nfields;
61859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nPPF = *nppf;
61959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFiles = *nfiles;
62059599516SKenneth E. Jansen   MPI_Comm_rank(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->myrank));
62159599516SKenneth E. Jansen   MPI_Comm_size(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->numprocs));
62259599516SKenneth E. Jansen 
62359599516SKenneth E. Jansen 
62459599516SKenneth E. Jansen   if( *nfiles > 1 ) { // split the ranks according to each mpiio file
62559599516SKenneth E. Jansen 
62659599516SKenneth E. Jansen     if ( s_assign_local_comm == 0) { // call mpi_comm_split for the first (and only) time
62759599516SKenneth E. Jansen 
62859599516SKenneth E. Jansen       if (PhastaIOActiveFiles[i]->myrank == 0) printf("Building subcommunicator\n");
62959599516SKenneth E. Jansen 
63059599516SKenneth E. Jansen       int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
63159599516SKenneth E. Jansen       MPI_Comm_split(MPI_COMM_WORLD,
63259599516SKenneth E. Jansen           color,
63359599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->myrank,
63459599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_comm));
63559599516SKenneth E. Jansen       MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
63659599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_numprocs));
63759599516SKenneth E. Jansen       MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
63859599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->local_myrank));
63959599516SKenneth E. Jansen 
64059599516SKenneth E. Jansen       // back up now these variables so that we do not need to call comm_split again
64159599516SKenneth E. Jansen       s_local_comm = PhastaIOActiveFiles[i]->local_comm;
64259599516SKenneth E. Jansen       s_local_size = PhastaIOActiveFiles[i]->local_numprocs;
64359599516SKenneth E. Jansen       s_local_rank = PhastaIOActiveFiles[i]->local_myrank;
64459599516SKenneth E. Jansen       s_assign_local_comm = 1;
64559599516SKenneth E. Jansen     }
64659599516SKenneth E. Jansen     else { // recycle the subcommunicator
64759599516SKenneth E. Jansen       if (PhastaIOActiveFiles[i]->myrank == 0) printf("Recycling subcommunicator\n");
64859599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_comm = s_local_comm;
64959599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_numprocs = s_local_size;
65059599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->local_myrank = s_local_rank;
65159599516SKenneth E. Jansen     }
65259599516SKenneth E. Jansen   }
65359599516SKenneth E. Jansen   else { // *nfiles == 1 here - no need to call mpi_comm_split here
65459599516SKenneth E. Jansen 
65559599516SKenneth E. Jansen     if (PhastaIOActiveFiles[i]->myrank == 0) printf("Bypassing subcommunicator\n");
65659599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_comm = MPI_COMM_WORLD;
65759599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_numprocs = PhastaIOActiveFiles[i]->numprocs;
65859599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->local_myrank = PhastaIOActiveFiles[i]->myrank;
65959599516SKenneth E. Jansen 
66059599516SKenneth E. Jansen   }
66159599516SKenneth E. Jansen 
66259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nppp =
66359599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
66459599516SKenneth E. Jansen 
66559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
66659599516SKenneth E. Jansen     (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
66759599516SKenneth E. Jansen     (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
66859599516SKenneth E. Jansen 
66959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_offset_table =
6702dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
67159599516SKenneth E. Jansen 
67259599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_read_table =
6732dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
67459599516SKenneth E. Jansen 
67559599516SKenneth E. Jansen   for (j=0; j<*nfields; j++)
67659599516SKenneth E. Jansen   {
67759599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[j] =
6782dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
67959599516SKenneth E. Jansen 
68059599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_read_table[j] =
6812dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
68259599516SKenneth E. Jansen   }
68359599516SKenneth E. Jansen   *filehandle = i;
68459599516SKenneth E. Jansen 
68559599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
68659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
68759599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
68859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
68959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
69059599516SKenneth E. Jansen 
69159599516SKenneth E. Jansen   // Time monitoring
69259599516SKenneth E. Jansen   endTimer(&timer_end);
69359599516SKenneth E. Jansen   printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
69459599516SKenneth E. Jansen 
69559599516SKenneth E. Jansen   phprintf_0("Info initphmpiio: quiting function");
69659599516SKenneth E. Jansen 
69759599516SKenneth E. Jansen   return i;
69859599516SKenneth E. Jansen }
69959599516SKenneth E. Jansen 
70059599516SKenneth E. Jansen /**
70159599516SKenneth E. Jansen  * Destruct the file struct and free buffers allocated in init function.
70259599516SKenneth E. Jansen  */
70359599516SKenneth E. Jansen void finalizephmpiio( int *fileDescriptor )
70459599516SKenneth E. Jansen {
70559599516SKenneth E. Jansen   double timer_start, timer_end;
70659599516SKenneth E. Jansen   startTimer(&timer_start);
70759599516SKenneth E. Jansen 
70859599516SKenneth E. Jansen   int i, j;
70959599516SKenneth E. Jansen   i = *fileDescriptor;
71059599516SKenneth E. Jansen   //PhastaIONextActiveIndex--;
71159599516SKenneth E. Jansen 
71259599516SKenneth E. Jansen   /* //free the offset table for this phasta file */
71359599516SKenneth E. Jansen   //for(j=0; j<MAX_FIELDS_NUMBER; j++) //Danger: undefined behavior for my_*_table.[j] not allocated or not initialized to NULL
71459599516SKenneth E. Jansen   for(j=0; j<PhastaIOActiveFiles[i]->nFields; j++)
71559599516SKenneth E. Jansen   {
71659599516SKenneth E. Jansen     free( PhastaIOActiveFiles[i]->my_offset_table[j]);
71759599516SKenneth E. Jansen     free( PhastaIOActiveFiles[i]->my_read_table[j]);
71859599516SKenneth E. Jansen   }
71959599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->my_offset_table );
72059599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->my_read_table );
72159599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->master_header );
72259599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->double_chunk );
72359599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->int_chunk );
72459599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->read_double_chunk );
72559599516SKenneth E. Jansen   free ( PhastaIOActiveFiles[i]->read_int_chunk );
72659599516SKenneth E. Jansen 
727946a6cddSCameron Smith   if( PhastaIOActiveFiles[i]->nFiles > 1 && s_assign_local_comm ) { // the comm was split
728946a6cddSCameron Smith     if (PhastaIOActiveFiles[i]->myrank == 0) printf("Freeing subcommunicator\n");
729946a6cddSCameron Smith     s_assign_local_comm = 0;
730946a6cddSCameron Smith     MPI_Comm_free(&(PhastaIOActiveFiles[i]->local_comm));
731946a6cddSCameron Smith   }
732946a6cddSCameron Smith 
73359599516SKenneth E. Jansen   free( PhastaIOActiveFiles[i]);
73459599516SKenneth E. Jansen 
73559599516SKenneth E. Jansen   endTimer(&timer_end);
73659599516SKenneth E. Jansen   printPerf("finalizempiio", timer_start, timer_end, 0, 0, "");
73759599516SKenneth E. Jansen 
73859599516SKenneth E. Jansen   PhastaIONextActiveIndex--;
73959599516SKenneth E. Jansen }
74059599516SKenneth E. Jansen 
74159599516SKenneth E. Jansen 
74259599516SKenneth E. Jansen /**
74359599516SKenneth E. Jansen  * Special init for M2N in order to create a subcommunicator for the reduced solution (requires PRINT_PERF to be false for now)
74459599516SKenneth E. Jansen  * Initialize the file struct members and allocate space for file struct buffers.
74559599516SKenneth E. Jansen  *
74659599516SKenneth E. Jansen  * Note: this function is only called when we are using new format. Old POSIX
74759599516SKenneth E. Jansen  * format should skip this routine and call openfile() directly instead.
74859599516SKenneth E. Jansen  */
74959599516SKenneth E. Jansen int initphmpiiosub( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[],MPI_Comm my_local_comm)
75059599516SKenneth E. Jansen {
75159599516SKenneth E. Jansen   // we init irank again in case query not called (e.g. syncIO write case)
75259599516SKenneth E. Jansen 
75359599516SKenneth E. Jansen   MPI_Comm_rank(my_local_comm, &irank);
75459599516SKenneth E. Jansen   MPI_Comm_size(my_local_comm, &mysize);
75559599516SKenneth E. Jansen 
75659599516SKenneth E. Jansen   phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d", irank, MasterHeaderSize);
75759599516SKenneth E. Jansen 
75859599516SKenneth E. Jansen   double timer_start, timer_end;
75959599516SKenneth E. Jansen   startTimer(&timer_start);
76059599516SKenneth E. Jansen 
76159599516SKenneth E. Jansen   char* imode = StringStripper( mode );
76259599516SKenneth E. Jansen 
76359599516SKenneth E. Jansen   // Note: if it's read, we presume query was called prior to init and
76459599516SKenneth E. Jansen   // MasterHeaderSize is already set to correct value from parsing header
76559599516SKenneth E. Jansen   // otherwise it's write then it needs some computation to be set
76659599516SKenneth E. Jansen   if ( cscompare( "read", imode ) ) {
76759599516SKenneth E. Jansen     // do nothing
76859599516SKenneth E. Jansen   }
76959599516SKenneth E. Jansen   else if( cscompare( "write", imode ) ) {
77059599516SKenneth E. Jansen     MasterHeaderSize =  computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
77159599516SKenneth E. Jansen   }
77259599516SKenneth E. Jansen   else {
77359599516SKenneth E. Jansen     printf("Error initphmpiio: can't recognize the mode %s", imode);
77459599516SKenneth E. Jansen     exit(1);
77559599516SKenneth E. Jansen   }
77659599516SKenneth E. Jansen   free ( imode );
77759599516SKenneth E. Jansen 
77859599516SKenneth E. Jansen   phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d", irank, MasterHeaderSize);
77959599516SKenneth E. Jansen 
78059599516SKenneth E. Jansen   int i, j;
78159599516SKenneth E. Jansen 
78259599516SKenneth E. Jansen   if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
78359599516SKenneth E. Jansen     printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
78459599516SKenneth E. Jansen     endTimer(&timer_end);
78559599516SKenneth E. Jansen     printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
78659599516SKenneth E. Jansen     return MAX_PHASTA_FILES_EXCEEDED;
78759599516SKenneth E. Jansen   }
78859599516SKenneth E. Jansen   //		else if( PhastaIONextActiveIndex == 0 )  //Hang in debug mode on Intrepid
78959599516SKenneth E. Jansen   //		{
79059599516SKenneth E. Jansen   //			for( i = 0; i < MAX_PHASTA_FILES; i++ );
79159599516SKenneth E. Jansen   //			{
79259599516SKenneth E. Jansen   //				PhastaIOActiveFiles[i] = NULL;
79359599516SKenneth E. Jansen   //			}
79459599516SKenneth E. Jansen   //		}
79559599516SKenneth E. Jansen 
79659599516SKenneth E. Jansen 
79759599516SKenneth E. Jansen   PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1,  sizeof( phastaio_file_t) );
79859599516SKenneth E. Jansen 
79959599516SKenneth E. Jansen   i = PhastaIONextActiveIndex;
80059599516SKenneth E. Jansen   PhastaIONextActiveIndex++;
80159599516SKenneth E. Jansen 
80259599516SKenneth E. Jansen   //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
80359599516SKenneth E. Jansen 
80459599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize;  // what does this mean??? TODO
80559599516SKenneth E. Jansen 
80659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->Wrong_Endian = false;
80759599516SKenneth E. Jansen 
80859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFields = *nfields;
80959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nPPF = *nppf;
81059599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nFiles = *nfiles;
81159599516SKenneth E. Jansen   MPI_Comm_rank(my_local_comm, &(PhastaIOActiveFiles[i]->myrank));
81259599516SKenneth E. Jansen   MPI_Comm_size(my_local_comm, &(PhastaIOActiveFiles[i]->numprocs));
81359599516SKenneth E. Jansen 
81459599516SKenneth E. Jansen   int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
81559599516SKenneth E. Jansen   MPI_Comm_split(my_local_comm,
81659599516SKenneth E. Jansen       color,
81759599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->myrank,
81859599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_comm));
81959599516SKenneth E. Jansen   MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
82059599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_numprocs));
82159599516SKenneth E. Jansen   MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
82259599516SKenneth E. Jansen       &(PhastaIOActiveFiles[i]->local_myrank));
82359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->nppp =
82459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
82559599516SKenneth E. Jansen 
82659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
82759599516SKenneth E. Jansen     (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
82859599516SKenneth E. Jansen     (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
82959599516SKenneth E. Jansen 
83059599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_offset_table =
8312dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
83259599516SKenneth E. Jansen 
83359599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->my_read_table =
8342dd307a1SCameron Smith     ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
83559599516SKenneth E. Jansen 
83659599516SKenneth E. Jansen   for (j=0; j<*nfields; j++)
83759599516SKenneth E. Jansen   {
83859599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[j] =
8392dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
84059599516SKenneth E. Jansen 
84159599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_read_table[j] =
8422dd307a1SCameron Smith       ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
84359599516SKenneth E. Jansen   }
84459599516SKenneth E. Jansen   *filehandle = i;
84559599516SKenneth E. Jansen 
84659599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
84759599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
84859599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
84959599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
85059599516SKenneth E. Jansen   PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
85159599516SKenneth E. Jansen 
85259599516SKenneth E. Jansen   // Time monitoring
85359599516SKenneth E. Jansen   endTimer(&timer_end);
85459599516SKenneth E. Jansen   printPerf("initphmpiiosub", timer_start, timer_end, 0, 0, "");
85559599516SKenneth E. Jansen 
85659599516SKenneth E. Jansen   phprintf_0("Info initphmpiiosub: quiting function");
85759599516SKenneth E. Jansen 
85859599516SKenneth E. Jansen   return i;
85959599516SKenneth E. Jansen }
86059599516SKenneth E. Jansen 
86159599516SKenneth E. Jansen 
86259599516SKenneth E. Jansen 
86359599516SKenneth E. Jansen /** open file for both POSIX and MPI-IO syncIO format.
86459599516SKenneth E. Jansen  *
86559599516SKenneth E. Jansen  * If it's old POSIX format, simply call posix fopen().
86659599516SKenneth E. Jansen  *
86759599516SKenneth E. Jansen  * If it's MPI-IO foramt:
86859599516SKenneth E. Jansen  * in "read" mode, it builds the header table that points to the offset of
86959599516SKenneth E. Jansen  * fields for parts;
87059599516SKenneth E. Jansen  * in "write" mode, it opens the file with MPI-IO open routine.
87159599516SKenneth E. Jansen  */
872103be424SCameron Smith void openfile(const char filename[], const char mode[], int*  fileDescriptor )
87359599516SKenneth E. Jansen {
87459599516SKenneth E. Jansen   phprintf_0("Info: entering openfile");
87559599516SKenneth E. Jansen 
87659599516SKenneth E. Jansen   double timer_start, timer_end;
87759599516SKenneth E. Jansen   startTimer(&timer_start);
87859599516SKenneth E. Jansen 
87959599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 )
88059599516SKenneth E. Jansen   {
88159599516SKenneth E. Jansen     FILE* file=NULL ;
88259599516SKenneth E. Jansen     *fileDescriptor = 0;
88359599516SKenneth E. Jansen     char* fname = StringStripper( filename );
88459599516SKenneth E. Jansen     char* imode = StringStripper( mode );
88559599516SKenneth E. Jansen 
88659599516SKenneth E. Jansen     if ( cscompare( "read", imode ) ) file = fopen(fname, "rb" );
88759599516SKenneth E. Jansen     else if( cscompare( "write", imode ) ) file = fopen(fname, "wb" );
88859599516SKenneth E. Jansen     else if( cscompare( "append", imode ) ) file = fopen(fname, "ab" );
88959599516SKenneth E. Jansen 
89059599516SKenneth E. Jansen     if ( !file ){
89159599516SKenneth E. Jansen       fprintf(stderr,"Error openfile: unable to open file %s",fname ) ;
89259599516SKenneth E. Jansen     } else {
89359599516SKenneth E. Jansen       fileArray.push_back( file );
89459599516SKenneth E. Jansen       byte_order.push_back( false );
89559599516SKenneth E. Jansen       header_type.push_back( sizeof(int) );
89659599516SKenneth E. Jansen       *fileDescriptor = fileArray.size();
89759599516SKenneth E. Jansen     }
89859599516SKenneth E. Jansen     free (fname);
89959599516SKenneth E. Jansen     free (imode);
90059599516SKenneth E. Jansen   }
90159599516SKenneth E. Jansen   else // else it would be parallel I/O, opposed to posix io
90259599516SKenneth E. Jansen   {
90359599516SKenneth E. Jansen     char* fname = StringStripper( filename );
90459599516SKenneth E. Jansen     char* imode = StringStripper( mode );
90559599516SKenneth E. Jansen     int rc;
90659599516SKenneth E. Jansen     int i = *fileDescriptor;
90759599516SKenneth E. Jansen     checkFileDescriptor("openfile",&i);
90859599516SKenneth E. Jansen     char* token;
90959599516SKenneth E. Jansen 
91059599516SKenneth E. Jansen     if ( cscompare( "read", imode ) )
91159599516SKenneth E. Jansen     {
91259599516SKenneth E. Jansen       //	      if (PhastaIOActiveFiles[i]->myrank == 0)
91359599516SKenneth E. Jansen       //                printf("\n **********\nRead open ... ... regular version\n");
91459599516SKenneth E. Jansen 
91559599516SKenneth E. Jansen       rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
91659599516SKenneth E. Jansen           fname,
91759599516SKenneth E. Jansen           MPI_MODE_RDONLY,
91859599516SKenneth E. Jansen           MPI_INFO_NULL,
91959599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->file_handle) );
92059599516SKenneth E. Jansen 
92159599516SKenneth E. Jansen       if(rc)
92259599516SKenneth E. Jansen       {
92359599516SKenneth E. Jansen         *fileDescriptor = UNABLE_TO_OPEN_FILE;
92459599516SKenneth E. Jansen         printf("Error openfile: Unable to open file %s! File descriptor = %d\n",fname,*fileDescriptor);
92559599516SKenneth E. Jansen         endTimer(&timer_end);
92659599516SKenneth E. Jansen         printPerf("openfile", timer_start, timer_end, 0, 0, "");
92759599516SKenneth E. Jansen         return;
92859599516SKenneth E. Jansen       }
92959599516SKenneth E. Jansen 
93059599516SKenneth E. Jansen       MPI_Status read_tag_status;
93159599516SKenneth E. Jansen       char read_out_tag[MAX_FIELDS_NAME_LENGTH];
93259599516SKenneth E. Jansen       int j;
93359599516SKenneth E. Jansen       int magic_number;
93459599516SKenneth E. Jansen 
93559599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
93659599516SKenneth E. Jansen         MPI_File_read_at( PhastaIOActiveFiles[i]->file_handle,
93759599516SKenneth E. Jansen             0,
93859599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header,
93959599516SKenneth E. Jansen             MasterHeaderSize,
94059599516SKenneth E. Jansen             MPI_CHAR,
94159599516SKenneth E. Jansen             &read_tag_status );
94259599516SKenneth E. Jansen       }
94359599516SKenneth E. Jansen 
94459599516SKenneth E. Jansen       MPI_Bcast( PhastaIOActiveFiles[i]->master_header,
94559599516SKenneth E. Jansen           MasterHeaderSize,
94659599516SKenneth E. Jansen           MPI_CHAR,
94759599516SKenneth E. Jansen           0,
94859599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->local_comm );
94959599516SKenneth E. Jansen 
95059599516SKenneth E. Jansen       memcpy( read_out_tag,
95159599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->master_header,
95259599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
95359599516SKenneth E. Jansen 
95459599516SKenneth E. Jansen       if ( cscompare ("MPI_IO_Tag",read_out_tag) )
95559599516SKenneth E. Jansen       {
95659599516SKenneth E. Jansen         // Test endianess ...
95759599516SKenneth E. Jansen         memcpy ( &magic_number,
95859599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
95959599516SKenneth E. Jansen             sizeof(int) );                                                   // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
96059599516SKenneth E. Jansen 
96159599516SKenneth E. Jansen         if ( magic_number != ENDIAN_TEST_NUMBER )
96259599516SKenneth E. Jansen         {
96359599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->Wrong_Endian = true;
96459599516SKenneth E. Jansen         }
96559599516SKenneth E. Jansen 
96659599516SKenneth E. Jansen         memcpy( read_out_tag,
96759599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH+1, // TODO: WHY +1???
96859599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH );
96959599516SKenneth E. Jansen 
97059599516SKenneth E. Jansen         // Read in # fields ...
97159599516SKenneth E. Jansen         token = strtok ( read_out_tag, ":" );
97259599516SKenneth E. Jansen         token = strtok( NULL," ,;<>" );
97359599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->nFields = atoi( token );
97459599516SKenneth E. Jansen 
9752dd307a1SCameron Smith         unsigned long **header_table;
9762dd307a1SCameron Smith         header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
97759599516SKenneth E. Jansen 
97859599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
97959599516SKenneth E. Jansen         {
9802dd307a1SCameron Smith           header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
98159599516SKenneth E. Jansen         }
98259599516SKenneth E. Jansen 
98359599516SKenneth E. Jansen         // Read in the offset table ...
98459599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
98559599516SKenneth E. Jansen         {
98659599516SKenneth E. Jansen           if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
98759599516SKenneth E. Jansen             memcpy( header_table[j],
98859599516SKenneth E. Jansen                 PhastaIOActiveFiles[i]->master_header +
98959599516SKenneth E. Jansen                 VERSION_INFO_HEADER_SIZE +
9902dd307a1SCameron Smith                 j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
9912dd307a1SCameron Smith                 PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
99259599516SKenneth E. Jansen           }
99359599516SKenneth E. Jansen 
99459599516SKenneth E. Jansen           MPI_Scatter( header_table[j],
99559599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->nppp,
99659599516SKenneth E. Jansen               MPI_LONG_LONG_INT,
99759599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->my_read_table[j],
99859599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->nppp,
99959599516SKenneth E. Jansen               MPI_LONG_LONG_INT,
100059599516SKenneth E. Jansen               0,
100159599516SKenneth E. Jansen               PhastaIOActiveFiles[i]->local_comm );
100259599516SKenneth E. Jansen 
100359599516SKenneth E. Jansen           // Swap byte order if endianess is different ...
100459599516SKenneth E. Jansen           if ( PhastaIOActiveFiles[i]->Wrong_Endian ) {
100559599516SKenneth E. Jansen             SwapArrayByteOrder( PhastaIOActiveFiles[i]->my_read_table[j],
10064187599aSCameron Smith                 sizeof(unsigned long),
100759599516SKenneth E. Jansen                 PhastaIOActiveFiles[i]->nppp );
100859599516SKenneth E. Jansen           }
100959599516SKenneth E. Jansen         }
101059599516SKenneth E. Jansen 
101159599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
101259599516SKenneth E. Jansen           free ( header_table[j] );
101359599516SKenneth E. Jansen         }
101459599516SKenneth E. Jansen         free (header_table);
101559599516SKenneth E. Jansen 
101659599516SKenneth E. Jansen       } // end of if MPI_IO_TAG
101759599516SKenneth E. Jansen       else //else not valid MPI file
101859599516SKenneth E. Jansen       {
101959599516SKenneth E. Jansen         *fileDescriptor = NOT_A_MPI_FILE;
102059599516SKenneth 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);
102159599516SKenneth E. Jansen         //Printing MasterHeaderSize is useful to test a compiler bug on Intrepid BGP
102259599516SKenneth E. Jansen         endTimer(&timer_end);
102359599516SKenneth E. Jansen         printPerf("openfile", timer_start, timer_end, 0, 0, "");
102459599516SKenneth E. Jansen         return;
102559599516SKenneth E. Jansen       }
102659599516SKenneth E. Jansen     } // end of if "read"
102759599516SKenneth E. Jansen     else if( cscompare( "write", imode ) )
102859599516SKenneth E. Jansen     {
102959599516SKenneth E. Jansen       rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
103059599516SKenneth E. Jansen           fname,
103159599516SKenneth E. Jansen           MPI_MODE_WRONLY | MPI_MODE_CREATE,
103259599516SKenneth E. Jansen           MPI_INFO_NULL,
103359599516SKenneth E. Jansen           &(PhastaIOActiveFiles[i]->file_handle) );
103459599516SKenneth E. Jansen       if(rc)
103559599516SKenneth E. Jansen       {
103659599516SKenneth E. Jansen         *fileDescriptor = UNABLE_TO_OPEN_FILE;
103759599516SKenneth E. Jansen         return;
103859599516SKenneth E. Jansen       }
103959599516SKenneth E. Jansen     } // end of if "write"
104059599516SKenneth E. Jansen     free (fname);
104159599516SKenneth E. Jansen     free (imode);
104259599516SKenneth E. Jansen   } // end of if FileIndex != 0
104359599516SKenneth E. Jansen 
104459599516SKenneth E. Jansen   endTimer(&timer_end);
104559599516SKenneth E. Jansen   printPerf("openfile", timer_start, timer_end, 0, 0, "");
104659599516SKenneth E. Jansen }
104759599516SKenneth E. Jansen 
104859599516SKenneth E. Jansen /** close file for both POSIX and MPI-IO syncIO format.
104959599516SKenneth E. Jansen  *
105059599516SKenneth E. Jansen  * If it's old POSIX format, simply call posix fclose().
105159599516SKenneth E. Jansen  *
105259599516SKenneth E. Jansen  * If it's MPI-IO foramt:
105359599516SKenneth E. Jansen  * in "read" mode, it simply close file with MPI-IO close routine.
105459599516SKenneth E. Jansen  * in "write" mode, rank 0 in each group will re-assemble the master header and
105559599516SKenneth E. Jansen  * offset table and write to the beginning of file, then close the file.
105659599516SKenneth E. Jansen  */
1057103be424SCameron Smith void closefile( int* fileDescriptor, const char mode[] )
105859599516SKenneth E. Jansen {
105959599516SKenneth E. Jansen   double timer_start, timer_end;
106059599516SKenneth E. Jansen   startTimer(&timer_start);
106159599516SKenneth E. Jansen 
106259599516SKenneth E. Jansen   int i = *fileDescriptor;
106359599516SKenneth E. Jansen   checkFileDescriptor("closefile",&i);
106459599516SKenneth E. Jansen 
106559599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
106659599516SKenneth E. Jansen     char* imode = StringStripper( mode );
106759599516SKenneth E. Jansen 
106859599516SKenneth E. Jansen     if( cscompare( "write", imode )
106959599516SKenneth E. Jansen         || cscompare( "append", imode ) ) {
107059599516SKenneth E. Jansen       fflush( fileArray[ *fileDescriptor - 1 ] );
107159599516SKenneth E. Jansen     }
107259599516SKenneth E. Jansen 
107359599516SKenneth E. Jansen     fclose( fileArray[ *fileDescriptor - 1 ] );
107459599516SKenneth E. Jansen     free (imode);
107559599516SKenneth E. Jansen   }
107659599516SKenneth E. Jansen   else {
107759599516SKenneth E. Jansen     char* imode = StringStripper( mode );
107859599516SKenneth E. Jansen 
107959599516SKenneth E. Jansen     //write master header here:
108059599516SKenneth E. Jansen     if ( cscompare( "write", imode ) ) {
108159599516SKenneth E. Jansen       //	      if ( PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields < 2*ONE_MEGABYTE/8 )  //SHOULD BE CHECKED
108259599516SKenneth E. Jansen       //		MasterHeaderSize = 4*ONE_MEGABYTE;
108359599516SKenneth E. Jansen       //	      else
108459599516SKenneth E. Jansen       //		MasterHeaderSize = 4*ONE_MEGABYTE + PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields * 8 - 2*ONE_MEGABYTE;
108559599516SKenneth E. Jansen 
108659599516SKenneth E. Jansen       MasterHeaderSize = computeMHSize( PhastaIOActiveFiles[i]->nFields, PhastaIOActiveFiles[i]->nPPF, LATEST_WRITE_VERSION);
108759599516SKenneth E. Jansen       phprintf_0("Info closefile: myrank = %d, MasterHeaderSize = %d\n", PhastaIOActiveFiles[i]->myrank, MasterHeaderSize);
108859599516SKenneth E. Jansen 
108959599516SKenneth E. Jansen       MPI_Status write_header_status;
109059599516SKenneth E. Jansen       char mpi_tag[MAX_FIELDS_NAME_LENGTH];
109159599516SKenneth E. Jansen       char version[MAX_FIELDS_NAME_LENGTH/4];
109259599516SKenneth E. Jansen       char mhsize[MAX_FIELDS_NAME_LENGTH/4];
109359599516SKenneth E. Jansen       int magic_number = ENDIAN_TEST_NUMBER;
109459599516SKenneth E. Jansen 
109559599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 )
109659599516SKenneth E. Jansen       {
109759599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
109859599516SKenneth E. Jansen         sprintf(mpi_tag, "MPI_IO_Tag : ");
109959599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header,
110059599516SKenneth E. Jansen             mpi_tag,
110159599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
110259599516SKenneth E. Jansen 
110359599516SKenneth E. Jansen         bzero((void*)version,MAX_FIELDS_NAME_LENGTH/4);
110459599516SKenneth E. Jansen         // this version is "1", print version in ASCII
110559599516SKenneth E. Jansen         sprintf(version, "version : %d",1);
110659599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/2,
110759599516SKenneth E. Jansen             version,
110859599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4);
110959599516SKenneth E. Jansen 
111059599516SKenneth E. Jansen         // master header size is computed using the formula above
111159599516SKenneth E. Jansen         bzero((void*)mhsize,MAX_FIELDS_NAME_LENGTH/4);
111259599516SKenneth E. Jansen         sprintf(mhsize, "mhsize : ");
111359599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/4*3,
111459599516SKenneth E. Jansen             mhsize,
111559599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH/4);
111659599516SKenneth E. Jansen 
111759599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
111859599516SKenneth E. Jansen         sprintf(mpi_tag,
111959599516SKenneth E. Jansen             "\nnFields : %d\n",
112059599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nFields);
112159599516SKenneth E. Jansen         memcpy(PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH,
112259599516SKenneth E. Jansen             mpi_tag,
112359599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
112459599516SKenneth E. Jansen 
112559599516SKenneth E. Jansen         bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
112659599516SKenneth E. Jansen         sprintf(mpi_tag, "\nnPPF : %d\n", PhastaIOActiveFiles[i]->nPPF);
112759599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+
112859599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nFields *
112959599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH +
113059599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH * 2,
113159599516SKenneth E. Jansen             mpi_tag,
113259599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH);
113359599516SKenneth E. Jansen 
113459599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
113559599516SKenneth E. Jansen             &magic_number,
113659599516SKenneth E. Jansen             sizeof(int));
113759599516SKenneth E. Jansen 
113859599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("mhsize : ") -1 + MAX_FIELDS_NAME_LENGTH/4*3,
113959599516SKenneth E. Jansen             &MasterHeaderSize,
114059599516SKenneth E. Jansen             sizeof(int));
114159599516SKenneth E. Jansen       }
114259599516SKenneth E. Jansen 
114359599516SKenneth E. Jansen       int j = 0;
11442dd307a1SCameron Smith       unsigned long **header_table;
11452dd307a1SCameron Smith       header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
114659599516SKenneth E. Jansen 
114759599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
11482dd307a1SCameron Smith         header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
114959599516SKenneth E. Jansen       }
115059599516SKenneth E. Jansen 
115159599516SKenneth E. Jansen       //if( irank == 0 ) printf("gonna mpi_gather, myrank = %d\n", irank);
115259599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
115359599516SKenneth E. Jansen         MPI_Gather( PhastaIOActiveFiles[i]->my_offset_table[j],
115459599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nppp,
115559599516SKenneth E. Jansen             MPI_LONG_LONG_INT,
115659599516SKenneth E. Jansen             header_table[j],
115759599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->nppp,
115859599516SKenneth E. Jansen             MPI_LONG_LONG_INT,
115959599516SKenneth E. Jansen             0,
116059599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->local_comm );
116159599516SKenneth E. Jansen       }
116259599516SKenneth E. Jansen 
116359599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
116459599516SKenneth E. Jansen 
116559599516SKenneth E. Jansen         //if( irank == 0 ) printf("gonna memcpy for every procs, myrank = %d\n", irank);
116659599516SKenneth E. Jansen         for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
116759599516SKenneth E. Jansen           memcpy ( PhastaIOActiveFiles[i]->master_header +
116859599516SKenneth E. Jansen               VERSION_INFO_HEADER_SIZE +
11692dd307a1SCameron Smith               j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
117059599516SKenneth E. Jansen               header_table[j],
11712dd307a1SCameron Smith               PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
117259599516SKenneth E. Jansen         }
117359599516SKenneth E. Jansen 
117459599516SKenneth E. Jansen         //if( irank == 0 ) printf("gonna file_write_at(), myrank = %d\n", irank);
117559599516SKenneth E. Jansen         MPI_File_write_at( PhastaIOActiveFiles[i]->file_handle,
117659599516SKenneth E. Jansen             0,
117759599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->master_header,
117859599516SKenneth E. Jansen             MasterHeaderSize,
117959599516SKenneth E. Jansen             MPI_CHAR,
118059599516SKenneth E. Jansen             &write_header_status );
118159599516SKenneth E. Jansen       }
118259599516SKenneth E. Jansen 
118359599516SKenneth E. Jansen       ////free(PhastaIOActiveFiles[i]->master_header);
118459599516SKenneth E. Jansen 
118559599516SKenneth E. Jansen       for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
118659599516SKenneth E. Jansen         free ( header_table[j] );
118759599516SKenneth E. Jansen       }
118859599516SKenneth E. Jansen       free (header_table);
118959599516SKenneth E. Jansen     }
119059599516SKenneth E. Jansen 
119159599516SKenneth E. Jansen     //if( irank == 0 ) printf("gonna file_close(), myrank = %d\n", irank);
119259599516SKenneth E. Jansen     MPI_File_close( &( PhastaIOActiveFiles[i]->file_handle ) );
119359599516SKenneth E. Jansen     free ( imode );
119459599516SKenneth E. Jansen   }
119559599516SKenneth E. Jansen 
119659599516SKenneth E. Jansen   endTimer(&timer_end);
119759599516SKenneth E. Jansen   printPerf("closefile_", timer_start, timer_end, 0, 0, "");
119859599516SKenneth E. Jansen }
119959599516SKenneth E. Jansen 
12003872e963SCameron Smith int readHeader( FILE* f, const char phrase[],
12013872e963SCameron Smith     int* params, int numParams, const char* iotype) {
12023872e963SCameron Smith   isBinary(iotype);
12033872e963SCameron Smith   return readHeader(f,phrase,params,numParams);
12043872e963SCameron Smith }
12053872e963SCameron Smith 
1206103be424SCameron Smith void readheader(
1207103be424SCameron Smith     int* fileDescriptor,
120859599516SKenneth E. Jansen     const  char keyphrase[],
120959599516SKenneth E. Jansen     void* valueArray,
121059599516SKenneth E. Jansen     int*  nItems,
121159599516SKenneth E. Jansen     const char  datatype[],
121259599516SKenneth E. Jansen     const char  iotype[] )
121359599516SKenneth E. Jansen {
121459599516SKenneth E. Jansen   double timer_start, timer_end;
1215d3337298SCameron Smith 
121659599516SKenneth E. Jansen   startTimer(&timer_start);
121759599516SKenneth E. Jansen 
121859599516SKenneth E. Jansen   int i = *fileDescriptor;
121959599516SKenneth E. Jansen   checkFileDescriptor("readheader",&i);
122059599516SKenneth E. Jansen 
122159599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
122259599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
122359599516SKenneth E. Jansen     FILE* fileObject;
122459599516SKenneth E. Jansen     int* valueListInt;
122559599516SKenneth E. Jansen 
122659599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
122759599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
122859599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
122959599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
123059599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
123159599516SKenneth E. Jansen       endTimer(&timer_end);
123259599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
123359599516SKenneth E. Jansen       return;
123459599516SKenneth E. Jansen     }
123559599516SKenneth E. Jansen 
1236961a4ff6SCameron Smith     LastHeaderKey[filePtr] = std::string(keyphrase);
123759599516SKenneth E. Jansen     LastHeaderNotFound = false;
123859599516SKenneth E. Jansen 
123959599516SKenneth E. Jansen     fileObject = fileArray[ filePtr ] ;
124059599516SKenneth E. Jansen     Wrong_Endian = byte_order[ filePtr ];
124159599516SKenneth E. Jansen 
124259599516SKenneth E. Jansen     isBinary( iotype );
124359599516SKenneth E. Jansen     typeSize( datatype );   //redundant call, just avoid a compiler warning.
124459599516SKenneth E. Jansen 
124559599516SKenneth E. Jansen     // right now we are making the assumption that we will only write integers
124659599516SKenneth E. Jansen     // on the header line.
124759599516SKenneth E. Jansen 
124859599516SKenneth E. Jansen     valueListInt = static_cast< int* >( valueArray );
124959599516SKenneth E. Jansen     int ierr = readHeader( fileObject ,
125059599516SKenneth E. Jansen         keyphrase,
125159599516SKenneth E. Jansen         valueListInt,
125259599516SKenneth E. Jansen         *nItems ) ;
125359599516SKenneth E. Jansen 
125459599516SKenneth E. Jansen     byte_order[ filePtr ] = Wrong_Endian ;
125559599516SKenneth E. Jansen 
125659599516SKenneth E. Jansen     if ( ierr ) LastHeaderNotFound = true;
125759599516SKenneth E. Jansen 
125859599516SKenneth E. Jansen     //return ; // don't return, go to the end to print perf
125959599516SKenneth E. Jansen   }
126059599516SKenneth E. Jansen   else {
126159599516SKenneth E. Jansen     int* valueListInt;
126259599516SKenneth E. Jansen     valueListInt = static_cast <int*>(valueArray);
1263400e9fc0SCameron Smith     char* token = NULL;
126459599516SKenneth E. Jansen     bool FOUND = false ;
126559599516SKenneth E. Jansen     isBinary( iotype );
126659599516SKenneth E. Jansen 
126759599516SKenneth E. Jansen     MPI_Status read_offset_status;
126859599516SKenneth E. Jansen     char read_out_tag[MAX_FIELDS_NAME_LENGTH];
1269400e9fc0SCameron Smith     memset(read_out_tag, '\0', MAX_FIELDS_NAME_LENGTH);
127059599516SKenneth E. Jansen     char readouttag[MAX_FIELDS_NUMBER][MAX_FIELDS_NAME_LENGTH];
127159599516SKenneth E. Jansen     int j;
127259599516SKenneth E. Jansen 
127359599516SKenneth E. Jansen     int string_length = strlen( keyphrase );
127459599516SKenneth E. Jansen     char* buffer = (char*) malloc ( string_length+1 );
127559599516SKenneth E. Jansen 
127659599516SKenneth E. Jansen     strcpy ( buffer, keyphrase );
127759599516SKenneth E. Jansen     buffer[ string_length ] = '\0';
127859599516SKenneth E. Jansen 
127959599516SKenneth E. Jansen     char* st2 = strtok ( buffer, "@" );
128059599516SKenneth E. Jansen     st2 = strtok (NULL, "@");
128159599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = atoi(st2);
128259599516SKenneth E. Jansen     if ( char* p = strpbrk(buffer, "@") )
128359599516SKenneth E. Jansen       *p = '\0';
128459599516SKenneth E. Jansen 
128559599516SKenneth E. Jansen     // Check if the user has input the right GPid
128659599516SKenneth E. Jansen     if ( ( PhastaIOActiveFiles[i]->GPid <=
128759599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->myrank *
128859599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->nppp )||
128959599516SKenneth E. Jansen         ( PhastaIOActiveFiles[i]->GPid >
129059599516SKenneth E. Jansen           ( PhastaIOActiveFiles[i]->myrank + 1 ) *
129159599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->nppp ) )
129259599516SKenneth E. Jansen     {
129359599516SKenneth E. Jansen       *fileDescriptor = NOT_A_MPI_FILE;
129459599516SKenneth 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);
129559599516SKenneth E. Jansen       // It is possible atoi could not generate a clear integer from st2 because of additional garbage character in keyphrase
129659599516SKenneth E. Jansen       endTimer(&timer_end);
129759599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
129859599516SKenneth E. Jansen       return;
129959599516SKenneth E. Jansen     }
130059599516SKenneth E. Jansen 
130159599516SKenneth E. Jansen     // Find the field we want ...
130259599516SKenneth E. Jansen     //for ( j = 0; j<MAX_FIELDS_NUMBER; j++ )
130359599516SKenneth E. Jansen     for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
130459599516SKenneth E. Jansen     {
130559599516SKenneth E. Jansen       memcpy( readouttag[j],
130659599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->master_header + j*MAX_FIELDS_NAME_LENGTH+MAX_FIELDS_NAME_LENGTH*2+1,
130759599516SKenneth E. Jansen           MAX_FIELDS_NAME_LENGTH-1 );
130859599516SKenneth E. Jansen     }
130959599516SKenneth E. Jansen 
131059599516SKenneth E. Jansen     for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
131159599516SKenneth E. Jansen     {
131259599516SKenneth E. Jansen       token = strtok ( readouttag[j], ":" );
131359599516SKenneth E. Jansen 
131459599516SKenneth E. Jansen       //if ( cscompare( buffer, token ) )
131559599516SKenneth E. Jansen       if ( cscompare( token , buffer ) && cscompare( buffer, token ) )
131659599516SKenneth 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").
131759599516SKenneth 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).
131859599516SKenneth 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.
131959599516SKenneth E. Jansen       {
132059599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->read_field_count = j;
132159599516SKenneth E. Jansen         FOUND = true;
132259599516SKenneth E. Jansen         //printf("buffer: %s | token: %s | j: %d\n",buffer,token,j);
132359599516SKenneth E. Jansen         break;
132459599516SKenneth E. Jansen       }
132559599516SKenneth E. Jansen     }
132659599516SKenneth E. Jansen     free(buffer);
132759599516SKenneth E. Jansen 
132859599516SKenneth E. Jansen     if (!FOUND)
132959599516SKenneth E. Jansen     {
133059599516SKenneth E. Jansen       //if(irank==0) printf("Warning readheader: Not found %s \n",keyphrase); //PhastaIOActiveFiles[i]->myrank is certainly initialized here.
133159599516SKenneth E. Jansen       if(PhastaIOActiveFiles[i]->myrank == 0) printf("WARNING readheader: Not found %s\n",keyphrase);
133259599516SKenneth E. Jansen       endTimer(&timer_end);
133359599516SKenneth E. Jansen       printPerf("readheader", timer_start, timer_end, 0, 0, "");
133459599516SKenneth E. Jansen       return;
133559599516SKenneth E. Jansen     }
133659599516SKenneth E. Jansen 
133759599516SKenneth E. Jansen     // Find the part we want ...
133859599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->read_part_count = PhastaIOActiveFiles[i]->GPid -
133959599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->myrank * PhastaIOActiveFiles[i]->nppp - 1;
134059599516SKenneth E. Jansen 
134159599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset =
134259599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->my_read_table[PhastaIOActiveFiles[i]->read_field_count][PhastaIOActiveFiles[i]->read_part_count];
134359599516SKenneth E. Jansen 
134459599516SKenneth E. Jansen     // 	  printf("****Rank %d offset is %d\n",PhastaIOActiveFiles[i]->myrank,PhastaIOActiveFiles[i]->my_offset);
134559599516SKenneth E. Jansen 
134659599516SKenneth E. Jansen     // Read each datablock header here ...
134759599516SKenneth E. Jansen 
134859599516SKenneth E. Jansen     MPI_File_read_at_all( PhastaIOActiveFiles[i]->file_handle,
134959599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->my_offset+1,
135059599516SKenneth E. Jansen         read_out_tag,
135159599516SKenneth E. Jansen         MAX_FIELDS_NAME_LENGTH-1,
135259599516SKenneth E. Jansen         MPI_CHAR,
135359599516SKenneth E. Jansen         &read_offset_status );
135459599516SKenneth E. Jansen     token = strtok ( read_out_tag, ":" );
135559599516SKenneth E. Jansen 
135659599516SKenneth E. Jansen     // 	  printf("&&&&Rank %d read_out_tag is %s\n",PhastaIOActiveFiles[i]->myrank,read_out_tag);
135759599516SKenneth E. Jansen 
135859599516SKenneth 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.
135959599516SKenneth E. Jansen     {
136059599516SKenneth E. Jansen       FOUND = true ;
136159599516SKenneth E. Jansen       token = strtok( NULL, " ,;<>" );
136259599516SKenneth E. Jansen       for( j=0; j < *nItems && ( token = strtok( NULL," ,;<>") ); j++ )
136359599516SKenneth E. Jansen         valueListInt[j] = atoi( token );
136459599516SKenneth E. Jansen 
136559599516SKenneth E. Jansen       if ( j < *nItems )
136659599516SKenneth E. Jansen       {
136759599516SKenneth E. Jansen         fprintf( stderr, "Expected # of ints not found for: %s\n", keyphrase );
136859599516SKenneth E. Jansen       }
136959599516SKenneth E. Jansen     }
137059599516SKenneth E. Jansen     else {
137159599516SKenneth E. Jansen       //if(irank==0)
137259599516SKenneth E. Jansen       if(PhastaIOActiveFiles[i]->myrank == 0)
137359599516SKenneth E. Jansen         // If we enter this if, there is a problem with the name of some fields
137459599516SKenneth E. Jansen       {
137559599516SKenneth E. Jansen         printf("Error readheader: Unexpected mismatch between keyphrase = %s and token = %s\n",keyphrase,token);
137659599516SKenneth E. Jansen       }
137759599516SKenneth E. Jansen     }
137859599516SKenneth E. Jansen   }
137959599516SKenneth E. Jansen 
138059599516SKenneth E. Jansen   endTimer(&timer_end);
138159599516SKenneth E. Jansen   printPerf("readheader", timer_start, timer_end, 0, 0, "");
138259599516SKenneth E. Jansen 
138359599516SKenneth E. Jansen }
138459599516SKenneth E. Jansen 
13853956dcfeSCameron Smith void readDataBlock(
13863956dcfeSCameron Smith     FILE* fileObject,
13873956dcfeSCameron Smith     void* valueArray,
13883956dcfeSCameron Smith     int nItems,
13893956dcfeSCameron Smith     const char  datatype[],
13903956dcfeSCameron Smith     const char  iotype[] )
13913956dcfeSCameron Smith {
13923956dcfeSCameron Smith   isBinary(iotype);
13933956dcfeSCameron Smith   size_t type_size = typeSize( datatype );
13943956dcfeSCameron Smith   if ( binary_format ) {
13953956dcfeSCameron Smith     char junk = '\0';
13963956dcfeSCameron Smith     fread( valueArray, type_size, nItems, fileObject );
13973956dcfeSCameron Smith     fread( &junk, sizeof(char), 1 , fileObject );
13983956dcfeSCameron Smith     if ( Wrong_Endian ) SwapArrayByteOrder( valueArray, type_size, nItems );
13993956dcfeSCameron Smith   } else {
14003956dcfeSCameron Smith     char* ts1 = StringStripper( datatype );
14013956dcfeSCameron Smith     if ( cscompare( "integer", ts1 ) ) {
14023956dcfeSCameron Smith       for( int n=0; n < nItems ; n++ )
14033956dcfeSCameron Smith         fscanf(fileObject, "%d\n",(int*)((int*)valueArray+n) );
14043956dcfeSCameron Smith     } else if ( cscompare( "double", ts1 ) ) {
14053956dcfeSCameron Smith       for( int n=0; n < nItems ; n++ )
14063956dcfeSCameron Smith         fscanf(fileObject, "%lf\n",(double*)((double*)valueArray+n) );
14073956dcfeSCameron Smith     }
14083956dcfeSCameron Smith     free (ts1);
14093956dcfeSCameron Smith   }
14103956dcfeSCameron Smith }
14113956dcfeSCameron Smith 
1412103be424SCameron Smith void readdatablock(
1413103be424SCameron Smith     int*  fileDescriptor,
141459599516SKenneth E. Jansen     const char keyphrase[],
141559599516SKenneth E. Jansen     void* valueArray,
141659599516SKenneth E. Jansen     int*  nItems,
141759599516SKenneth E. Jansen     const char  datatype[],
141859599516SKenneth E. Jansen     const char  iotype[] )
141959599516SKenneth E. Jansen {
142059599516SKenneth E. Jansen   //if(irank == 0) printf("entering readdatablock()\n");
14212dd307a1SCameron Smith   unsigned long data_size = 0;
142259599516SKenneth E. Jansen   double timer_start, timer_end;
142359599516SKenneth E. Jansen   startTimer(&timer_start);
142459599516SKenneth E. Jansen 
142559599516SKenneth E. Jansen   int i = *fileDescriptor;
142659599516SKenneth E. Jansen   checkFileDescriptor("readdatablock",&i);
142759599516SKenneth E. Jansen 
142859599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
142959599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
143059599516SKenneth E. Jansen     FILE* fileObject;
143159599516SKenneth E. Jansen 
143259599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
143359599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
143459599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before\n") ;
143559599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
143659599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
143759599516SKenneth E. Jansen       endTimer(&timer_end);
143859599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
143959599516SKenneth E. Jansen       return;
144059599516SKenneth E. Jansen     }
144159599516SKenneth E. Jansen 
144259599516SKenneth E. Jansen     // error check..
144359599516SKenneth E. Jansen     // since we require that a consistant header always preceed the data block
144459599516SKenneth E. Jansen     // let us check to see that it is actually the case.
144559599516SKenneth E. Jansen 
1446961a4ff6SCameron Smith     if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
144759599516SKenneth E. Jansen       fprintf(stderr, "Header not consistant with data block\n");
1448961a4ff6SCameron Smith       fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
144959599516SKenneth E. Jansen       fprintf(stderr, "DataBlock: %s\n ", keyphrase );
145059599516SKenneth E. Jansen       fprintf(stderr, "Please recheck read sequence \n");
145159599516SKenneth E. Jansen       if( Strict_Error ) {
145259599516SKenneth E. Jansen         fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
145359599516SKenneth E. Jansen         endTimer(&timer_end);
145459599516SKenneth E. Jansen         printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
145559599516SKenneth E. Jansen         return;
145659599516SKenneth E. Jansen       }
145759599516SKenneth E. Jansen     }
145859599516SKenneth E. Jansen 
145959599516SKenneth E. Jansen     if ( LastHeaderNotFound ) {
146059599516SKenneth E. Jansen       endTimer(&timer_end);
146159599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
146259599516SKenneth E. Jansen       return;
146359599516SKenneth E. Jansen     }
146459599516SKenneth E. Jansen     fileObject = fileArray[ filePtr ];
146559599516SKenneth E. Jansen     Wrong_Endian = byte_order[ filePtr ];
1466bae968b9SCameron Smith     LastHeaderKey.erase(filePtr);
14673956dcfeSCameron Smith     readDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
146859599516SKenneth E. Jansen 
146959599516SKenneth E. Jansen     //return;
147059599516SKenneth E. Jansen   }
147159599516SKenneth E. Jansen   else {
147259599516SKenneth E. Jansen     // 	  printf("read data block\n");
147359599516SKenneth E. Jansen     MPI_Status read_data_status;
147459599516SKenneth E. Jansen     size_t type_size = typeSize( datatype );
147559599516SKenneth E. Jansen     int nUnits = *nItems;
147659599516SKenneth E. Jansen     isBinary( iotype );
147759599516SKenneth E. Jansen 
147859599516SKenneth E. Jansen     // read datablock then
147959599516SKenneth E. Jansen     //MR CHANGE
148059599516SKenneth E. Jansen     //             if ( cscompare ( datatype, "double"))
148159599516SKenneth E. Jansen     char* ts2 = StringStripper( datatype );
148259599516SKenneth E. Jansen     if ( cscompare ( "double" , ts2))
148359599516SKenneth E. Jansen       //MR CHANGE END
148459599516SKenneth E. Jansen     {
148559599516SKenneth E. Jansen 
148659599516SKenneth E. Jansen       MPI_File_read_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
148759599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
148859599516SKenneth E. Jansen           valueArray,
148959599516SKenneth E. Jansen           nUnits,
149059599516SKenneth E. Jansen           MPI_DOUBLE );
149159599516SKenneth E. Jansen       MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
149259599516SKenneth E. Jansen           valueArray,
149359599516SKenneth E. Jansen           &read_data_status );
149459599516SKenneth E. Jansen       data_size=8*nUnits;
149559599516SKenneth E. Jansen 
149659599516SKenneth E. Jansen     }
149759599516SKenneth E. Jansen     //MR CHANGE
149859599516SKenneth E. Jansen     //             else if ( cscompare ( datatype, "integer"))
149959599516SKenneth E. Jansen     else if ( cscompare ( "integer" , ts2))
150059599516SKenneth E. Jansen       //MR CHANGE END
150159599516SKenneth E. Jansen     {
150259599516SKenneth E. Jansen       MPI_File_read_at_all_begin(PhastaIOActiveFiles[i]->file_handle,
150359599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
150459599516SKenneth E. Jansen           valueArray,
150559599516SKenneth E. Jansen           nUnits,
150659599516SKenneth E. Jansen           MPI_INT );
150759599516SKenneth E. Jansen       MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
150859599516SKenneth E. Jansen           valueArray,
150959599516SKenneth E. Jansen           &read_data_status );
151059599516SKenneth E. Jansen       data_size=4*nUnits;
151159599516SKenneth E. Jansen     }
151259599516SKenneth E. Jansen     else
151359599516SKenneth E. Jansen     {
151459599516SKenneth E. Jansen       *fileDescriptor = DATA_TYPE_ILLEGAL;
151559599516SKenneth E. Jansen       printf("readdatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
151659599516SKenneth E. Jansen       endTimer(&timer_end);
151759599516SKenneth E. Jansen       printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
151859599516SKenneth E. Jansen       return;
151959599516SKenneth E. Jansen     }
152059599516SKenneth E. Jansen     free(ts2);
152159599516SKenneth E. Jansen 
152259599516SKenneth E. Jansen 
152359599516SKenneth E. Jansen     // 	  printf("%d Read finishe\n",PhastaIOActiveFiles[i]->myrank);
152459599516SKenneth E. Jansen 
152559599516SKenneth E. Jansen     // Swap data byte order if endianess is different ...
152659599516SKenneth E. Jansen     if ( PhastaIOActiveFiles[i]->Wrong_Endian )
152759599516SKenneth E. Jansen     {
152859599516SKenneth E. Jansen       SwapArrayByteOrder( valueArray, type_size, nUnits );
152959599516SKenneth E. Jansen     }
153059599516SKenneth E. Jansen   }
153159599516SKenneth E. Jansen 
153259599516SKenneth E. Jansen   endTimer(&timer_end);
153359599516SKenneth E. Jansen   char extra_msg[1024];
153459599516SKenneth E. Jansen   memset(extra_msg, '\0', 1024);
153559599516SKenneth E. Jansen   char* key = StringStripper(keyphrase);
153659599516SKenneth E. Jansen   sprintf(extra_msg, " field is %s ", key);
153759599516SKenneth E. Jansen   printPerf("readdatablock", timer_start, timer_end, data_size, 1, extra_msg);
153859599516SKenneth E. Jansen   free(key);
153959599516SKenneth E. Jansen 
154059599516SKenneth E. Jansen }
154159599516SKenneth E. Jansen 
1542103be424SCameron Smith void writeHeader(
1543103be424SCameron Smith     FILE* f,
1544ea868eb1SCameron Smith     const char keyphrase[],
1545ea868eb1SCameron Smith     const void* valueArray,
1546ea868eb1SCameron Smith     const int nItems,
1547ea868eb1SCameron Smith     const int ndataItems,
1548ea868eb1SCameron Smith     const char datatype[],
1549ea868eb1SCameron Smith     const char iotype[])
1550ea868eb1SCameron Smith {
1551ea868eb1SCameron Smith   isBinary( iotype );
1552ea868eb1SCameron Smith 
1553ea868eb1SCameron Smith   const int _newline =
1554ea868eb1SCameron Smith     ( ndataItems > 0 ) ? sizeof( char ) : 0;
1555ea868eb1SCameron Smith   int size_of_nextblock =
1556ea868eb1SCameron Smith     ( binary_format ) ? typeSize(datatype) * ndataItems + _newline : ndataItems;
1557ea868eb1SCameron Smith 
1558ea868eb1SCameron Smith   fprintf( f, "%s : < %d > ", keyphrase, size_of_nextblock );
1559ea868eb1SCameron Smith   for( int i = 0; i < nItems; i++ )
1560ea868eb1SCameron Smith     fprintf( f, "%d ", *((int*)((int*)valueArray+i)));
1561ea868eb1SCameron Smith   fprintf( f, "\n");
1562ea868eb1SCameron Smith }
1563ea868eb1SCameron Smith 
1564103be424SCameron Smith void writeheader(
1565103be424SCameron Smith     const int* fileDescriptor,
156659599516SKenneth E. Jansen     const char keyphrase[],
156759599516SKenneth E. Jansen     const void* valueArray,
156859599516SKenneth E. Jansen     const int* nItems,
156959599516SKenneth E. Jansen     const int* ndataItems,
157059599516SKenneth E. Jansen     const char datatype[],
157159599516SKenneth E. Jansen     const char iotype[])
157259599516SKenneth E. Jansen {
157359599516SKenneth E. Jansen 
157459599516SKenneth E. Jansen   //if(irank == 0) printf("entering writeheader()\n");
157559599516SKenneth E. Jansen 
157659599516SKenneth E. Jansen   double timer_start, timer_end;
157759599516SKenneth E. Jansen   startTimer(&timer_start);
157859599516SKenneth E. Jansen 
157959599516SKenneth E. Jansen   int i = *fileDescriptor;
158059599516SKenneth E. Jansen   checkFileDescriptor("writeheader",&i);
158159599516SKenneth E. Jansen 
158259599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
158359599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
158459599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
158559599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
158659599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
158759599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
158859599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
158959599516SKenneth E. Jansen       endTimer(&timer_end);
159059599516SKenneth E. Jansen       printPerf("writeheader", timer_start, timer_end, 0, 0, "");
159159599516SKenneth E. Jansen       return;
159259599516SKenneth E. Jansen     }
159359599516SKenneth E. Jansen 
1594961a4ff6SCameron Smith     LastHeaderKey[filePtr] = std::string(keyphrase);
159559599516SKenneth E. Jansen     DataSize = *ndataItems;
1596ea868eb1SCameron Smith     FILE* fileObject = fileArray[ filePtr ] ;
1597ea868eb1SCameron Smith     header_type[ filePtr ] = typeSize( datatype );
1598ea868eb1SCameron Smith     writeHeader(fileObject,keyphrase,valueArray,*nItems,
1599ea868eb1SCameron Smith         *ndataItems,datatype,iotype);
160059599516SKenneth E. Jansen   }
160159599516SKenneth E. Jansen   else { // else it's parallel I/O
160259599516SKenneth E. Jansen     DataSize = *ndataItems;
160359599516SKenneth E. Jansen     size_t type_size = typeSize( datatype );
160459599516SKenneth E. Jansen     isBinary( iotype );
160559599516SKenneth E. Jansen     char mpi_tag[MAX_FIELDS_NAME_LENGTH];
160659599516SKenneth E. Jansen 
160759599516SKenneth E. Jansen     int string_length = strlen( keyphrase );
160859599516SKenneth E. Jansen     char* buffer = (char*) malloc ( string_length+1 );
160959599516SKenneth E. Jansen 
161059599516SKenneth E. Jansen     strcpy ( buffer, keyphrase);
161159599516SKenneth E. Jansen     buffer[ string_length ] = '\0';
161259599516SKenneth E. Jansen 
161359599516SKenneth E. Jansen     char* st2 = strtok ( buffer, "@" );
161459599516SKenneth E. Jansen     st2 = strtok (NULL, "@");
161559599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = atoi(st2);
161659599516SKenneth E. Jansen 
161759599516SKenneth E. Jansen     if ( char* p = strpbrk(buffer, "@") )
161859599516SKenneth E. Jansen       *p = '\0';
161959599516SKenneth E. Jansen 
162059599516SKenneth E. Jansen     bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
162159599516SKenneth E. Jansen     sprintf(mpi_tag, "\n%s : %d\n", buffer, PhastaIOActiveFiles[i]->field_count);
16222dd307a1SCameron Smith     unsigned long offset_value;
162359599516SKenneth E. Jansen 
162459599516SKenneth E. Jansen     int temp = *ndataItems;
16252dd307a1SCameron Smith     unsigned long number_of_items = (unsigned long)temp;
162659599516SKenneth E. Jansen     MPI_Barrier(PhastaIOActiveFiles[i]->local_comm);
162759599516SKenneth E. Jansen 
162859599516SKenneth E. Jansen     MPI_Scan( &number_of_items,
162959599516SKenneth E. Jansen         &offset_value,
163059599516SKenneth E. Jansen         1,
163159599516SKenneth E. Jansen         MPI_LONG_LONG_INT,
163259599516SKenneth E. Jansen         MPI_SUM,
163359599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_comm );
163459599516SKenneth E. Jansen 
163559599516SKenneth E. Jansen     offset_value = (offset_value - number_of_items) * type_size;
163659599516SKenneth E. Jansen 
163759599516SKenneth E. Jansen     offset_value += PhastaIOActiveFiles[i]->local_myrank *
163859599516SKenneth E. Jansen       DB_HEADER_SIZE +
163959599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->next_start_address;
164059599516SKenneth E. Jansen     // This offset is the starting address of each datablock header...
164159599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset = offset_value;
164259599516SKenneth E. Jansen 
164359599516SKenneth E. Jansen     // Write in my offset table ...
164459599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->my_offset_table[PhastaIOActiveFiles[i]->field_count][PhastaIOActiveFiles[i]->part_count] =
164559599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->my_offset;
164659599516SKenneth E. Jansen 
164759599516SKenneth E. Jansen     // Update the next-start-address ...
164859599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->next_start_address = offset_value +
164959599516SKenneth E. Jansen       number_of_items * type_size +
165059599516SKenneth E. Jansen       DB_HEADER_SIZE;
165159599516SKenneth E. Jansen     MPI_Bcast( &(PhastaIOActiveFiles[i]->next_start_address),
165259599516SKenneth E. Jansen         1,
165359599516SKenneth E. Jansen         MPI_LONG_LONG_INT,
165459599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_numprocs-1,
165559599516SKenneth E. Jansen         PhastaIOActiveFiles[i]->local_comm );
165659599516SKenneth E. Jansen 
165759599516SKenneth E. Jansen     // Prepare datablock header ...
165859599516SKenneth E. Jansen     int _newline = (*ndataItems>0)?sizeof(char):0;
165959599516SKenneth E. Jansen     unsigned int size_of_nextblock = type_size * (*ndataItems) + _newline;
166059599516SKenneth E. Jansen 
166159599516SKenneth E. Jansen     //char datablock_header[255];
166259599516SKenneth E. Jansen     //bzero((void*)datablock_header,255);
166359599516SKenneth E. Jansen     char datablock_header[DB_HEADER_SIZE];
166459599516SKenneth E. Jansen     bzero((void*)datablock_header,DB_HEADER_SIZE);
166559599516SKenneth E. Jansen 
166659599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->GPid = PhastaIOActiveFiles[i]->nppp*PhastaIOActiveFiles[i]->myrank+PhastaIOActiveFiles[i]->part_count;
166759599516SKenneth E. Jansen     sprintf( datablock_header,
166859599516SKenneth E. Jansen         "\n%s : < %u >",
166959599516SKenneth E. Jansen         keyphrase,
167059599516SKenneth E. Jansen         size_of_nextblock );
167159599516SKenneth E. Jansen 
167259599516SKenneth E. Jansen     for ( int j = 0; j < *nItems; j++ )
167359599516SKenneth E. Jansen     {
167459599516SKenneth E. Jansen       sprintf( datablock_header,
167559599516SKenneth E. Jansen           "%s %d ",
167659599516SKenneth E. Jansen           datablock_header,
167759599516SKenneth E. Jansen           *((int*)((int*)valueArray+j)));
167859599516SKenneth E. Jansen     }
167959599516SKenneth E. Jansen     sprintf( datablock_header,
168059599516SKenneth E. Jansen         "%s\n ",
168159599516SKenneth E. Jansen         datablock_header );
168259599516SKenneth E. Jansen 
168359599516SKenneth E. Jansen     // Write datablock header ...
168459599516SKenneth E. Jansen     //MR CHANGE
168559599516SKenneth E. Jansen     // 	if ( cscompare(datatype,"double") )
168659599516SKenneth E. Jansen     char* ts1 = StringStripper( datatype );
168759599516SKenneth E. Jansen     if ( cscompare("double",ts1) )
168859599516SKenneth E. Jansen       //MR CHANGE END
168959599516SKenneth E. Jansen     {
169059599516SKenneth E. Jansen       free ( PhastaIOActiveFiles[i]->double_chunk );
169159599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->double_chunk = ( double * )malloc( (sizeof( double )*number_of_items+ DB_HEADER_SIZE));
169259599516SKenneth E. Jansen 
169359599516SKenneth E. Jansen       double * aa = ( double * )datablock_header;
169459599516SKenneth E. Jansen       memcpy(PhastaIOActiveFiles[i]->double_chunk, aa, DB_HEADER_SIZE);
169559599516SKenneth E. Jansen     }
169659599516SKenneth E. Jansen     //MR CHANGE
169759599516SKenneth E. Jansen     // 	if  ( cscompare(datatype,"integer") )
169859599516SKenneth E. Jansen     else if ( cscompare("integer",ts1) )
169959599516SKenneth E. Jansen       //MR CHANGE END
170059599516SKenneth E. Jansen     {
170159599516SKenneth E. Jansen       free ( PhastaIOActiveFiles[i]->int_chunk );
170259599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->int_chunk = ( int * )malloc( (sizeof( int )*number_of_items+ DB_HEADER_SIZE));
170359599516SKenneth E. Jansen 
170459599516SKenneth E. Jansen       int * aa = ( int * )datablock_header;
170559599516SKenneth E. Jansen       memcpy(PhastaIOActiveFiles[i]->int_chunk, aa, DB_HEADER_SIZE);
170659599516SKenneth E. Jansen     }
170759599516SKenneth E. Jansen     else {
170859599516SKenneth E. Jansen       //             *fileDescriptor = DATA_TYPE_ILLEGAL;
170959599516SKenneth E. Jansen       printf("writeheader - DATA_TYPE_ILLEGAL - %s\n",datatype);
171059599516SKenneth E. Jansen       endTimer(&timer_end);
171159599516SKenneth E. Jansen       printPerf("writeheader", timer_start, timer_end, 0, 0, "");
171259599516SKenneth E. Jansen       return;
171359599516SKenneth E. Jansen     }
171459599516SKenneth E. Jansen     free(ts1);
171559599516SKenneth E. Jansen 
171659599516SKenneth E. Jansen     PhastaIOActiveFiles[i]->part_count++;
171759599516SKenneth E. Jansen     if ( PhastaIOActiveFiles[i]->part_count == PhastaIOActiveFiles[i]->nppp ) {
171859599516SKenneth E. Jansen       //A new field will be written
171959599516SKenneth E. Jansen       if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
172059599516SKenneth E. Jansen         memcpy( PhastaIOActiveFiles[i]->master_header +
172159599516SKenneth E. Jansen             PhastaIOActiveFiles[i]->field_count *
172259599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH +
172359599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH * 2,
172459599516SKenneth E. Jansen             mpi_tag,
172559599516SKenneth E. Jansen             MAX_FIELDS_NAME_LENGTH-1);
172659599516SKenneth E. Jansen       }
172759599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->field_count++;
172859599516SKenneth E. Jansen       PhastaIOActiveFiles[i]->part_count=0;
172959599516SKenneth E. Jansen     }
173059599516SKenneth E. Jansen     free(buffer);
173159599516SKenneth E. Jansen   }
173259599516SKenneth E. Jansen 
173359599516SKenneth E. Jansen   endTimer(&timer_end);
173459599516SKenneth E. Jansen   printPerf("writeheader", timer_start, timer_end, 0, 0, "");
173559599516SKenneth E. Jansen }
173659599516SKenneth E. Jansen 
1737103be424SCameron Smith void writeDataBlock(
1738103be424SCameron Smith     FILE* f,
1739ea868eb1SCameron Smith     const void* valueArray,
1740ea868eb1SCameron Smith     const int   nItems,
1741ea868eb1SCameron Smith     const char  datatype[],
1742ea868eb1SCameron Smith     const char  iotype[]  )
1743ea868eb1SCameron Smith {
1744ea868eb1SCameron Smith   isBinary( iotype );
1745ea868eb1SCameron Smith   size_t type_size = typeSize( datatype );
1746ea868eb1SCameron Smith   if ( binary_format ) {
1747ea868eb1SCameron Smith     fwrite( valueArray, type_size, nItems, f );
1748ea868eb1SCameron Smith     fprintf( f,"\n");
1749ea868eb1SCameron Smith   } else {
1750ea868eb1SCameron Smith     char* ts1 = StringStripper( datatype );
1751ea868eb1SCameron Smith     if ( cscompare( "integer", ts1 ) ) {
1752*be3da47bSCameron Smith       const int* vals = (int*) valueArray;
1753ea868eb1SCameron Smith       for( int n=0; n < nItems ; n++ )
1754*be3da47bSCameron Smith         fprintf(f,"%d\n",vals[n]);
1755ea868eb1SCameron Smith     } else if ( cscompare( "double", ts1 ) ) {
1756*be3da47bSCameron Smith       const double* vals = (double*) valueArray;
1757ea868eb1SCameron Smith       for( int n=0; n < nItems ; n++ )
1758*be3da47bSCameron Smith         fprintf(f,"%f\n",vals[n]);
1759ea868eb1SCameron Smith     }
1760ea868eb1SCameron Smith     free (ts1);
1761ea868eb1SCameron Smith   }
1762ea868eb1SCameron Smith }
1763ea868eb1SCameron Smith 
1764103be424SCameron Smith void writedatablock(
1765103be424SCameron Smith     const int* fileDescriptor,
176659599516SKenneth E. Jansen     const char keyphrase[],
176759599516SKenneth E. Jansen     const void* valueArray,
176859599516SKenneth E. Jansen     const int* nItems,
176959599516SKenneth E. Jansen     const char datatype[],
177059599516SKenneth E. Jansen     const char iotype[] )
177159599516SKenneth E. Jansen {
177259599516SKenneth E. Jansen   //if(irank == 0) printf("entering writedatablock()\n");
177359599516SKenneth E. Jansen 
17742dd307a1SCameron Smith   unsigned long data_size = 0;
177559599516SKenneth E. Jansen   double timer_start, timer_end;
177659599516SKenneth E. Jansen   startTimer(&timer_start);
177759599516SKenneth E. Jansen 
177859599516SKenneth E. Jansen   int i = *fileDescriptor;
177959599516SKenneth E. Jansen   checkFileDescriptor("writedatablock",&i);
178059599516SKenneth E. Jansen 
178159599516SKenneth E. Jansen   if ( PhastaIONextActiveIndex == 0 ) {
178259599516SKenneth E. Jansen     int filePtr = *fileDescriptor - 1;
178359599516SKenneth E. Jansen 
178459599516SKenneth E. Jansen     if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
178559599516SKenneth E. Jansen       fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
178659599516SKenneth E. Jansen       fprintf(stderr,"openfile_ function has to be called before \n") ;
178759599516SKenneth E. Jansen       fprintf(stderr,"acessing the file\n ") ;
178859599516SKenneth E. Jansen       fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
178959599516SKenneth E. Jansen       endTimer(&timer_end);
179059599516SKenneth E. Jansen       printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
179159599516SKenneth E. Jansen       return;
179259599516SKenneth E. Jansen     }
179359599516SKenneth E. Jansen     // since we require that a consistant header always preceed the data block
179459599516SKenneth E. Jansen     // let us check to see that it is actually the case.
179559599516SKenneth E. Jansen 
1796961a4ff6SCameron Smith     if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
179759599516SKenneth E. Jansen       fprintf(stderr, "Header not consistant with data block\n");
1798961a4ff6SCameron Smith       fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
179959599516SKenneth E. Jansen       fprintf(stderr, "DataBlock: %s\n ", keyphrase );
180059599516SKenneth E. Jansen       fprintf(stderr, "Please recheck write sequence \n");
180159599516SKenneth E. Jansen       if( Strict_Error ) {
180259599516SKenneth E. Jansen         fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
180359599516SKenneth E. Jansen         endTimer(&timer_end);
180459599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
180559599516SKenneth E. Jansen         return;
180659599516SKenneth E. Jansen       }
180759599516SKenneth E. Jansen     }
180859599516SKenneth E. Jansen 
180959599516SKenneth E. Jansen     FILE* fileObject =  fileArray[ filePtr ] ;
181059599516SKenneth E. Jansen     size_t type_size=typeSize( datatype );
181159599516SKenneth E. Jansen     isBinary( iotype );
181259599516SKenneth E. Jansen 
1813bae968b9SCameron Smith     LastHeaderKey.erase(filePtr);
1814bae968b9SCameron Smith 
181559599516SKenneth E. Jansen     if ( header_type[filePtr] != (int)type_size ) {
181659599516SKenneth E. Jansen       fprintf(stderr,"header and datablock differ on typeof data in the block for\n");
181759599516SKenneth E. Jansen       fprintf(stderr,"keyphrase : %s\n", keyphrase);
181859599516SKenneth E. Jansen       if( Strict_Error ) {
181959599516SKenneth E. Jansen         fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
182059599516SKenneth E. Jansen         endTimer(&timer_end);
182159599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
182259599516SKenneth E. Jansen         return;
182359599516SKenneth E. Jansen       }
182459599516SKenneth E. Jansen     }
182559599516SKenneth E. Jansen 
182659599516SKenneth E. Jansen     int nUnits = *nItems;
182759599516SKenneth E. Jansen 
182859599516SKenneth E. Jansen     if ( nUnits != DataSize ) {
182959599516SKenneth E. Jansen       fprintf(stderr,"header and datablock differ on number of data items for\n");
183059599516SKenneth E. Jansen       fprintf(stderr,"keyphrase : %s\n", keyphrase);
183159599516SKenneth E. Jansen       if( Strict_Error ) {
183259599516SKenneth E. Jansen         fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
183359599516SKenneth E. Jansen         endTimer(&timer_end);
183459599516SKenneth E. Jansen         printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
183559599516SKenneth E. Jansen         return;
183659599516SKenneth E. Jansen       }
183759599516SKenneth E. Jansen     }
1838ea868eb1SCameron Smith     writeDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
183959599516SKenneth E. Jansen   }
184059599516SKenneth E. Jansen   else {  // syncIO case
184159599516SKenneth E. Jansen     MPI_Status write_data_status;
184259599516SKenneth E. Jansen     isBinary( iotype );
184359599516SKenneth E. Jansen     int nUnits = *nItems;
184459599516SKenneth E. Jansen 
184559599516SKenneth E. Jansen     //MR CHANGE
184659599516SKenneth E. Jansen     // 	if ( cscompare(datatype,"double") )
184759599516SKenneth E. Jansen     char* ts1 = StringStripper( datatype );
184859599516SKenneth E. Jansen     if ( cscompare("double",ts1) )
184959599516SKenneth E. Jansen       //MR CHANGE END
185059599516SKenneth E. Jansen     {
185159599516SKenneth E. Jansen       memcpy((PhastaIOActiveFiles[i]->double_chunk+DB_HEADER_SIZE/sizeof(double)), valueArray, nUnits*sizeof(double));
185259599516SKenneth E. Jansen       MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
185359599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset,
185459599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->double_chunk,
185559599516SKenneth E. Jansen           //BLOCK_SIZE/sizeof(double),
185659599516SKenneth E. Jansen           nUnits+DB_HEADER_SIZE/sizeof(double),
185759599516SKenneth E. Jansen           MPI_DOUBLE );
185859599516SKenneth E. Jansen       MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
185959599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->double_chunk,
186059599516SKenneth E. Jansen           &write_data_status );
186159599516SKenneth E. Jansen       data_size=8*nUnits;
186259599516SKenneth E. Jansen     }
186359599516SKenneth E. Jansen     //MR CHANGE
186459599516SKenneth E. Jansen     // 	else if ( cscompare ( datatype, "integer"))
186559599516SKenneth E. Jansen     else if ( cscompare("integer",ts1) )
186659599516SKenneth E. Jansen       //MR CHANGE END
186759599516SKenneth E. Jansen     {
186859599516SKenneth E. Jansen       memcpy((PhastaIOActiveFiles[i]->int_chunk+DB_HEADER_SIZE/sizeof(int)), valueArray, nUnits*sizeof(int));
186959599516SKenneth E. Jansen       MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
187059599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->my_offset,
187159599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->int_chunk,
187259599516SKenneth E. Jansen           nUnits+DB_HEADER_SIZE/sizeof(int),
187359599516SKenneth E. Jansen           MPI_INT );
187459599516SKenneth E. Jansen       MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
187559599516SKenneth E. Jansen           PhastaIOActiveFiles[i]->int_chunk,
187659599516SKenneth E. Jansen           &write_data_status );
187759599516SKenneth E. Jansen       data_size=4*nUnits;
187859599516SKenneth E. Jansen     }
187959599516SKenneth E. Jansen     else {
188059599516SKenneth E. Jansen       printf("Error: writedatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
188159599516SKenneth E. Jansen       endTimer(&timer_end);
188259599516SKenneth E. Jansen       printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
188359599516SKenneth E. Jansen       return;
188459599516SKenneth E. Jansen     }
188559599516SKenneth E. Jansen     free(ts1);
188659599516SKenneth E. Jansen   }
188759599516SKenneth E. Jansen 
188859599516SKenneth E. Jansen   endTimer(&timer_end);
188959599516SKenneth E. Jansen   char extra_msg[1024];
189059599516SKenneth E. Jansen   memset(extra_msg, '\0', 1024);
189159599516SKenneth E. Jansen   char* key = StringStripper(keyphrase);
189259599516SKenneth E. Jansen   sprintf(extra_msg, " field is %s ", key);
189359599516SKenneth E. Jansen   printPerf("writedatablock", timer_start, timer_end, data_size, 1, extra_msg);
189459599516SKenneth E. Jansen   free(key);
189559599516SKenneth E. Jansen 
189659599516SKenneth E. Jansen }
189759599516SKenneth E. Jansen 
1898103be424SCameron Smith void SwapArrayByteOrder( void* array, int nbytes, int nItems )
189959599516SKenneth E. Jansen {
190059599516SKenneth E. Jansen   /* This swaps the byte order for the array of nItems each
190159599516SKenneth E. Jansen      of size nbytes , This will be called only locally  */
190259599516SKenneth E. Jansen   int i,j;
190359599516SKenneth E. Jansen   unsigned char* ucDst = (unsigned char*)array;
190459599516SKenneth E. Jansen   for(i=0; i < nItems; i++) {
190559599516SKenneth E. Jansen     for(j=0; j < (nbytes/2); j++)
190659599516SKenneth E. Jansen       std::swap( ucDst[j] , ucDst[(nbytes - 1) - j] );
190759599516SKenneth E. Jansen     ucDst += nbytes;
190859599516SKenneth E. Jansen   }
190959599516SKenneth E. Jansen }
191059599516SKenneth E. Jansen 
1911103be424SCameron Smith void writestring( int* fileDescriptor, const char inString[] )
191259599516SKenneth E. Jansen {
191359599516SKenneth E. Jansen   int filePtr = *fileDescriptor - 1;
191459599516SKenneth E. Jansen   FILE* fileObject = fileArray[filePtr] ;
191559599516SKenneth E. Jansen   fprintf(fileObject,"%s",inString );
191659599516SKenneth E. Jansen   return;
191759599516SKenneth E. Jansen }
191859599516SKenneth E. Jansen 
1919103be424SCameron Smith void Gather_Headers( int* fileDescriptor, vector< string >& headers )
192059599516SKenneth E. Jansen {
192159599516SKenneth E. Jansen   FILE* fileObject;
192259599516SKenneth E. Jansen   char Line[1024];
192359599516SKenneth E. Jansen 
192459599516SKenneth E. Jansen   fileObject = fileArray[ (*fileDescriptor)-1 ];
192559599516SKenneth E. Jansen 
192659599516SKenneth E. Jansen   while( !feof(fileObject) ) {
192759599516SKenneth E. Jansen     fgets( Line, 1024, fileObject);
192859599516SKenneth E. Jansen     if ( Line[0] == '#' ) {
192959599516SKenneth E. Jansen       headers.push_back( Line );
193059599516SKenneth E. Jansen     } else {
193159599516SKenneth E. Jansen       break;
193259599516SKenneth E. Jansen     }
193359599516SKenneth E. Jansen   }
193459599516SKenneth E. Jansen   rewind( fileObject );
193559599516SKenneth E. Jansen   clearerr( fileObject );
193659599516SKenneth E. Jansen }
193759599516SKenneth E. Jansen 
1938103be424SCameron Smith void isWrong( void ) {
1939103be424SCameron Smith   (Wrong_Endian) ? fprintf(stdout,"YES\n") : fprintf(stdout,"NO\n");
1940103be424SCameron Smith }
194159599516SKenneth E. Jansen 
1942103be424SCameron Smith void togglestrictmode( void ) { Strict_Error = !Strict_Error; }
1943103be424SCameron Smith 
1944103be424SCameron Smith int isLittleEndian( void )
194559599516SKenneth E. Jansen {
194659599516SKenneth E. Jansen   // this function returns a 1 if the current running architecture is
194759599516SKenneth E. Jansen   // LittleEndian Byte Ordered, else it returns a zero
194859599516SKenneth E. Jansen   union  {
194959599516SKenneth E. Jansen     long a;
195059599516SKenneth E. Jansen     char c[sizeof( long )];
195159599516SKenneth E. Jansen   } endianUnion;
195259599516SKenneth E. Jansen   endianUnion.a = 1 ;
195359599516SKenneth E. Jansen   if ( endianUnion.c[sizeof(long)-1] != 1 ) return 1 ;
195459599516SKenneth E. Jansen   else return 0;
195559599516SKenneth E. Jansen }
195659599516SKenneth E. Jansen 
195759599516SKenneth E. Jansen namespace PHASTA {
195859599516SKenneth E. Jansen   const char* const PhastaIO_traits<int>::type_string = "integer";
195959599516SKenneth E. Jansen   const char* const PhastaIO_traits<double>::type_string = "double";
196059599516SKenneth E. Jansen }
1961