159599516SKenneth E. Jansen /*
259599516SKenneth E. Jansen *
359599516SKenneth E. Jansen * This is the SyncIO library that uses MPI-IO collective functions to
459599516SKenneth E. Jansen * implement a flexible I/O checkpoint solution for a large number of
559599516SKenneth E. Jansen * processors.
659599516SKenneth E. Jansen *
759599516SKenneth E. Jansen * Previous developer: Ning Liu (liun2@cs.rpi.edu)
859599516SKenneth E. Jansen * Jing Fu (fuj@cs.rpi.edu)
959599516SKenneth E. Jansen * Current developers: Michel Rasquin (Michel.Rasquin@colorado.edu),
1059599516SKenneth E. Jansen * Ben Matthews (benjamin.a.matthews@colorado.edu)
1159599516SKenneth E. Jansen *
1259599516SKenneth E. Jansen */
1359599516SKenneth E. Jansen
1459599516SKenneth E. Jansen #include <map>
1559599516SKenneth E. Jansen #include <vector>
1659599516SKenneth E. Jansen #include <string>
17a746c838SCameron Smith #include <cstdarg>
1859599516SKenneth E. Jansen #include <string.h>
1959599516SKenneth E. Jansen #include <ctype.h>
2059599516SKenneth E. Jansen #include <stdlib.h>
2159599516SKenneth E. Jansen #include <stdio.h>
2259599516SKenneth E. Jansen #include <math.h>
23d3337298SCameron Smith #include <sstream>
2459599516SKenneth E. Jansen #include "phastaIO.h"
258f9016f6SCameron Smith #include "phiotmrc.h"
264600f8abSCameron Smith #include "phiompi.h"
2759599516SKenneth E. Jansen #include "mpi.h"
28f42e0444SCameron Smith #include "phiostats.h"
298f9016f6SCameron Smith #include "phiotimer.h"
3097a07b0aSCameron Smith #include <assert.h>
3159599516SKenneth E. Jansen
324600f8abSCameron Smith /* OS-specific things try to stay here */
334600f8abSCameron Smith #include <sys/stat.h>
344600f8abSCameron Smith #include <unistd.h>
354600f8abSCameron Smith #include <errno.h>
364600f8abSCameron Smith
374600f8abSCameron Smith
3859599516SKenneth E. Jansen #define VERSION_INFO_HEADER_SIZE 8192
3959599516SKenneth E. Jansen #define DB_HEADER_SIZE 1024
4059599516SKenneth E. Jansen #define ONE_MEGABYTE 1048576
4159599516SKenneth E. Jansen #define TWO_MEGABYTE 2097152
4259599516SKenneth E. Jansen #define ENDIAN_TEST_NUMBER 12180 // Troy's Zip Code!!
4359599516SKenneth E. Jansen #define MAX_PHASTA_FILES 64
4459599516SKenneth E. Jansen #define MAX_PHASTA_FILE_NAME_LENGTH 1024
4559599516SKenneth 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
4659599516SKenneth E. Jansen // -3 for MPI_IO_Tag, nFields and nppf, -4 for extra security (former nFiles)
4759599516SKenneth E. Jansen #define MAX_FIELDS_NAME_LENGTH 128
4859599516SKenneth E. Jansen #define DefaultMHSize (4*ONE_MEGABYTE)
4959599516SKenneth E. Jansen //#define DefaultMHSize (8350) //For test
5059599516SKenneth E. Jansen #define LATEST_WRITE_VERSION 1
5159599516SKenneth E. Jansen #define inv1024sq 953.674316406e-9 // = 1/1024/1024
5259599516SKenneth E. Jansen int MasterHeaderSize = -1;
5359599516SKenneth E. Jansen
5459599516SKenneth E. Jansen bool PRINT_PERF = false; // default print no perf results
5559599516SKenneth E. Jansen int irank = -1; // global rank, should never be manually manipulated
5659599516SKenneth E. Jansen int mysize = -1;
5759599516SKenneth E. Jansen
5859599516SKenneth E. Jansen // Static variables are bad but used here to store the subcommunicator and associated variables
5959599516SKenneth 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)
6059599516SKenneth E. Jansen static int s_assign_local_comm = 0;
6159599516SKenneth E. Jansen static MPI_Comm s_local_comm;
6259599516SKenneth E. Jansen static int s_local_size = -1;
6359599516SKenneth E. Jansen static int s_local_rank = -1;
6459599516SKenneth E. Jansen
6559599516SKenneth E. Jansen // the following defines are for debug printf
6659599516SKenneth E. Jansen #define PHASTAIO_DEBUG 0 //default to not print any debugging info
6759599516SKenneth E. Jansen
phprintf(const char * fmt,...)68fe03b87fSCameron Smith void phprintf(const char* fmt, ...) {
69be3da47bSCameron Smith (void)fmt;
7059599516SKenneth E. Jansen #if PHASTAIO_DEBUG
71fe03b87fSCameron Smith char format[1024];
72fe03b87fSCameron Smith snprintf(format, sizeof(format), "phastaIO debug: %s", fmt);
73fe03b87fSCameron Smith va_list ap;
74fe03b87fSCameron Smith va_start(ap,fmt);
75a746c838SCameron Smith vprintf(format,ap);
76fe03b87fSCameron Smith va_end(ap);
7759599516SKenneth E. Jansen #endif
78fe03b87fSCameron Smith }
79fe03b87fSCameron Smith
phprintf_0(const char * fmt,...)80fe03b87fSCameron Smith void phprintf_0(const char* fmt, ...) {
81be3da47bSCameron Smith (void)fmt;
82fe03b87fSCameron Smith #if PHASTAIO_DEBUG
83fe03b87fSCameron Smith int rank = 0;
84fe03b87fSCameron Smith MPI_Comm_rank(MPI_COMM_WORLD, &rank);
85fe03b87fSCameron Smith if(rank == 0){
86fe03b87fSCameron Smith char format[1024];
87fe03b87fSCameron Smith snprintf(format, sizeof(format), "phastaIO debug: irank=0 %s", fmt);
88fe03b87fSCameron Smith va_list ap;
89a746c838SCameron Smith va_start(ap,fmt);
90fe03b87fSCameron Smith vprintf(format, ap);
91fe03b87fSCameron Smith va_end(ap);
92fe03b87fSCameron Smith }
93fe03b87fSCameron Smith #endif
94fe03b87fSCameron Smith }
9559599516SKenneth E. Jansen
9659599516SKenneth E. Jansen enum PhastaIO_Errors
9759599516SKenneth E. Jansen {
9859599516SKenneth E. Jansen MAX_PHASTA_FILES_EXCEEDED = -1,
9959599516SKenneth E. Jansen UNABLE_TO_OPEN_FILE = -2,
10059599516SKenneth E. Jansen NOT_A_MPI_FILE = -3,
10159599516SKenneth E. Jansen GPID_EXCEEDED = -4,
10210d56689SCameron Smith DATA_TYPE_ILLEGAL = -5
10359599516SKenneth E. Jansen };
10459599516SKenneth E. Jansen
10559599516SKenneth E. Jansen using namespace std;
10659599516SKenneth E. Jansen
10759599516SKenneth E. Jansen namespace{
10859599516SKenneth E. Jansen
109961a4ff6SCameron Smith map<int, std::string> LastHeaderKey;
11059599516SKenneth E. Jansen vector< FILE* > fileArray;
11159599516SKenneth E. Jansen vector< bool > byte_order;
11259599516SKenneth E. Jansen vector< int > header_type;
11359599516SKenneth E. Jansen int DataSize=0;
11459599516SKenneth E. Jansen bool LastHeaderNotFound = false;
11559599516SKenneth E. Jansen bool Wrong_Endian = false ;
11659599516SKenneth E. Jansen bool Strict_Error = false ;
11759599516SKenneth E. Jansen bool binary_format = true;
11859599516SKenneth E. Jansen
11959599516SKenneth E. Jansen /***********************************************************************/
12059599516SKenneth E. Jansen /***************** NEW PHASTA IO CODE STARTS HERE **********************/
12159599516SKenneth E. Jansen /***********************************************************************/
12259599516SKenneth E. Jansen
12359599516SKenneth E. Jansen typedef struct
12459599516SKenneth E. Jansen {
12559599516SKenneth E. Jansen char filename[MAX_PHASTA_FILE_NAME_LENGTH]; /* defafults to 1024 */
1262dd307a1SCameron Smith unsigned long my_offset;
1272dd307a1SCameron Smith unsigned long next_start_address;
1282dd307a1SCameron Smith unsigned long **my_offset_table;
1292dd307a1SCameron Smith unsigned long **my_read_table;
13059599516SKenneth E. Jansen
13159599516SKenneth E. Jansen double * double_chunk;
13259599516SKenneth E. Jansen double * read_double_chunk;
13359599516SKenneth E. Jansen
13459599516SKenneth E. Jansen int field_count;
13559599516SKenneth E. Jansen int part_count;
13659599516SKenneth E. Jansen int read_field_count;
13759599516SKenneth E. Jansen int read_part_count;
13859599516SKenneth E. Jansen int GPid;
13959599516SKenneth E. Jansen int start_id;
14059599516SKenneth E. Jansen
14159599516SKenneth E. Jansen int mhsize;
14259599516SKenneth E. Jansen
14359599516SKenneth E. Jansen int myrank;
14459599516SKenneth E. Jansen int numprocs;
14559599516SKenneth E. Jansen int local_myrank;
14659599516SKenneth E. Jansen int local_numprocs;
14759599516SKenneth E. Jansen
14859599516SKenneth E. Jansen int nppp;
14959599516SKenneth E. Jansen int nPPF;
15059599516SKenneth E. Jansen int nFiles;
15159599516SKenneth E. Jansen int nFields;
15259599516SKenneth E. Jansen
15359599516SKenneth E. Jansen int * int_chunk;
15459599516SKenneth E. Jansen int * read_int_chunk;
15559599516SKenneth E. Jansen
15659599516SKenneth E. Jansen int Wrong_Endian; /* default to false */
15759599516SKenneth E. Jansen char * master_header;
15859599516SKenneth E. Jansen MPI_File file_handle;
15959599516SKenneth E. Jansen MPI_Comm local_comm;
16059599516SKenneth E. Jansen } phastaio_file_t;
16159599516SKenneth E. Jansen
16259599516SKenneth E. Jansen typedef struct
16359599516SKenneth E. Jansen {
16459599516SKenneth E. Jansen int nppf, nfields;
16559599516SKenneth E. Jansen char * masterHeader;
16659599516SKenneth E. Jansen }serial_file;
16759599516SKenneth E. Jansen
16859599516SKenneth E. Jansen serial_file *SerialFile;
16959599516SKenneth E. Jansen phastaio_file_t *PhastaIOActiveFiles[MAX_PHASTA_FILES];
17059599516SKenneth E. Jansen int PhastaIONextActiveIndex = 0; /* indicates next index to allocate */
17159599516SKenneth E. Jansen
17259599516SKenneth E. Jansen // the caller has the responsibility to delete the returned string
17359599516SKenneth E. Jansen // TODO: StringStipper("nbc value? ") returns NULL?
StringStripper(const char istring[])174103be424SCameron Smith char* StringStripper( const char istring[] ) {
17559599516SKenneth E. Jansen int length = strlen( istring );
17659599516SKenneth E. Jansen char* dest = (char *)malloc( length + 1 );
17759599516SKenneth E. Jansen strcpy( dest, istring );
17859599516SKenneth E. Jansen dest[ length ] = '\0';
17959599516SKenneth E. Jansen if ( char* p = strpbrk( dest, " ") )
18059599516SKenneth E. Jansen *p = '\0';
18159599516SKenneth E. Jansen return dest;
18259599516SKenneth E. Jansen }
18359599516SKenneth E. Jansen
cscompare(const char teststring[],const char targetstring[])184103be424SCameron Smith inline int cscompare( const char teststring[], const char targetstring[] ) {
18559599516SKenneth E. Jansen char* s1 = const_cast<char*>(teststring);
18659599516SKenneth E. Jansen char* s2 = const_cast<char*>(targetstring);
18759599516SKenneth E. Jansen
18859599516SKenneth E. Jansen while( *s1 == ' ') s1++;
18959599516SKenneth E. Jansen while( *s2 == ' ') s2++;
19059599516SKenneth E. Jansen while( ( *s1 )
19159599516SKenneth E. Jansen && ( *s2 )
19259599516SKenneth E. Jansen && ( *s2 != '?')
19359599516SKenneth E. Jansen && ( tolower( *s1 )==tolower( *s2 ) ) ) {
19459599516SKenneth E. Jansen s1++;
19559599516SKenneth E. Jansen s2++;
19659599516SKenneth E. Jansen while( *s1 == ' ') s1++;
19759599516SKenneth E. Jansen while( *s2 == ' ') s2++;
19859599516SKenneth E. Jansen }
19959599516SKenneth E. Jansen if ( !( *s1 ) || ( *s1 == '?') ) return 1;
20059599516SKenneth E. Jansen else return 0;
20159599516SKenneth E. Jansen }
20259599516SKenneth E. Jansen
isBinary(const char iotype[])203103be424SCameron Smith inline void isBinary( const char iotype[] ) {
20459599516SKenneth E. Jansen char* fname = StringStripper( iotype );
20559599516SKenneth E. Jansen if ( cscompare( fname, "binary" ) ) binary_format = true;
20659599516SKenneth E. Jansen else binary_format = false;
20759599516SKenneth E. Jansen free (fname);
20859599516SKenneth E. Jansen
20959599516SKenneth E. Jansen }
21059599516SKenneth E. Jansen
typeSize(const char typestring[])211103be424SCameron Smith inline size_t typeSize( const char typestring[] ) {
21259599516SKenneth E. Jansen char* ts1 = StringStripper( typestring );
21359599516SKenneth E. Jansen if ( cscompare( "integer", ts1 ) ) {
21459599516SKenneth E. Jansen free (ts1);
21559599516SKenneth E. Jansen return sizeof(int);
21659599516SKenneth E. Jansen } else if ( cscompare( "double", ts1 ) ) {
21759599516SKenneth E. Jansen free (ts1);
21859599516SKenneth E. Jansen return sizeof( double );
21959599516SKenneth E. Jansen } else {
22059599516SKenneth E. Jansen free (ts1);
22159599516SKenneth E. Jansen fprintf(stderr,"unknown type : %s\n",ts1);
22259599516SKenneth E. Jansen return 0;
22359599516SKenneth E. Jansen }
22459599516SKenneth E. Jansen }
22559599516SKenneth E. Jansen
readHeader(FILE * fileObject,const char phrase[],int * params,int expect)226103be424SCameron Smith int readHeader(
227103be424SCameron Smith FILE* fileObject,
22859599516SKenneth E. Jansen const char phrase[],
22959599516SKenneth E. Jansen int* params,
23059599516SKenneth E. Jansen int expect ) {
23159599516SKenneth E. Jansen char* text_header;
23259599516SKenneth E. Jansen char* token;
23397a07b0aSCameron Smith char Line[1024] = "\0";
23459599516SKenneth E. Jansen char junk;
23559599516SKenneth E. Jansen bool FOUND = false ;
23659599516SKenneth E. Jansen int real_length;
23759599516SKenneth E. Jansen int skip_size, integer_value;
23859599516SKenneth E. Jansen int rewind_count=0;
23959599516SKenneth E. Jansen
24059599516SKenneth E. Jansen if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
24159599516SKenneth E. Jansen rewind( fileObject );
24259599516SKenneth E. Jansen clearerr( fileObject );
24359599516SKenneth E. Jansen rewind_count++;
24459599516SKenneth E. Jansen fgets( Line, 1024, fileObject );
24559599516SKenneth E. Jansen }
24659599516SKenneth E. Jansen
24759599516SKenneth E. Jansen while( !FOUND && ( rewind_count < 2 ) ) {
24859599516SKenneth E. Jansen if ( ( Line[0] != '\n' ) && ( real_length = strcspn( Line, "#" )) ) {
24959599516SKenneth E. Jansen text_header = (char *)malloc( real_length + 1 );
25059599516SKenneth E. Jansen strncpy( text_header, Line, real_length );
25159599516SKenneth E. Jansen text_header[ real_length ] =static_cast<char>(NULL);
25259599516SKenneth E. Jansen token = strtok ( text_header, ":" );
25397a07b0aSCameron Smith assert(token);
25497a07b0aSCameron Smith if( cscompare( phrase , token ) ) {
25559599516SKenneth E. Jansen FOUND = true ;
25659599516SKenneth E. Jansen token = strtok( NULL, " ,;<>" );
25797a07b0aSCameron Smith assert(token);
25859599516SKenneth E. Jansen skip_size = atoi( token );
25959599516SKenneth E. Jansen int i;
26059599516SKenneth E. Jansen for( i=0; i < expect && ( token = strtok( NULL," ,;<>") ); i++) {
26197a07b0aSCameron Smith assert(token);
26259599516SKenneth E. Jansen params[i] = atoi( token );
26359599516SKenneth E. Jansen }
26459599516SKenneth E. Jansen if ( i < expect ) {
265103be424SCameron Smith fprintf(stderr,
266103be424SCameron Smith "Aloha Expected # of ints not found for: %s\n",phrase );
26759599516SKenneth E. Jansen }
26859599516SKenneth E. Jansen } else if ( cscompare(token,"byteorder magic number") ) {
26959599516SKenneth E. Jansen if ( binary_format ) {
27059599516SKenneth E. Jansen fread((void*)&integer_value,sizeof(int),1,fileObject);
27159599516SKenneth E. Jansen fread( &junk, sizeof(char), 1 , fileObject );
27259599516SKenneth E. Jansen if ( 362436 != integer_value ) Wrong_Endian = true;
27359599516SKenneth E. Jansen } else{
27459599516SKenneth E. Jansen fscanf(fileObject, "%d\n", &integer_value );
27559599516SKenneth E. Jansen }
27659599516SKenneth E. Jansen } else {
27759599516SKenneth E. Jansen /* some other header, so just skip over */
27859599516SKenneth E. Jansen token = strtok( NULL, " ,;<>" );
27997a07b0aSCameron Smith assert(token);
28059599516SKenneth E. Jansen skip_size = atoi( token );
28159599516SKenneth E. Jansen if ( binary_format)
28259599516SKenneth E. Jansen fseek( fileObject, skip_size, SEEK_CUR );
28359599516SKenneth E. Jansen else
28459599516SKenneth E. Jansen for( int gama=0; gama < skip_size; gama++ )
28559599516SKenneth E. Jansen fgets( Line, 1024, fileObject );
28659599516SKenneth E. Jansen }
28759599516SKenneth E. Jansen free (text_header);
28859599516SKenneth E. Jansen } // end of if before while loop
28959599516SKenneth E. Jansen
29059599516SKenneth E. Jansen if ( !FOUND )
29159599516SKenneth E. Jansen if( !fgets( Line, 1024, fileObject ) && feof( fileObject ) ) {
29259599516SKenneth E. Jansen rewind( fileObject );
29359599516SKenneth E. Jansen clearerr( fileObject );
29459599516SKenneth E. Jansen rewind_count++;
29559599516SKenneth E. Jansen fgets( Line, 1024, fileObject );
29659599516SKenneth E. Jansen }
29759599516SKenneth E. Jansen }
29859599516SKenneth E. Jansen
29959599516SKenneth E. Jansen if ( !FOUND ) {
30059599516SKenneth E. Jansen //fprintf(stderr, "Error: Could not find: %s\n", phrase);
30159599516SKenneth E. Jansen if(irank == 0) printf("WARNING: Could not find: %s\n", phrase);
30259599516SKenneth E. Jansen return 1;
30359599516SKenneth E. Jansen }
30459599516SKenneth E. Jansen return 0;
30559599516SKenneth E. Jansen }
30659599516SKenneth E. Jansen
30759599516SKenneth E. Jansen } // end unnamed namespace
30859599516SKenneth E. Jansen
30959599516SKenneth E. Jansen
31059599516SKenneth E. Jansen // begin of publicly visible functions
31159599516SKenneth E. Jansen
31259599516SKenneth E. Jansen /**
31359599516SKenneth E. Jansen * This function takes a long long pointer and assign (start) phiotmrc value to it
31459599516SKenneth E. Jansen */
startTimer(double * start)31559599516SKenneth E. Jansen void startTimer(double* start) {
31659599516SKenneth E. Jansen if( !PRINT_PERF ) return;
31759599516SKenneth E. Jansen MPI_Barrier(MPI_COMM_WORLD);
31859599516SKenneth E. Jansen *start = phiotmrc();
31959599516SKenneth E. Jansen }
32059599516SKenneth E. Jansen
32159599516SKenneth E. Jansen /**
32259599516SKenneth E. Jansen * This function takes a long long pointer and assign (end) phiotmrc value to it
32359599516SKenneth E. Jansen */
endTimer(double * end)32459599516SKenneth E. Jansen void endTimer(double* end) {
32559599516SKenneth E. Jansen if( !PRINT_PERF ) return;
32659599516SKenneth E. Jansen *end = phiotmrc();
32759599516SKenneth E. Jansen MPI_Barrier(MPI_COMM_WORLD);
32859599516SKenneth E. Jansen }
32959599516SKenneth E. Jansen
33059599516SKenneth E. Jansen /**
33159599516SKenneth E. Jansen * choose to print some performance results (or not) according to
33259599516SKenneth E. Jansen * the PRINT_PERF macro
33359599516SKenneth E. Jansen */
printPerf(const char * func_name,double start,double end,unsigned long datasize,int printdatainfo,const char * extra_msg)33459599516SKenneth E. Jansen void printPerf(
33559599516SKenneth E. Jansen const char* func_name,
33659599516SKenneth E. Jansen double start,
33759599516SKenneth E. Jansen double end,
3382dd307a1SCameron Smith unsigned long datasize,
33959599516SKenneth E. Jansen int printdatainfo,
34059599516SKenneth E. Jansen const char* extra_msg) {
34159599516SKenneth E. Jansen if( !PRINT_PERF ) return;
3422dd307a1SCameron Smith unsigned long data_size = datasize;
34359599516SKenneth E. Jansen double time = end - start;
3442dd307a1SCameron Smith unsigned long isizemin,isizemax,isizetot;
34559599516SKenneth E. Jansen double sizemin,sizemax,sizeavg,sizetot,rate;
34659599516SKenneth E. Jansen double tmin, tmax, tavg, ttot;
34759599516SKenneth E. Jansen
34859599516SKenneth E. Jansen MPI_Allreduce(&time, &tmin,1, MPI_DOUBLE, MPI_MIN, MPI_COMM_WORLD);
34959599516SKenneth E. Jansen MPI_Allreduce(&time, &tmax,1, MPI_DOUBLE, MPI_MAX, MPI_COMM_WORLD);
35059599516SKenneth E. Jansen MPI_Allreduce(&time, &ttot,1, MPI_DOUBLE, MPI_SUM, MPI_COMM_WORLD);
35159599516SKenneth E. Jansen tavg = ttot/mysize;
35259599516SKenneth E. Jansen
35359599516SKenneth E. Jansen if(irank == 0) {
35459599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) printf("** 1PFPP ");
35559599516SKenneth E. Jansen else printf("** syncIO ");
35659599516SKenneth E. Jansen printf("%s(): Tmax = %f sec, Tmin = %f sec, Tavg = %f sec", func_name, tmax, tmin, tavg);
35759599516SKenneth E. Jansen }
35859599516SKenneth E. Jansen
35959599516SKenneth E. Jansen if(printdatainfo == 1) { // if printdatainfo ==1, compute I/O rate and block size
36059599516SKenneth E. Jansen MPI_Allreduce(&data_size,&isizemin,1,MPI_LONG_LONG_INT,MPI_MIN,MPI_COMM_WORLD);
36159599516SKenneth E. Jansen MPI_Allreduce(&data_size,&isizemax,1,MPI_LONG_LONG_INT,MPI_MAX,MPI_COMM_WORLD);
36259599516SKenneth E. Jansen MPI_Allreduce(&data_size,&isizetot,1,MPI_LONG_LONG_INT,MPI_SUM,MPI_COMM_WORLD);
36359599516SKenneth E. Jansen
36459599516SKenneth E. Jansen sizemin=(double)(isizemin*inv1024sq);
36559599516SKenneth E. Jansen sizemax=(double)(isizemax*inv1024sq);
36659599516SKenneth E. Jansen sizetot=(double)(isizetot*inv1024sq);
36759599516SKenneth E. Jansen sizeavg=(double)(1.0*sizetot/mysize);
36859599516SKenneth E. Jansen rate=(double)(1.0*sizetot/tmax);
36959599516SKenneth E. Jansen
37059599516SKenneth E. Jansen if( irank == 0) {
371103be424SCameron Smith printf(", Rate = %f MB/s [%s] \n \t\t\t"
372103be424SCameron Smith " block size: Min= %f MB; Max= %f MB; Avg= %f MB; Tot= %f MB\n",
373103be424SCameron Smith rate, extra_msg, sizemin,sizemax,sizeavg,sizetot);
37459599516SKenneth E. Jansen }
37559599516SKenneth E. Jansen }
37659599516SKenneth E. Jansen else {
37759599516SKenneth E. Jansen if(irank == 0) {
37859599516SKenneth E. Jansen printf(" \n");
37959599516SKenneth E. Jansen //printf(" (%s) \n", extra_msg);
38059599516SKenneth E. Jansen }
38159599516SKenneth E. Jansen }
38259599516SKenneth E. Jansen }
38359599516SKenneth E. Jansen
38459599516SKenneth E. Jansen /**
38559599516SKenneth E. Jansen * This function is normally called at the beginning of a read operation, before
38659599516SKenneth E. Jansen * init function.
38759599516SKenneth E. Jansen * This function (uses rank 0) reads out nfields, nppf, master header size,
38859599516SKenneth E. Jansen * endianess and allocates for masterHeader string.
38959599516SKenneth E. Jansen * These values are essential for following read operations. Rank 0 will bcast
39059599516SKenneth E. Jansen * these values to other ranks in the commm world
39159599516SKenneth E. Jansen *
39259599516SKenneth E. Jansen * If the file set is of old POSIX format, it would throw error and exit
39359599516SKenneth E. Jansen */
queryphmpiio(const char filename[],int * nfields,int * nppf)39459599516SKenneth E. Jansen void queryphmpiio(const char filename[],int *nfields, int *nppf)
39559599516SKenneth E. Jansen {
39659599516SKenneth E. Jansen MPI_Comm_rank(MPI_COMM_WORLD, &irank);
39759599516SKenneth E. Jansen MPI_Comm_size(MPI_COMM_WORLD, &mysize);
39859599516SKenneth E. Jansen
39959599516SKenneth E. Jansen if(irank == 0) {
40059599516SKenneth E. Jansen FILE * fileHandle;
40159599516SKenneth E. Jansen char* fname = StringStripper( filename );
40259599516SKenneth E. Jansen
4038f9016f6SCameron Smith PHASTAIO_OPENTIME(fileHandle = fopen (fname,"rb");)
40459599516SKenneth E. Jansen if (fileHandle == NULL ) {
40559599516SKenneth E. Jansen printf("\nError: File %s doesn't exist! Please check!\n",fname);
40659599516SKenneth E. Jansen }
40759599516SKenneth E. Jansen else {
40859599516SKenneth E. Jansen SerialFile =(serial_file *)calloc( 1, sizeof( serial_file) );
40959599516SKenneth E. Jansen int meta_size_limit = VERSION_INFO_HEADER_SIZE;
41059599516SKenneth E. Jansen SerialFile->masterHeader = (char *)malloc( meta_size_limit );
41159599516SKenneth E. Jansen fread(SerialFile->masterHeader, 1, meta_size_limit, fileHandle);
41259599516SKenneth E. Jansen
41359599516SKenneth E. Jansen char read_out_tag[MAX_FIELDS_NAME_LENGTH];
41459599516SKenneth E. Jansen char version[MAX_FIELDS_NAME_LENGTH/4];
41559599516SKenneth E. Jansen int mhsize;
41659599516SKenneth E. Jansen char * token;
41759599516SKenneth E. Jansen int magic_number;
41859599516SKenneth E. Jansen
41959599516SKenneth E. Jansen memcpy( read_out_tag,
42059599516SKenneth E. Jansen SerialFile->masterHeader,
42159599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH-1 );
42259599516SKenneth E. Jansen
42359599516SKenneth E. Jansen if ( cscompare ("MPI_IO_Tag",read_out_tag) ) {
42459599516SKenneth E. Jansen // Test endianess ...
42559599516SKenneth E. Jansen memcpy (&magic_number,
42659599516SKenneth E. Jansen SerialFile->masterHeader + sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
42759599516SKenneth E. Jansen sizeof(int) ); // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
42859599516SKenneth E. Jansen
42959599516SKenneth E. Jansen if ( magic_number != ENDIAN_TEST_NUMBER ) {
43059599516SKenneth E. Jansen printf("Endian is different!\n");
43159599516SKenneth E. Jansen // Will do swap later
43259599516SKenneth E. Jansen }
43359599516SKenneth E. Jansen
43459599516SKenneth E. Jansen // test version, old version, default masterheader size is 4M
43559599516SKenneth E. Jansen // newer version masterheader size is read from first line
43659599516SKenneth E. Jansen memcpy(version,
43759599516SKenneth E. Jansen SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/2,
43859599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH/4 - 1); //TODO: why -1?
43959599516SKenneth E. Jansen
44059599516SKenneth E. Jansen if( cscompare ("version",version) ) {
44159599516SKenneth E. Jansen // if there is "version" tag in the file, then it is newer format
44259599516SKenneth E. Jansen // read master header size from here, otherwise use default
44359599516SKenneth E. Jansen // Note: if version is "1", we know mhsize is at 3/4 place...
44459599516SKenneth E. Jansen
44559599516SKenneth E. Jansen token = strtok(version, ":");
44659599516SKenneth E. Jansen token = strtok(NULL, " ,;<>" );
44759599516SKenneth E. Jansen int iversion = atoi(token);
44859599516SKenneth E. Jansen
44959599516SKenneth E. Jansen if( iversion == 1) {
45059599516SKenneth E. Jansen memcpy( &mhsize,
45159599516SKenneth E. Jansen SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH/4*3 + sizeof("mhsize : ")-1,
45259599516SKenneth E. Jansen sizeof(int));
45359599516SKenneth E. Jansen if ( magic_number != ENDIAN_TEST_NUMBER )
45459599516SKenneth E. Jansen SwapArrayByteOrder(&mhsize, sizeof(int), 1);
45559599516SKenneth E. Jansen
45659599516SKenneth E. Jansen if( mhsize > DefaultMHSize ) {
45759599516SKenneth E. Jansen //if actual headersize is larger than default, let's re-read
45859599516SKenneth E. Jansen free(SerialFile->masterHeader);
45959599516SKenneth E. Jansen SerialFile->masterHeader = (char *)malloc(mhsize);
46059599516SKenneth E. Jansen fseek(fileHandle, 0, SEEK_SET); // reset the file stream position
46159599516SKenneth E. Jansen fread(SerialFile->masterHeader,1,mhsize,fileHandle);
46259599516SKenneth E. Jansen }
46359599516SKenneth E. Jansen }
46459599516SKenneth E. Jansen //TODO: check if this is a valid int??
46559599516SKenneth E. Jansen MasterHeaderSize = mhsize;
46659599516SKenneth E. Jansen }
46759599516SKenneth E. Jansen else { // else it's version 0's format w/o version tag, implicating MHSize=4M
46859599516SKenneth E. Jansen MasterHeaderSize = DefaultMHSize;
46959599516SKenneth E. Jansen }
47059599516SKenneth E. Jansen
47159599516SKenneth E. Jansen memcpy( read_out_tag,
47259599516SKenneth E. Jansen SerialFile->masterHeader+MAX_FIELDS_NAME_LENGTH+1,
47359599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH ); //TODO: why +1
47459599516SKenneth E. Jansen
47559599516SKenneth E. Jansen // Read in # fields ...
47659599516SKenneth E. Jansen token = strtok( read_out_tag, ":" );
47759599516SKenneth E. Jansen token = strtok( NULL," ,;<>" );
47859599516SKenneth E. Jansen *nfields = atoi( token );
47959599516SKenneth E. Jansen if ( *nfields > MAX_FIELDS_NUMBER) {
48059599516SKenneth E. Jansen printf("Error queryphmpiio: nfields is larger than MAX_FIELDS_NUMBER!\n");
48159599516SKenneth E. Jansen }
48259599516SKenneth E. Jansen SerialFile->nfields=*nfields; //TODO: sanity check of this int?
48359599516SKenneth E. Jansen
48459599516SKenneth E. Jansen memcpy( read_out_tag,
48559599516SKenneth E. Jansen SerialFile->masterHeader + MAX_FIELDS_NAME_LENGTH * 2
48659599516SKenneth E. Jansen + *nfields * MAX_FIELDS_NAME_LENGTH,
48759599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH);
48859599516SKenneth E. Jansen
48959599516SKenneth E. Jansen token = strtok( read_out_tag, ":" );
49059599516SKenneth E. Jansen token = strtok( NULL," ,;<>" );
49159599516SKenneth E. Jansen *nppf = atoi( token );
49259599516SKenneth E. Jansen SerialFile->nppf=*nppf; //TODO: sanity check of int
49359599516SKenneth E. Jansen } // end of if("MPI_IO_TAG")
49459599516SKenneth E. Jansen else {
49559599516SKenneth 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);
49659599516SKenneth E. Jansen exit(1);
49759599516SKenneth E. Jansen }
4988f9016f6SCameron Smith PHASTAIO_CLOSETIME(fclose(fileHandle);)
49959599516SKenneth E. Jansen free(SerialFile->masterHeader);
50059599516SKenneth E. Jansen free(SerialFile);
50159599516SKenneth E. Jansen } //end of else
50259599516SKenneth E. Jansen free(fname);
50359599516SKenneth E. Jansen }
50459599516SKenneth E. Jansen
50559599516SKenneth E. Jansen // Bcast value to every one
50659599516SKenneth E. Jansen MPI_Bcast( nfields, 1, MPI_INT, 0, MPI_COMM_WORLD );
50759599516SKenneth E. Jansen MPI_Bcast( nppf, 1, MPI_INT, 0, MPI_COMM_WORLD );
50859599516SKenneth E. Jansen MPI_Bcast( &MasterHeaderSize, 1, MPI_INT, 0, MPI_COMM_WORLD );
509a746c838SCameron Smith phprintf("Info queryphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
51059599516SKenneth E. Jansen }
51159599516SKenneth E. Jansen
51259599516SKenneth E. Jansen /**
51359599516SKenneth E. Jansen * This function computes the right master header size (round to size of 2^n).
51459599516SKenneth E. Jansen * This is only needed for file format version 1 in "write" mode.
51559599516SKenneth E. Jansen */
computeMHSize(int nfields,int nppf,int version)51659599516SKenneth E. Jansen int computeMHSize(int nfields, int nppf, int version) {
51710d56689SCameron Smith int mhsize=0;
51859599516SKenneth E. Jansen if(version == 1) {
51959599516SKenneth 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
52059599516SKenneth E. Jansen int meta_info_size = VERSION_INFO_HEADER_SIZE;
5212dd307a1SCameron Smith int actual_size = nfields * nppf * sizeof(unsigned long) + meta_info_size;
52259599516SKenneth E. Jansen //printf("actual_size = %d, offset table size = %d\n", actual_size, nfields * nppf * sizeof(long long));
52359599516SKenneth E. Jansen if (actual_size > DefaultMHSize) {
52459599516SKenneth E. Jansen mhsize = (int) ceil( (double) actual_size/DefaultMHSize); // it's rounded to ceiling of this value
52559599516SKenneth E. Jansen mhsize *= DefaultMHSize;
52659599516SKenneth E. Jansen }
52759599516SKenneth E. Jansen else {
52859599516SKenneth E. Jansen mhsize = DefaultMHSize;
52959599516SKenneth E. Jansen }
530b744cbcaSCameron Smith } else {
531b744cbcaSCameron Smith int rank = 0;
532b744cbcaSCameron Smith MPI_Comm_rank(MPI_COMM_WORLD, &rank);
533b744cbcaSCameron Smith if(!rank) {
534b744cbcaSCameron Smith fprintf(stderr,
535b744cbcaSCameron Smith "ERROR invalid version passed to %s... exiting\n", __func__);
536b744cbcaSCameron Smith exit(EXIT_FAILURE);
537b744cbcaSCameron Smith }
53859599516SKenneth E. Jansen }
53959599516SKenneth E. Jansen return mhsize;
54059599516SKenneth E. Jansen }
54159599516SKenneth E. Jansen
54259599516SKenneth E. Jansen /**
54359599516SKenneth E. Jansen * Computes correct color of a rank according to number of files.
54459599516SKenneth E. Jansen */
computeColor(int myrank,int numprocs,int nfiles)54559599516SKenneth E. Jansen extern "C" int computeColor( int myrank, int numprocs, int nfiles) {
54659599516SKenneth E. Jansen int color =
54759599516SKenneth E. Jansen (int)(myrank / (numprocs / nfiles));
54859599516SKenneth E. Jansen return color;
54959599516SKenneth E. Jansen }
55059599516SKenneth E. Jansen
55159599516SKenneth E. Jansen
55259599516SKenneth E. Jansen /**
55359599516SKenneth E. Jansen * Check the file descriptor.
55459599516SKenneth E. Jansen */
checkFileDescriptor(const char fctname[],int * fileDescriptor)555103be424SCameron Smith void checkFileDescriptor(const char fctname[], int* fileDescriptor ) {
55659599516SKenneth E. Jansen if ( *fileDescriptor < 0 ) {
55759599516SKenneth E. Jansen printf("Error: File descriptor = %d in %s\n",*fileDescriptor,fctname);
55859599516SKenneth E. Jansen exit(1);
55959599516SKenneth E. Jansen }
56059599516SKenneth E. Jansen }
56159599516SKenneth E. Jansen
56259599516SKenneth E. Jansen /**
56359599516SKenneth E. Jansen * Initialize the file struct members and allocate space for file struct
56459599516SKenneth E. Jansen * buffers.
56559599516SKenneth E. Jansen *
56659599516SKenneth E. Jansen * Note: this function is only called when we are using new format. Old POSIX
56759599516SKenneth E. Jansen * format should skip this routine and call openfile() directly instead.
56859599516SKenneth E. Jansen */
initphmpiio(int * nfields,int * nppf,int * nfiles,int * filehandle,const char mode[])56959599516SKenneth E. Jansen int initphmpiio( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[])
57059599516SKenneth E. Jansen {
57159599516SKenneth E. Jansen // we init irank again in case query not called (e.g. syncIO write case)
57259599516SKenneth E. Jansen MPI_Comm_rank(MPI_COMM_WORLD, &irank);
57359599516SKenneth E. Jansen MPI_Comm_size(MPI_COMM_WORLD, &mysize);
57459599516SKenneth E. Jansen
575a746c838SCameron Smith phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d, nfields %d, nppf %d, nfiles %d\n", irank, MasterHeaderSize, *nfields, *nppf, *nfiles);
57659599516SKenneth E. Jansen
57759599516SKenneth E. Jansen double timer_start, timer_end;
57859599516SKenneth E. Jansen startTimer(&timer_start);
57959599516SKenneth E. Jansen
58059599516SKenneth E. Jansen char* imode = StringStripper( mode );
58159599516SKenneth E. Jansen
58259599516SKenneth E. Jansen // Note: if it's read, we presume query was called prior to init and
58359599516SKenneth E. Jansen // MasterHeaderSize is already set to correct value from parsing header
58459599516SKenneth E. Jansen // otherwise it's write then it needs some computation to be set
58559599516SKenneth E. Jansen if ( cscompare( "read", imode ) ) {
58659599516SKenneth E. Jansen // do nothing
58759599516SKenneth E. Jansen }
58859599516SKenneth E. Jansen else if( cscompare( "write", imode ) ) {
58959599516SKenneth E. Jansen MasterHeaderSize = computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
59059599516SKenneth E. Jansen }
59159599516SKenneth E. Jansen else {
59259599516SKenneth E. Jansen printf("Error initphmpiio: can't recognize the mode %s", imode);
59359599516SKenneth E. Jansen exit(1);
59459599516SKenneth E. Jansen }
59559599516SKenneth E. Jansen free ( imode );
59659599516SKenneth E. Jansen
597a746c838SCameron Smith phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
59859599516SKenneth E. Jansen
59959599516SKenneth E. Jansen int i, j;
60059599516SKenneth E. Jansen
60159599516SKenneth E. Jansen if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
60259599516SKenneth E. Jansen printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
60359599516SKenneth E. Jansen endTimer(&timer_end);
60459599516SKenneth E. Jansen printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
60559599516SKenneth E. Jansen return MAX_PHASTA_FILES_EXCEEDED;
60659599516SKenneth E. Jansen }
60759599516SKenneth E. Jansen // else if( PhastaIONextActiveIndex == 0 ) //Hang in debug mode on Intrepid
60859599516SKenneth E. Jansen // {
60959599516SKenneth E. Jansen // for( i = 0; i < MAX_PHASTA_FILES; i++ );
61059599516SKenneth E. Jansen // {
61159599516SKenneth E. Jansen // PhastaIOActiveFiles[i] = NULL;
61259599516SKenneth E. Jansen // }
61359599516SKenneth E. Jansen // }
61459599516SKenneth E. Jansen
61559599516SKenneth E. Jansen
61659599516SKenneth E. Jansen PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1, sizeof( phastaio_file_t) );
61759599516SKenneth E. Jansen
61859599516SKenneth E. Jansen i = PhastaIONextActiveIndex;
61959599516SKenneth E. Jansen PhastaIONextActiveIndex++;
62059599516SKenneth E. Jansen
62159599516SKenneth E. Jansen //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
62259599516SKenneth E. Jansen
62359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize; // what does this mean??? TODO
62459599516SKenneth E. Jansen
62559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->Wrong_Endian = false;
62659599516SKenneth E. Jansen
62759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFields = *nfields;
62859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nPPF = *nppf;
62959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFiles = *nfiles;
63059599516SKenneth E. Jansen MPI_Comm_rank(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->myrank));
63159599516SKenneth E. Jansen MPI_Comm_size(MPI_COMM_WORLD, &(PhastaIOActiveFiles[i]->numprocs));
63259599516SKenneth E. Jansen
63359599516SKenneth E. Jansen
63459599516SKenneth E. Jansen if( *nfiles > 1 ) { // split the ranks according to each mpiio file
63559599516SKenneth E. Jansen
63659599516SKenneth E. Jansen if ( s_assign_local_comm == 0) { // call mpi_comm_split for the first (and only) time
63759599516SKenneth E. Jansen
63859599516SKenneth E. Jansen if (PhastaIOActiveFiles[i]->myrank == 0) printf("Building subcommunicator\n");
63959599516SKenneth E. Jansen
64059599516SKenneth E. Jansen int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
64159599516SKenneth E. Jansen MPI_Comm_split(MPI_COMM_WORLD,
64259599516SKenneth E. Jansen color,
64359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->myrank,
64459599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_comm));
64559599516SKenneth E. Jansen MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
64659599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_numprocs));
64759599516SKenneth E. Jansen MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
64859599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_myrank));
64959599516SKenneth E. Jansen
65059599516SKenneth E. Jansen // back up now these variables so that we do not need to call comm_split again
65159599516SKenneth E. Jansen s_local_comm = PhastaIOActiveFiles[i]->local_comm;
65259599516SKenneth E. Jansen s_local_size = PhastaIOActiveFiles[i]->local_numprocs;
65359599516SKenneth E. Jansen s_local_rank = PhastaIOActiveFiles[i]->local_myrank;
65459599516SKenneth E. Jansen s_assign_local_comm = 1;
65559599516SKenneth E. Jansen }
65659599516SKenneth E. Jansen else { // recycle the subcommunicator
65759599516SKenneth E. Jansen if (PhastaIOActiveFiles[i]->myrank == 0) printf("Recycling subcommunicator\n");
65859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm = s_local_comm;
65959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_numprocs = s_local_size;
66059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_myrank = s_local_rank;
66159599516SKenneth E. Jansen }
66259599516SKenneth E. Jansen }
66359599516SKenneth E. Jansen else { // *nfiles == 1 here - no need to call mpi_comm_split here
66459599516SKenneth E. Jansen
66559599516SKenneth E. Jansen if (PhastaIOActiveFiles[i]->myrank == 0) printf("Bypassing subcommunicator\n");
66659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm = MPI_COMM_WORLD;
66759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_numprocs = PhastaIOActiveFiles[i]->numprocs;
66859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_myrank = PhastaIOActiveFiles[i]->myrank;
66959599516SKenneth E. Jansen
67059599516SKenneth E. Jansen }
67159599516SKenneth E. Jansen
67259599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp =
67359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
67459599516SKenneth E. Jansen
67559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
67659599516SKenneth E. Jansen (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
67759599516SKenneth E. Jansen (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
67859599516SKenneth E. Jansen
67959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset_table =
6802dd307a1SCameron Smith ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
68159599516SKenneth E. Jansen
68259599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table =
6832dd307a1SCameron Smith ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
68459599516SKenneth E. Jansen
68559599516SKenneth E. Jansen for (j=0; j<*nfields; j++)
68659599516SKenneth E. Jansen {
68759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset_table[j] =
6882dd307a1SCameron Smith ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
68959599516SKenneth E. Jansen
69059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table[j] =
6912dd307a1SCameron Smith ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
69259599516SKenneth E. Jansen }
69359599516SKenneth E. Jansen *filehandle = i;
69459599516SKenneth E. Jansen
69559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
69659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
69759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
69859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
69959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
70059599516SKenneth E. Jansen
70159599516SKenneth E. Jansen // Time monitoring
70259599516SKenneth E. Jansen endTimer(&timer_end);
70359599516SKenneth E. Jansen printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
70459599516SKenneth E. Jansen
70559599516SKenneth E. Jansen phprintf_0("Info initphmpiio: quiting function");
70659599516SKenneth E. Jansen
70759599516SKenneth E. Jansen return i;
70859599516SKenneth E. Jansen }
70959599516SKenneth E. Jansen
71059599516SKenneth E. Jansen /**
71159599516SKenneth E. Jansen * Destruct the file struct and free buffers allocated in init function.
71259599516SKenneth E. Jansen */
finalizephmpiio(int * fileDescriptor)71359599516SKenneth E. Jansen void finalizephmpiio( int *fileDescriptor )
71459599516SKenneth E. Jansen {
71559599516SKenneth E. Jansen double timer_start, timer_end;
71659599516SKenneth E. Jansen startTimer(&timer_start);
71759599516SKenneth E. Jansen
71859599516SKenneth E. Jansen int i, j;
71959599516SKenneth E. Jansen i = *fileDescriptor;
72059599516SKenneth E. Jansen //PhastaIONextActiveIndex--;
72159599516SKenneth E. Jansen
72259599516SKenneth E. Jansen /* //free the offset table for this phasta file */
72359599516SKenneth E. Jansen //for(j=0; j<MAX_FIELDS_NUMBER; j++) //Danger: undefined behavior for my_*_table.[j] not allocated or not initialized to NULL
72459599516SKenneth E. Jansen for(j=0; j<PhastaIOActiveFiles[i]->nFields; j++)
72559599516SKenneth E. Jansen {
72659599516SKenneth E. Jansen free( PhastaIOActiveFiles[i]->my_offset_table[j]);
72759599516SKenneth E. Jansen free( PhastaIOActiveFiles[i]->my_read_table[j]);
72859599516SKenneth E. Jansen }
72959599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->my_offset_table );
73059599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->my_read_table );
73159599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->master_header );
73259599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->double_chunk );
73359599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->int_chunk );
73459599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->read_double_chunk );
73559599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->read_int_chunk );
73659599516SKenneth E. Jansen
737946a6cddSCameron Smith if( PhastaIOActiveFiles[i]->nFiles > 1 && s_assign_local_comm ) { // the comm was split
738946a6cddSCameron Smith if (PhastaIOActiveFiles[i]->myrank == 0) printf("Freeing subcommunicator\n");
739946a6cddSCameron Smith s_assign_local_comm = 0;
740946a6cddSCameron Smith MPI_Comm_free(&(PhastaIOActiveFiles[i]->local_comm));
741946a6cddSCameron Smith }
742946a6cddSCameron Smith
74359599516SKenneth E. Jansen free( PhastaIOActiveFiles[i]);
74459599516SKenneth E. Jansen
74559599516SKenneth E. Jansen endTimer(&timer_end);
74659599516SKenneth E. Jansen printPerf("finalizempiio", timer_start, timer_end, 0, 0, "");
74759599516SKenneth E. Jansen
74859599516SKenneth E. Jansen PhastaIONextActiveIndex--;
74959599516SKenneth E. Jansen }
75059599516SKenneth E. Jansen
75159599516SKenneth E. Jansen
75259599516SKenneth E. Jansen /**
75359599516SKenneth E. Jansen * Special init for M2N in order to create a subcommunicator for the reduced solution (requires PRINT_PERF to be false for now)
75459599516SKenneth E. Jansen * Initialize the file struct members and allocate space for file struct buffers.
75559599516SKenneth E. Jansen *
75659599516SKenneth E. Jansen * Note: this function is only called when we are using new format. Old POSIX
75759599516SKenneth E. Jansen * format should skip this routine and call openfile() directly instead.
75859599516SKenneth E. Jansen */
initphmpiiosub(int * nfields,int * nppf,int * nfiles,int * filehandle,const char mode[],MPI_Comm my_local_comm)75959599516SKenneth E. Jansen int initphmpiiosub( int *nfields, int *nppf, int *nfiles, int *filehandle, const char mode[],MPI_Comm my_local_comm)
76059599516SKenneth E. Jansen {
76159599516SKenneth E. Jansen // we init irank again in case query not called (e.g. syncIO write case)
76259599516SKenneth E. Jansen
76359599516SKenneth E. Jansen MPI_Comm_rank(my_local_comm, &irank);
76459599516SKenneth E. Jansen MPI_Comm_size(my_local_comm, &mysize);
76559599516SKenneth E. Jansen
766a746c838SCameron Smith phprintf("Info initphmpiio: entering function, myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
76759599516SKenneth E. Jansen
76859599516SKenneth E. Jansen double timer_start, timer_end;
76959599516SKenneth E. Jansen startTimer(&timer_start);
77059599516SKenneth E. Jansen
77159599516SKenneth E. Jansen char* imode = StringStripper( mode );
77259599516SKenneth E. Jansen
77359599516SKenneth E. Jansen // Note: if it's read, we presume query was called prior to init and
77459599516SKenneth E. Jansen // MasterHeaderSize is already set to correct value from parsing header
77559599516SKenneth E. Jansen // otherwise it's write then it needs some computation to be set
77659599516SKenneth E. Jansen if ( cscompare( "read", imode ) ) {
77759599516SKenneth E. Jansen // do nothing
77859599516SKenneth E. Jansen }
77959599516SKenneth E. Jansen else if( cscompare( "write", imode ) ) {
78059599516SKenneth E. Jansen MasterHeaderSize = computeMHSize(*nfields, *nppf, LATEST_WRITE_VERSION);
78159599516SKenneth E. Jansen }
78259599516SKenneth E. Jansen else {
78359599516SKenneth E. Jansen printf("Error initphmpiio: can't recognize the mode %s", imode);
78459599516SKenneth E. Jansen exit(1);
78559599516SKenneth E. Jansen }
78659599516SKenneth E. Jansen free ( imode );
78759599516SKenneth E. Jansen
788a746c838SCameron Smith phprintf("Info initphmpiio: myrank = %d, MasterHeaderSize = %d\n", irank, MasterHeaderSize);
78959599516SKenneth E. Jansen
79059599516SKenneth E. Jansen int i, j;
79159599516SKenneth E. Jansen
79259599516SKenneth E. Jansen if( PhastaIONextActiveIndex == MAX_PHASTA_FILES ) {
79359599516SKenneth E. Jansen printf("Error initphmpiio: PhastaIONextActiveIndex = MAX_PHASTA_FILES");
79459599516SKenneth E. Jansen endTimer(&timer_end);
79559599516SKenneth E. Jansen printPerf("initphmpiio", timer_start, timer_end, 0, 0, "");
79659599516SKenneth E. Jansen return MAX_PHASTA_FILES_EXCEEDED;
79759599516SKenneth E. Jansen }
79859599516SKenneth E. Jansen // else if( PhastaIONextActiveIndex == 0 ) //Hang in debug mode on Intrepid
79959599516SKenneth E. Jansen // {
80059599516SKenneth E. Jansen // for( i = 0; i < MAX_PHASTA_FILES; i++ );
80159599516SKenneth E. Jansen // {
80259599516SKenneth E. Jansen // PhastaIOActiveFiles[i] = NULL;
80359599516SKenneth E. Jansen // }
80459599516SKenneth E. Jansen // }
80559599516SKenneth E. Jansen
80659599516SKenneth E. Jansen
80759599516SKenneth E. Jansen PhastaIOActiveFiles[PhastaIONextActiveIndex] = (phastaio_file_t *)calloc( 1, sizeof( phastaio_file_t) );
80859599516SKenneth E. Jansen
80959599516SKenneth E. Jansen i = PhastaIONextActiveIndex;
81059599516SKenneth E. Jansen PhastaIONextActiveIndex++;
81159599516SKenneth E. Jansen
81259599516SKenneth E. Jansen //PhastaIOActiveFiles[i]->next_start_address = 2*TWO_MEGABYTE;
81359599516SKenneth E. Jansen
81459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->next_start_address = MasterHeaderSize; // what does this mean??? TODO
81559599516SKenneth E. Jansen
81659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->Wrong_Endian = false;
81759599516SKenneth E. Jansen
81859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFields = *nfields;
81959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nPPF = *nppf;
82059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFiles = *nfiles;
82159599516SKenneth E. Jansen MPI_Comm_rank(my_local_comm, &(PhastaIOActiveFiles[i]->myrank));
82259599516SKenneth E. Jansen MPI_Comm_size(my_local_comm, &(PhastaIOActiveFiles[i]->numprocs));
82359599516SKenneth E. Jansen
82459599516SKenneth E. Jansen int color = computeColor(PhastaIOActiveFiles[i]->myrank, PhastaIOActiveFiles[i]->numprocs, PhastaIOActiveFiles[i]->nFiles);
82559599516SKenneth E. Jansen MPI_Comm_split(my_local_comm,
82659599516SKenneth E. Jansen color,
82759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->myrank,
82859599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_comm));
82959599516SKenneth E. Jansen MPI_Comm_size(PhastaIOActiveFiles[i]->local_comm,
83059599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_numprocs));
83159599516SKenneth E. Jansen MPI_Comm_rank(PhastaIOActiveFiles[i]->local_comm,
83259599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->local_myrank));
83359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp =
83459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nPPF/PhastaIOActiveFiles[i]->local_numprocs;
83559599516SKenneth E. Jansen
83659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->start_id = PhastaIOActiveFiles[i]->nPPF *
83759599516SKenneth E. Jansen (int)(PhastaIOActiveFiles[i]->myrank/PhastaIOActiveFiles[i]->local_numprocs) +
83859599516SKenneth E. Jansen (PhastaIOActiveFiles[i]->local_myrank * PhastaIOActiveFiles[i]->nppp);
83959599516SKenneth E. Jansen
84059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset_table =
8412dd307a1SCameron Smith ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
84259599516SKenneth E. Jansen
84359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table =
8442dd307a1SCameron Smith ( unsigned long ** ) calloc( MAX_FIELDS_NUMBER , sizeof( unsigned long *) );
84559599516SKenneth E. Jansen
84659599516SKenneth E. Jansen for (j=0; j<*nfields; j++)
84759599516SKenneth E. Jansen {
84859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset_table[j] =
8492dd307a1SCameron Smith ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
85059599516SKenneth E. Jansen
85159599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table[j] =
8522dd307a1SCameron Smith ( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nppp , sizeof( unsigned long) );
85359599516SKenneth E. Jansen }
85459599516SKenneth E. Jansen *filehandle = i;
85559599516SKenneth E. Jansen
85659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header = (char *)calloc(MasterHeaderSize,sizeof( char ));
85759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->double_chunk = (double *)calloc(1,sizeof( double ));
85859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->int_chunk = (int *)calloc(1,sizeof( int ));
85959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_double_chunk = (double *)calloc(1,sizeof( double ));
86059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_int_chunk = (int *)calloc(1,sizeof( int ));
86159599516SKenneth E. Jansen
86259599516SKenneth E. Jansen // Time monitoring
86359599516SKenneth E. Jansen endTimer(&timer_end);
86459599516SKenneth E. Jansen printPerf("initphmpiiosub", timer_start, timer_end, 0, 0, "");
86559599516SKenneth E. Jansen
86659599516SKenneth E. Jansen phprintf_0("Info initphmpiiosub: quiting function");
86759599516SKenneth E. Jansen
86859599516SKenneth E. Jansen return i;
86959599516SKenneth E. Jansen }
87059599516SKenneth E. Jansen
8714600f8abSCameron Smith namespace {
87259599516SKenneth E. Jansen
8734600f8abSCameron Smith enum {
8744600f8abSCameron Smith DIR_MODE = S_IRWXU | S_IRWXG | S_IROTH | S_IXOTH
8754600f8abSCameron Smith };
8764600f8abSCameron Smith
my_mkdir(std::string name)8775f4c514aSCameron Smith bool my_mkdir(std::string name) {
8785f4c514aSCameron Smith if(name.empty())
8795f4c514aSCameron Smith return true;
8804600f8abSCameron Smith errno = 0;
8815f4c514aSCameron Smith int err = mkdir(name.c_str(), DIR_MODE);
8824600f8abSCameron Smith if ((err == -1) && (errno == EEXIST)) {
8834600f8abSCameron Smith errno = 0;
8844600f8abSCameron Smith err = 0;
8854600f8abSCameron Smith return false;
8864600f8abSCameron Smith }
8874600f8abSCameron Smith assert(!err);
8884600f8abSCameron Smith return true;
8894600f8abSCameron Smith }
8904600f8abSCameron Smith
8914600f8abSCameron Smith enum {
8924600f8abSCameron Smith DIR_FANOUT = 2048
8934600f8abSCameron Smith };
8944600f8abSCameron Smith
getSubDirPrefix()8954600f8abSCameron Smith std::string getSubDirPrefix() {
8964600f8abSCameron Smith if (phio_peers() <= DIR_FANOUT)
8974600f8abSCameron Smith return string("");
8984600f8abSCameron Smith int self = phio_self();
8994600f8abSCameron Smith int subSelf = self % DIR_FANOUT;
9004600f8abSCameron Smith int subGroup = self / DIR_FANOUT;
9014600f8abSCameron Smith std::stringstream ss;
9024600f8abSCameron Smith ss << subGroup << '/';
9034600f8abSCameron Smith return ss.str();
9044600f8abSCameron Smith }
9054600f8abSCameron Smith }
90659599516SKenneth E. Jansen
90759599516SKenneth E. Jansen /** open file for both POSIX and MPI-IO syncIO format.
90859599516SKenneth E. Jansen *
90959599516SKenneth E. Jansen * If it's old POSIX format, simply call posix fopen().
91059599516SKenneth E. Jansen *
91159599516SKenneth E. Jansen * If it's MPI-IO foramt:
91259599516SKenneth E. Jansen * in "read" mode, it builds the header table that points to the offset of
91359599516SKenneth E. Jansen * fields for parts;
91459599516SKenneth E. Jansen * in "write" mode, it opens the file with MPI-IO open routine.
91559599516SKenneth E. Jansen */
openfile(const char filename[],const char mode[],int * fileDescriptor)916103be424SCameron Smith void openfile(const char filename[], const char mode[], int* fileDescriptor )
91759599516SKenneth E. Jansen {
91859599516SKenneth E. Jansen phprintf_0("Info: entering openfile");
91959599516SKenneth E. Jansen
92059599516SKenneth E. Jansen double timer_start, timer_end;
92159599516SKenneth E. Jansen startTimer(&timer_start);
92259599516SKenneth E. Jansen
92359599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 )
92459599516SKenneth E. Jansen {
92559599516SKenneth E. Jansen FILE* file=NULL ;
92659599516SKenneth E. Jansen *fileDescriptor = 0;
92759599516SKenneth E. Jansen char* fname = StringStripper( filename );
92859599516SKenneth E. Jansen char* imode = StringStripper( mode );
92959599516SKenneth E. Jansen
9304600f8abSCameron Smith std::string posixname = getSubDirPrefix();
9314600f8abSCameron Smith if (!phio_self())
9325f4c514aSCameron Smith my_mkdir(posixname);
9334600f8abSCameron Smith phio_barrier();
9344600f8abSCameron Smith posixname += string(fname);
9354600f8abSCameron Smith if ( cscompare( "read", imode ) )
9368f9016f6SCameron Smith PHASTAIO_OPENTIME(file = fopen(posixname.c_str(), "rb" );)
9374600f8abSCameron Smith else if( cscompare( "write", imode ) )
9388f9016f6SCameron Smith PHASTAIO_OPENTIME(file = fopen(posixname.c_str(), "wb" );)
9394600f8abSCameron Smith else if( cscompare( "append", imode ) )
9408f9016f6SCameron Smith PHASTAIO_OPENTIME(file = fopen(posixname.c_str(), "ab" );)
94159599516SKenneth E. Jansen
94259599516SKenneth E. Jansen if ( !file ){
94367701978SCameron Smith fprintf(stderr,"Error openfile: unable to open file %s\n",fname);
94459599516SKenneth E. Jansen } else {
94559599516SKenneth E. Jansen fileArray.push_back( file );
94659599516SKenneth E. Jansen byte_order.push_back( false );
94759599516SKenneth E. Jansen header_type.push_back( sizeof(int) );
94859599516SKenneth E. Jansen *fileDescriptor = fileArray.size();
94959599516SKenneth E. Jansen }
95059599516SKenneth E. Jansen free (fname);
95159599516SKenneth E. Jansen free (imode);
95259599516SKenneth E. Jansen }
95359599516SKenneth E. Jansen else // else it would be parallel I/O, opposed to posix io
95459599516SKenneth E. Jansen {
95559599516SKenneth E. Jansen char* fname = StringStripper( filename );
95659599516SKenneth E. Jansen char* imode = StringStripper( mode );
95759599516SKenneth E. Jansen int rc;
95859599516SKenneth E. Jansen int i = *fileDescriptor;
95959599516SKenneth E. Jansen checkFileDescriptor("openfile",&i);
96059599516SKenneth E. Jansen char* token;
96159599516SKenneth E. Jansen
96259599516SKenneth E. Jansen if ( cscompare( "read", imode ) )
96359599516SKenneth E. Jansen {
96459599516SKenneth E. Jansen // if (PhastaIOActiveFiles[i]->myrank == 0)
96559599516SKenneth E. Jansen // printf("\n **********\nRead open ... ... regular version\n");
96659599516SKenneth E. Jansen
9678f9016f6SCameron Smith PHASTAIO_OPENTIME(
96859599516SKenneth E. Jansen rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
96959599516SKenneth E. Jansen fname,
97059599516SKenneth E. Jansen MPI_MODE_RDONLY,
97159599516SKenneth E. Jansen MPI_INFO_NULL,
97259599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->file_handle) );
9738f9016f6SCameron Smith )
97459599516SKenneth E. Jansen
97559599516SKenneth E. Jansen if(rc)
97659599516SKenneth E. Jansen {
97759599516SKenneth E. Jansen *fileDescriptor = UNABLE_TO_OPEN_FILE;
9784d60bba2SCameron Smith int error_string_length;
9794d60bba2SCameron Smith char error_string[4096];
9804d60bba2SCameron Smith MPI_Error_string(rc, error_string, &error_string_length);
9814d60bba2SCameron Smith fprintf(stderr, "Error openfile: Unable to open file %s! MPI reports \"%s\"\n",fname,error_string);
98259599516SKenneth E. Jansen endTimer(&timer_end);
98359599516SKenneth E. Jansen printPerf("openfile", timer_start, timer_end, 0, 0, "");
98459599516SKenneth E. Jansen return;
98559599516SKenneth E. Jansen }
98659599516SKenneth E. Jansen
98759599516SKenneth E. Jansen MPI_Status read_tag_status;
98859599516SKenneth E. Jansen char read_out_tag[MAX_FIELDS_NAME_LENGTH];
98959599516SKenneth E. Jansen int j;
99059599516SKenneth E. Jansen int magic_number;
99159599516SKenneth E. Jansen
99259599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
99359599516SKenneth E. Jansen MPI_File_read_at( PhastaIOActiveFiles[i]->file_handle,
99459599516SKenneth E. Jansen 0,
99559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header,
99659599516SKenneth E. Jansen MasterHeaderSize,
99759599516SKenneth E. Jansen MPI_CHAR,
99859599516SKenneth E. Jansen &read_tag_status );
99959599516SKenneth E. Jansen }
100059599516SKenneth E. Jansen
100159599516SKenneth E. Jansen MPI_Bcast( PhastaIOActiveFiles[i]->master_header,
100259599516SKenneth E. Jansen MasterHeaderSize,
100359599516SKenneth E. Jansen MPI_CHAR,
100459599516SKenneth E. Jansen 0,
100559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm );
100659599516SKenneth E. Jansen
100759599516SKenneth E. Jansen memcpy( read_out_tag,
100859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header,
100959599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH-1 );
101059599516SKenneth E. Jansen
101159599516SKenneth E. Jansen if ( cscompare ("MPI_IO_Tag",read_out_tag) )
101259599516SKenneth E. Jansen {
101359599516SKenneth E. Jansen // Test endianess ...
101459599516SKenneth E. Jansen memcpy ( &magic_number,
101559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
101659599516SKenneth E. Jansen sizeof(int) ); // masterheader should look like "MPI_IO_Tag : 12180 " with 12180 in binary format
101759599516SKenneth E. Jansen
101859599516SKenneth E. Jansen if ( magic_number != ENDIAN_TEST_NUMBER )
101959599516SKenneth E. Jansen {
102059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->Wrong_Endian = true;
102159599516SKenneth E. Jansen }
102259599516SKenneth E. Jansen
102359599516SKenneth E. Jansen memcpy( read_out_tag,
102459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH+1, // TODO: WHY +1???
102559599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH );
102659599516SKenneth E. Jansen
102759599516SKenneth E. Jansen // Read in # fields ...
102859599516SKenneth E. Jansen token = strtok ( read_out_tag, ":" );
102959599516SKenneth E. Jansen token = strtok( NULL," ,;<>" );
103059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFields = atoi( token );
103159599516SKenneth E. Jansen
10322dd307a1SCameron Smith unsigned long **header_table;
10332dd307a1SCameron Smith header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
103459599516SKenneth E. Jansen
103559599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
103659599516SKenneth E. Jansen {
10372dd307a1SCameron Smith header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
103859599516SKenneth E. Jansen }
103959599516SKenneth E. Jansen
104059599516SKenneth E. Jansen // Read in the offset table ...
104159599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ )
104259599516SKenneth E. Jansen {
104359599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
104459599516SKenneth E. Jansen memcpy( header_table[j],
104559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header +
104659599516SKenneth E. Jansen VERSION_INFO_HEADER_SIZE +
10472dd307a1SCameron Smith j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
10482dd307a1SCameron Smith PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
104959599516SKenneth E. Jansen }
105059599516SKenneth E. Jansen
105159599516SKenneth E. Jansen MPI_Scatter( header_table[j],
105259599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp,
105359599516SKenneth E. Jansen MPI_LONG_LONG_INT,
105459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table[j],
105559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp,
105659599516SKenneth E. Jansen MPI_LONG_LONG_INT,
105759599516SKenneth E. Jansen 0,
105859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm );
105959599516SKenneth E. Jansen
106059599516SKenneth E. Jansen // Swap byte order if endianess is different ...
106159599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->Wrong_Endian ) {
106259599516SKenneth E. Jansen SwapArrayByteOrder( PhastaIOActiveFiles[i]->my_read_table[j],
10634187599aSCameron Smith sizeof(unsigned long),
106459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp );
106559599516SKenneth E. Jansen }
106659599516SKenneth E. Jansen }
106759599516SKenneth E. Jansen
106859599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
106959599516SKenneth E. Jansen free ( header_table[j] );
107059599516SKenneth E. Jansen }
107159599516SKenneth E. Jansen free (header_table);
107259599516SKenneth E. Jansen
107359599516SKenneth E. Jansen } // end of if MPI_IO_TAG
107459599516SKenneth E. Jansen else //else not valid MPI file
107559599516SKenneth E. Jansen {
107659599516SKenneth E. Jansen *fileDescriptor = NOT_A_MPI_FILE;
107759599516SKenneth 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);
107859599516SKenneth E. Jansen //Printing MasterHeaderSize is useful to test a compiler bug on Intrepid BGP
107959599516SKenneth E. Jansen endTimer(&timer_end);
108059599516SKenneth E. Jansen printPerf("openfile", timer_start, timer_end, 0, 0, "");
108159599516SKenneth E. Jansen return;
108259599516SKenneth E. Jansen }
108359599516SKenneth E. Jansen } // end of if "read"
108459599516SKenneth E. Jansen else if( cscompare( "write", imode ) )
108559599516SKenneth E. Jansen {
10868f9016f6SCameron Smith PHASTAIO_OPENTIME(
108759599516SKenneth E. Jansen rc = MPI_File_open( PhastaIOActiveFiles[i]->local_comm,
108859599516SKenneth E. Jansen fname,
108959599516SKenneth E. Jansen MPI_MODE_WRONLY | MPI_MODE_CREATE,
109059599516SKenneth E. Jansen MPI_INFO_NULL,
109159599516SKenneth E. Jansen &(PhastaIOActiveFiles[i]->file_handle) );
10928f9016f6SCameron Smith )
10934d60bba2SCameron Smith if(rc != MPI_SUCCESS)
109459599516SKenneth E. Jansen {
109559599516SKenneth E. Jansen *fileDescriptor = UNABLE_TO_OPEN_FILE;
10964d60bba2SCameron Smith int error_string_length;
10974d60bba2SCameron Smith char error_string[4096];
10984d60bba2SCameron Smith MPI_Error_string(rc, error_string, &error_string_length);
10994d60bba2SCameron Smith fprintf(stderr, "Error openfile: Unable to open file %s! MPI reports \"%s\"\n",fname,error_string);
110059599516SKenneth E. Jansen return;
110159599516SKenneth E. Jansen }
110259599516SKenneth E. Jansen } // end of if "write"
110359599516SKenneth E. Jansen free (fname);
110459599516SKenneth E. Jansen free (imode);
110559599516SKenneth E. Jansen } // end of if FileIndex != 0
110659599516SKenneth E. Jansen
110759599516SKenneth E. Jansen endTimer(&timer_end);
110859599516SKenneth E. Jansen printPerf("openfile", timer_start, timer_end, 0, 0, "");
110959599516SKenneth E. Jansen }
111059599516SKenneth E. Jansen
111159599516SKenneth E. Jansen /** close file for both POSIX and MPI-IO syncIO format.
111259599516SKenneth E. Jansen *
111359599516SKenneth E. Jansen * If it's old POSIX format, simply call posix fclose().
111459599516SKenneth E. Jansen *
111559599516SKenneth E. Jansen * If it's MPI-IO foramt:
111659599516SKenneth E. Jansen * in "read" mode, it simply close file with MPI-IO close routine.
111759599516SKenneth E. Jansen * in "write" mode, rank 0 in each group will re-assemble the master header and
111859599516SKenneth E. Jansen * offset table and write to the beginning of file, then close the file.
111959599516SKenneth E. Jansen */
closefile(int * fileDescriptor,const char mode[])1120103be424SCameron Smith void closefile( int* fileDescriptor, const char mode[] )
112159599516SKenneth E. Jansen {
112259599516SKenneth E. Jansen double timer_start, timer_end;
112359599516SKenneth E. Jansen startTimer(&timer_start);
112459599516SKenneth E. Jansen
112559599516SKenneth E. Jansen int i = *fileDescriptor;
112659599516SKenneth E. Jansen checkFileDescriptor("closefile",&i);
112759599516SKenneth E. Jansen
112859599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) {
112959599516SKenneth E. Jansen char* imode = StringStripper( mode );
113059599516SKenneth E. Jansen
113159599516SKenneth E. Jansen if( cscompare( "write", imode )
113259599516SKenneth E. Jansen || cscompare( "append", imode ) ) {
113359599516SKenneth E. Jansen fflush( fileArray[ *fileDescriptor - 1 ] );
113459599516SKenneth E. Jansen }
113559599516SKenneth E. Jansen
11368f9016f6SCameron Smith PHASTAIO_CLOSETIME(fclose( fileArray[ *fileDescriptor - 1 ] );)
113759599516SKenneth E. Jansen free (imode);
113859599516SKenneth E. Jansen }
113959599516SKenneth E. Jansen else {
114059599516SKenneth E. Jansen char* imode = StringStripper( mode );
114159599516SKenneth E. Jansen
114259599516SKenneth E. Jansen //write master header here:
114359599516SKenneth E. Jansen if ( cscompare( "write", imode ) ) {
114459599516SKenneth E. Jansen // if ( PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields < 2*ONE_MEGABYTE/8 ) //SHOULD BE CHECKED
114559599516SKenneth E. Jansen // MasterHeaderSize = 4*ONE_MEGABYTE;
114659599516SKenneth E. Jansen // else
114759599516SKenneth E. Jansen // MasterHeaderSize = 4*ONE_MEGABYTE + PhastaIOActiveFiles[i]->nPPF * PhastaIOActiveFiles[i]->nFields * 8 - 2*ONE_MEGABYTE;
114859599516SKenneth E. Jansen
114959599516SKenneth E. Jansen MasterHeaderSize = computeMHSize( PhastaIOActiveFiles[i]->nFields, PhastaIOActiveFiles[i]->nPPF, LATEST_WRITE_VERSION);
115059599516SKenneth E. Jansen phprintf_0("Info closefile: myrank = %d, MasterHeaderSize = %d\n", PhastaIOActiveFiles[i]->myrank, MasterHeaderSize);
115159599516SKenneth E. Jansen
115259599516SKenneth E. Jansen MPI_Status write_header_status;
115359599516SKenneth E. Jansen char mpi_tag[MAX_FIELDS_NAME_LENGTH];
115459599516SKenneth E. Jansen char version[MAX_FIELDS_NAME_LENGTH/4];
115559599516SKenneth E. Jansen char mhsize[MAX_FIELDS_NAME_LENGTH/4];
115659599516SKenneth E. Jansen int magic_number = ENDIAN_TEST_NUMBER;
115759599516SKenneth E. Jansen
115859599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->local_myrank == 0 )
115959599516SKenneth E. Jansen {
116059599516SKenneth E. Jansen bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
116159599516SKenneth E. Jansen sprintf(mpi_tag, "MPI_IO_Tag : ");
116259599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->master_header,
116359599516SKenneth E. Jansen mpi_tag,
116459599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH);
116559599516SKenneth E. Jansen
116659599516SKenneth E. Jansen bzero((void*)version,MAX_FIELDS_NAME_LENGTH/4);
116759599516SKenneth E. Jansen // this version is "1", print version in ASCII
116859599516SKenneth E. Jansen sprintf(version, "version : %d",1);
116959599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/2,
117059599516SKenneth E. Jansen version,
117159599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH/4);
117259599516SKenneth E. Jansen
117359599516SKenneth E. Jansen // master header size is computed using the formula above
117459599516SKenneth E. Jansen bzero((void*)mhsize,MAX_FIELDS_NAME_LENGTH/4);
117559599516SKenneth E. Jansen sprintf(mhsize, "mhsize : ");
117659599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->master_header + MAX_FIELDS_NAME_LENGTH/4*3,
117759599516SKenneth E. Jansen mhsize,
117859599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH/4);
117959599516SKenneth E. Jansen
118059599516SKenneth E. Jansen bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
118159599516SKenneth E. Jansen sprintf(mpi_tag,
118259599516SKenneth E. Jansen "\nnFields : %d\n",
118359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFields);
118459599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->master_header+MAX_FIELDS_NAME_LENGTH,
118559599516SKenneth E. Jansen mpi_tag,
118659599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH);
118759599516SKenneth E. Jansen
118859599516SKenneth E. Jansen bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
118959599516SKenneth E. Jansen sprintf(mpi_tag, "\nnPPF : %d\n", PhastaIOActiveFiles[i]->nPPF);
119059599516SKenneth E. Jansen memcpy( PhastaIOActiveFiles[i]->master_header+
119159599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nFields *
119259599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH +
119359599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH * 2,
119459599516SKenneth E. Jansen mpi_tag,
119559599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH);
119659599516SKenneth E. Jansen
119759599516SKenneth E. Jansen memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("MPI_IO_Tag : ")-1, //-1 sizeof returns the size of the string+1 for "\0"
119859599516SKenneth E. Jansen &magic_number,
119959599516SKenneth E. Jansen sizeof(int));
120059599516SKenneth E. Jansen
120159599516SKenneth E. Jansen memcpy( PhastaIOActiveFiles[i]->master_header+sizeof("mhsize : ") -1 + MAX_FIELDS_NAME_LENGTH/4*3,
120259599516SKenneth E. Jansen &MasterHeaderSize,
120359599516SKenneth E. Jansen sizeof(int));
120459599516SKenneth E. Jansen }
120559599516SKenneth E. Jansen
120659599516SKenneth E. Jansen int j = 0;
12072dd307a1SCameron Smith unsigned long **header_table;
12082dd307a1SCameron Smith header_table = ( unsigned long ** )calloc(PhastaIOActiveFiles[i]->nFields, sizeof(unsigned long *));
120959599516SKenneth E. Jansen
121059599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
12112dd307a1SCameron Smith header_table[j]=( unsigned long * ) calloc( PhastaIOActiveFiles[i]->nPPF , sizeof( unsigned long));
121259599516SKenneth E. Jansen }
121359599516SKenneth E. Jansen
121459599516SKenneth E. Jansen //if( irank == 0 ) printf("gonna mpi_gather, myrank = %d\n", irank);
121559599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
121659599516SKenneth E. Jansen MPI_Gather( PhastaIOActiveFiles[i]->my_offset_table[j],
121759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp,
121859599516SKenneth E. Jansen MPI_LONG_LONG_INT,
121959599516SKenneth E. Jansen header_table[j],
122059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp,
122159599516SKenneth E. Jansen MPI_LONG_LONG_INT,
122259599516SKenneth E. Jansen 0,
122359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm );
122459599516SKenneth E. Jansen }
122559599516SKenneth E. Jansen
122659599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
122759599516SKenneth E. Jansen
122859599516SKenneth E. Jansen //if( irank == 0 ) printf("gonna memcpy for every procs, myrank = %d\n", irank);
122959599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
123059599516SKenneth E. Jansen memcpy ( PhastaIOActiveFiles[i]->master_header +
123159599516SKenneth E. Jansen VERSION_INFO_HEADER_SIZE +
12322dd307a1SCameron Smith j * PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long),
123359599516SKenneth E. Jansen header_table[j],
12342dd307a1SCameron Smith PhastaIOActiveFiles[i]->nPPF * sizeof(unsigned long) );
123559599516SKenneth E. Jansen }
123659599516SKenneth E. Jansen
123759599516SKenneth E. Jansen //if( irank == 0 ) printf("gonna file_write_at(), myrank = %d\n", irank);
123859599516SKenneth E. Jansen MPI_File_write_at( PhastaIOActiveFiles[i]->file_handle,
123959599516SKenneth E. Jansen 0,
124059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header,
124159599516SKenneth E. Jansen MasterHeaderSize,
124259599516SKenneth E. Jansen MPI_CHAR,
124359599516SKenneth E. Jansen &write_header_status );
124459599516SKenneth E. Jansen }
124559599516SKenneth E. Jansen
124659599516SKenneth E. Jansen ////free(PhastaIOActiveFiles[i]->master_header);
124759599516SKenneth E. Jansen
124859599516SKenneth E. Jansen for ( j = 0; j < PhastaIOActiveFiles[i]->nFields; j++ ) {
124959599516SKenneth E. Jansen free ( header_table[j] );
125059599516SKenneth E. Jansen }
125159599516SKenneth E. Jansen free (header_table);
125259599516SKenneth E. Jansen }
125359599516SKenneth E. Jansen
125459599516SKenneth E. Jansen //if( irank == 0 ) printf("gonna file_close(), myrank = %d\n", irank);
12558f9016f6SCameron Smith PHASTAIO_CLOSETIME(
125659599516SKenneth E. Jansen MPI_File_close( &( PhastaIOActiveFiles[i]->file_handle ) );
12578f9016f6SCameron Smith )
125859599516SKenneth E. Jansen free ( imode );
125959599516SKenneth E. Jansen }
126059599516SKenneth E. Jansen
126159599516SKenneth E. Jansen endTimer(&timer_end);
126259599516SKenneth E. Jansen printPerf("closefile_", timer_start, timer_end, 0, 0, "");
126359599516SKenneth E. Jansen }
126459599516SKenneth E. Jansen
readHeader(FILE * f,const char phrase[],int * params,int numParams,const char * iotype)12653872e963SCameron Smith int readHeader( FILE* f, const char phrase[],
12663872e963SCameron Smith int* params, int numParams, const char* iotype) {
12673872e963SCameron Smith isBinary(iotype);
12683872e963SCameron Smith return readHeader(f,phrase,params,numParams);
12693872e963SCameron Smith }
12703872e963SCameron Smith
readheader(int * fileDescriptor,const char keyphrase[],void * valueArray,int * nItems,const char datatype[],const char iotype[])1271103be424SCameron Smith void readheader(
1272103be424SCameron Smith int* fileDescriptor,
127359599516SKenneth E. Jansen const char keyphrase[],
127459599516SKenneth E. Jansen void* valueArray,
127559599516SKenneth E. Jansen int* nItems,
127659599516SKenneth E. Jansen const char datatype[],
127759599516SKenneth E. Jansen const char iotype[] )
127859599516SKenneth E. Jansen {
127959599516SKenneth E. Jansen double timer_start, timer_end;
1280d3337298SCameron Smith
128159599516SKenneth E. Jansen startTimer(&timer_start);
128259599516SKenneth E. Jansen
128359599516SKenneth E. Jansen int i = *fileDescriptor;
128459599516SKenneth E. Jansen checkFileDescriptor("readheader",&i);
128559599516SKenneth E. Jansen
128659599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) {
128759599516SKenneth E. Jansen int filePtr = *fileDescriptor - 1;
128859599516SKenneth E. Jansen FILE* fileObject;
128959599516SKenneth E. Jansen int* valueListInt;
129059599516SKenneth E. Jansen
129159599516SKenneth E. Jansen if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
129259599516SKenneth E. Jansen fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
129359599516SKenneth E. Jansen fprintf(stderr,"openfile_ function has to be called before \n") ;
129459599516SKenneth E. Jansen fprintf(stderr,"acessing the file\n ") ;
129559599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
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
1301961a4ff6SCameron Smith LastHeaderKey[filePtr] = std::string(keyphrase);
130259599516SKenneth E. Jansen LastHeaderNotFound = false;
130359599516SKenneth E. Jansen
130459599516SKenneth E. Jansen fileObject = fileArray[ filePtr ] ;
130559599516SKenneth E. Jansen Wrong_Endian = byte_order[ filePtr ];
130659599516SKenneth E. Jansen
130759599516SKenneth E. Jansen isBinary( iotype );
130859599516SKenneth E. Jansen typeSize( datatype ); //redundant call, just avoid a compiler warning.
130959599516SKenneth E. Jansen
131059599516SKenneth E. Jansen // right now we are making the assumption that we will only write integers
131159599516SKenneth E. Jansen // on the header line.
131259599516SKenneth E. Jansen
131359599516SKenneth E. Jansen valueListInt = static_cast< int* >( valueArray );
131459599516SKenneth E. Jansen int ierr = readHeader( fileObject ,
131559599516SKenneth E. Jansen keyphrase,
131659599516SKenneth E. Jansen valueListInt,
131759599516SKenneth E. Jansen *nItems ) ;
131859599516SKenneth E. Jansen
131959599516SKenneth E. Jansen byte_order[ filePtr ] = Wrong_Endian ;
132059599516SKenneth E. Jansen
132159599516SKenneth E. Jansen if ( ierr ) LastHeaderNotFound = true;
132259599516SKenneth E. Jansen
132359599516SKenneth E. Jansen //return ; // don't return, go to the end to print perf
132459599516SKenneth E. Jansen }
132559599516SKenneth E. Jansen else {
132659599516SKenneth E. Jansen int* valueListInt;
132759599516SKenneth E. Jansen valueListInt = static_cast <int*>(valueArray);
1328400e9fc0SCameron Smith char* token = NULL;
132959599516SKenneth E. Jansen bool FOUND = false ;
133059599516SKenneth E. Jansen isBinary( iotype );
133159599516SKenneth E. Jansen
133259599516SKenneth E. Jansen MPI_Status read_offset_status;
133359599516SKenneth E. Jansen char read_out_tag[MAX_FIELDS_NAME_LENGTH];
1334400e9fc0SCameron Smith memset(read_out_tag, '\0', MAX_FIELDS_NAME_LENGTH);
133559599516SKenneth E. Jansen char readouttag[MAX_FIELDS_NUMBER][MAX_FIELDS_NAME_LENGTH];
133659599516SKenneth E. Jansen int j;
133759599516SKenneth E. Jansen
133859599516SKenneth E. Jansen int string_length = strlen( keyphrase );
133959599516SKenneth E. Jansen char* buffer = (char*) malloc ( string_length+1 );
134059599516SKenneth E. Jansen
134159599516SKenneth E. Jansen strcpy ( buffer, keyphrase );
134259599516SKenneth E. Jansen buffer[ string_length ] = '\0';
134359599516SKenneth E. Jansen
134459599516SKenneth E. Jansen char* st2 = strtok ( buffer, "@" );
134559599516SKenneth E. Jansen st2 = strtok (NULL, "@");
134659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->GPid = atoi(st2);
134759599516SKenneth E. Jansen if ( char* p = strpbrk(buffer, "@") )
134859599516SKenneth E. Jansen *p = '\0';
134959599516SKenneth E. Jansen
135059599516SKenneth E. Jansen // Check if the user has input the right GPid
135159599516SKenneth E. Jansen if ( ( PhastaIOActiveFiles[i]->GPid <=
135259599516SKenneth E. Jansen PhastaIOActiveFiles[i]->myrank *
135359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp )||
135459599516SKenneth E. Jansen ( PhastaIOActiveFiles[i]->GPid >
135559599516SKenneth E. Jansen ( PhastaIOActiveFiles[i]->myrank + 1 ) *
135659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->nppp ) )
135759599516SKenneth E. Jansen {
135859599516SKenneth E. Jansen *fileDescriptor = NOT_A_MPI_FILE;
135959599516SKenneth 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);
136059599516SKenneth E. Jansen // It is possible atoi could not generate a clear integer from st2 because of additional garbage character in keyphrase
136159599516SKenneth E. Jansen endTimer(&timer_end);
136259599516SKenneth E. Jansen printPerf("readheader", timer_start, timer_end, 0, 0, "");
136359599516SKenneth E. Jansen return;
136459599516SKenneth E. Jansen }
136559599516SKenneth E. Jansen
136659599516SKenneth E. Jansen // Find the field we want ...
136759599516SKenneth E. Jansen //for ( j = 0; j<MAX_FIELDS_NUMBER; j++ )
136859599516SKenneth E. Jansen for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
136959599516SKenneth E. Jansen {
137059599516SKenneth E. Jansen memcpy( readouttag[j],
137159599516SKenneth E. Jansen PhastaIOActiveFiles[i]->master_header + j*MAX_FIELDS_NAME_LENGTH+MAX_FIELDS_NAME_LENGTH*2+1,
137259599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH-1 );
137359599516SKenneth E. Jansen }
137459599516SKenneth E. Jansen
137559599516SKenneth E. Jansen for ( j = 0; j<PhastaIOActiveFiles[i]->nFields; j++ )
137659599516SKenneth E. Jansen {
137759599516SKenneth E. Jansen token = strtok ( readouttag[j], ":" );
137859599516SKenneth E. Jansen
137959599516SKenneth E. Jansen //if ( cscompare( buffer, token ) )
138059599516SKenneth E. Jansen if ( cscompare( token , buffer ) && cscompare( buffer, token ) )
138159599516SKenneth 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").
138259599516SKenneth 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).
138359599516SKenneth 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.
138459599516SKenneth E. Jansen {
138559599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_field_count = j;
138659599516SKenneth E. Jansen FOUND = true;
138759599516SKenneth E. Jansen //printf("buffer: %s | token: %s | j: %d\n",buffer,token,j);
138859599516SKenneth E. Jansen break;
138959599516SKenneth E. Jansen }
139059599516SKenneth E. Jansen }
139159599516SKenneth E. Jansen free(buffer);
139259599516SKenneth E. Jansen
139359599516SKenneth E. Jansen if (!FOUND)
139459599516SKenneth E. Jansen {
139559599516SKenneth E. Jansen //if(irank==0) printf("Warning readheader: Not found %s \n",keyphrase); //PhastaIOActiveFiles[i]->myrank is certainly initialized here.
139659599516SKenneth E. Jansen if(PhastaIOActiveFiles[i]->myrank == 0) printf("WARNING readheader: Not found %s\n",keyphrase);
139759599516SKenneth E. Jansen endTimer(&timer_end);
139859599516SKenneth E. Jansen printPerf("readheader", timer_start, timer_end, 0, 0, "");
139959599516SKenneth E. Jansen return;
140059599516SKenneth E. Jansen }
140159599516SKenneth E. Jansen
140259599516SKenneth E. Jansen // Find the part we want ...
140359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->read_part_count = PhastaIOActiveFiles[i]->GPid -
140459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->myrank * PhastaIOActiveFiles[i]->nppp - 1;
140559599516SKenneth E. Jansen
140659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset =
140759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_read_table[PhastaIOActiveFiles[i]->read_field_count][PhastaIOActiveFiles[i]->read_part_count];
140859599516SKenneth E. Jansen
140959599516SKenneth E. Jansen // printf("****Rank %d offset is %d\n",PhastaIOActiveFiles[i]->myrank,PhastaIOActiveFiles[i]->my_offset);
141059599516SKenneth E. Jansen
141159599516SKenneth E. Jansen // Read each datablock header here ...
141259599516SKenneth E. Jansen
141359599516SKenneth E. Jansen MPI_File_read_at_all( PhastaIOActiveFiles[i]->file_handle,
141459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset+1,
141559599516SKenneth E. Jansen read_out_tag,
141659599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH-1,
141759599516SKenneth E. Jansen MPI_CHAR,
141859599516SKenneth E. Jansen &read_offset_status );
141959599516SKenneth E. Jansen token = strtok ( read_out_tag, ":" );
142059599516SKenneth E. Jansen
142159599516SKenneth E. Jansen // printf("&&&&Rank %d read_out_tag is %s\n",PhastaIOActiveFiles[i]->myrank,read_out_tag);
142259599516SKenneth E. Jansen
142359599516SKenneth 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.
142459599516SKenneth E. Jansen {
142559599516SKenneth E. Jansen FOUND = true ;
142659599516SKenneth E. Jansen token = strtok( NULL, " ,;<>" );
142759599516SKenneth E. Jansen for( j=0; j < *nItems && ( token = strtok( NULL," ,;<>") ); j++ )
142859599516SKenneth E. Jansen valueListInt[j] = atoi( token );
142959599516SKenneth E. Jansen
143059599516SKenneth E. Jansen if ( j < *nItems )
143159599516SKenneth E. Jansen {
143259599516SKenneth E. Jansen fprintf( stderr, "Expected # of ints not found for: %s\n", keyphrase );
143359599516SKenneth E. Jansen }
143459599516SKenneth E. Jansen }
143559599516SKenneth E. Jansen else {
143659599516SKenneth E. Jansen //if(irank==0)
143759599516SKenneth E. Jansen if(PhastaIOActiveFiles[i]->myrank == 0)
143859599516SKenneth E. Jansen // If we enter this if, there is a problem with the name of some fields
143959599516SKenneth E. Jansen {
144059599516SKenneth E. Jansen printf("Error readheader: Unexpected mismatch between keyphrase = %s and token = %s\n",keyphrase,token);
144159599516SKenneth E. Jansen }
144259599516SKenneth E. Jansen }
144359599516SKenneth E. Jansen }
144459599516SKenneth E. Jansen
144559599516SKenneth E. Jansen endTimer(&timer_end);
144659599516SKenneth E. Jansen printPerf("readheader", timer_start, timer_end, 0, 0, "");
144759599516SKenneth E. Jansen
144859599516SKenneth E. Jansen }
144959599516SKenneth E. Jansen
readDataBlock(FILE * fileObject,void * valueArray,int nItems,const char datatype[],const char iotype[])14503956dcfeSCameron Smith void readDataBlock(
14513956dcfeSCameron Smith FILE* fileObject,
14523956dcfeSCameron Smith void* valueArray,
14533956dcfeSCameron Smith int nItems,
14543956dcfeSCameron Smith const char datatype[],
14553956dcfeSCameron Smith const char iotype[] )
14563956dcfeSCameron Smith {
14573956dcfeSCameron Smith isBinary(iotype);
14583956dcfeSCameron Smith size_t type_size = typeSize( datatype );
14598f9016f6SCameron Smith phastaioTime t0,t1;
1460389f1f6aSCameron Smith phastaio_time(&t0);
14613956dcfeSCameron Smith if ( binary_format ) {
14623956dcfeSCameron Smith char junk = '\0';
14633956dcfeSCameron Smith fread( valueArray, type_size, nItems, fileObject );
14643956dcfeSCameron Smith fread( &junk, sizeof(char), 1 , fileObject );
14653956dcfeSCameron Smith if ( Wrong_Endian ) SwapArrayByteOrder( valueArray, type_size, nItems );
14663956dcfeSCameron Smith } else {
14673956dcfeSCameron Smith char* ts1 = StringStripper( datatype );
14683956dcfeSCameron Smith if ( cscompare( "integer", ts1 ) ) {
14693956dcfeSCameron Smith for( int n=0; n < nItems ; n++ )
14703956dcfeSCameron Smith fscanf(fileObject, "%d\n",(int*)((int*)valueArray+n) );
14713956dcfeSCameron Smith } else if ( cscompare( "double", ts1 ) ) {
14723956dcfeSCameron Smith for( int n=0; n < nItems ; n++ )
14733956dcfeSCameron Smith fscanf(fileObject, "%lf\n",(double*)((double*)valueArray+n) );
14743956dcfeSCameron Smith }
14753956dcfeSCameron Smith free (ts1);
14763956dcfeSCameron Smith }
1477389f1f6aSCameron Smith phastaio_time(&t1);
1478389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1479389f1f6aSCameron Smith phastaio_addReadTime(elapsed);
1480f42e0444SCameron Smith phastaio_addReadBytes(nItems*type_size);
14813956dcfeSCameron Smith }
14823956dcfeSCameron Smith
readdatablock(int * fileDescriptor,const char keyphrase[],void * valueArray,int * nItems,const char datatype[],const char iotype[])1483103be424SCameron Smith void readdatablock(
1484103be424SCameron Smith int* fileDescriptor,
148559599516SKenneth E. Jansen const char keyphrase[],
148659599516SKenneth E. Jansen void* valueArray,
148759599516SKenneth E. Jansen int* nItems,
148859599516SKenneth E. Jansen const char datatype[],
148959599516SKenneth E. Jansen const char iotype[] )
149059599516SKenneth E. Jansen {
149159599516SKenneth E. Jansen //if(irank == 0) printf("entering readdatablock()\n");
14922dd307a1SCameron Smith unsigned long data_size = 0;
149359599516SKenneth E. Jansen double timer_start, timer_end;
149459599516SKenneth E. Jansen startTimer(&timer_start);
149559599516SKenneth E. Jansen
149659599516SKenneth E. Jansen int i = *fileDescriptor;
149759599516SKenneth E. Jansen checkFileDescriptor("readdatablock",&i);
149859599516SKenneth E. Jansen
149959599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) {
150059599516SKenneth E. Jansen int filePtr = *fileDescriptor - 1;
150159599516SKenneth E. Jansen FILE* fileObject;
150259599516SKenneth E. Jansen
150359599516SKenneth E. Jansen if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
150459599516SKenneth E. Jansen fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
150559599516SKenneth E. Jansen fprintf(stderr,"openfile_ function has to be called before\n") ;
150659599516SKenneth E. Jansen fprintf(stderr,"acessing the file\n ") ;
150759599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
150859599516SKenneth E. Jansen endTimer(&timer_end);
150959599516SKenneth E. Jansen printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
151059599516SKenneth E. Jansen return;
151159599516SKenneth E. Jansen }
151259599516SKenneth E. Jansen
151359599516SKenneth E. Jansen // error check..
151459599516SKenneth E. Jansen // since we require that a consistant header always preceed the data block
151559599516SKenneth E. Jansen // let us check to see that it is actually the case.
151659599516SKenneth E. Jansen
1517961a4ff6SCameron Smith if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
151859599516SKenneth E. Jansen fprintf(stderr, "Header not consistant with data block\n");
1519961a4ff6SCameron Smith fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
152059599516SKenneth E. Jansen fprintf(stderr, "DataBlock: %s\n ", keyphrase );
152159599516SKenneth E. Jansen fprintf(stderr, "Please recheck read sequence \n");
152259599516SKenneth E. Jansen if( Strict_Error ) {
152359599516SKenneth E. Jansen fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
152459599516SKenneth E. Jansen endTimer(&timer_end);
152559599516SKenneth E. Jansen printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
152659599516SKenneth E. Jansen return;
152759599516SKenneth E. Jansen }
152859599516SKenneth E. Jansen }
152959599516SKenneth E. Jansen
153059599516SKenneth E. Jansen if ( LastHeaderNotFound ) {
153159599516SKenneth E. Jansen endTimer(&timer_end);
153259599516SKenneth E. Jansen printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
153359599516SKenneth E. Jansen return;
153459599516SKenneth E. Jansen }
153559599516SKenneth E. Jansen fileObject = fileArray[ filePtr ];
153659599516SKenneth E. Jansen Wrong_Endian = byte_order[ filePtr ];
1537bae968b9SCameron Smith LastHeaderKey.erase(filePtr);
15383956dcfeSCameron Smith readDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
153959599516SKenneth E. Jansen
154059599516SKenneth E. Jansen //return;
154159599516SKenneth E. Jansen }
154259599516SKenneth E. Jansen else {
154359599516SKenneth E. Jansen // printf("read data block\n");
154459599516SKenneth E. Jansen MPI_Status read_data_status;
154559599516SKenneth E. Jansen size_t type_size = typeSize( datatype );
154659599516SKenneth E. Jansen int nUnits = *nItems;
154759599516SKenneth E. Jansen isBinary( iotype );
154859599516SKenneth E. Jansen
154959599516SKenneth E. Jansen // read datablock then
155059599516SKenneth E. Jansen //MR CHANGE
155159599516SKenneth E. Jansen // if ( cscompare ( datatype, "double"))
155259599516SKenneth E. Jansen char* ts2 = StringStripper( datatype );
155359599516SKenneth E. Jansen if ( cscompare ( "double" , ts2))
155459599516SKenneth E. Jansen //MR CHANGE END
155559599516SKenneth E. Jansen {
155659599516SKenneth E. Jansen
15578f9016f6SCameron Smith phastaioTime t0,t1;
1558389f1f6aSCameron Smith phastaio_time(&t0);
155959599516SKenneth E. Jansen MPI_File_read_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
156059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
156159599516SKenneth E. Jansen valueArray,
156259599516SKenneth E. Jansen nUnits,
156359599516SKenneth E. Jansen MPI_DOUBLE );
156459599516SKenneth E. Jansen MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
156559599516SKenneth E. Jansen valueArray,
156659599516SKenneth E. Jansen &read_data_status );
156759599516SKenneth E. Jansen data_size=8*nUnits;
1568389f1f6aSCameron Smith phastaio_time(&t1);
1569389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1570389f1f6aSCameron Smith phastaio_addReadTime(elapsed);
1571f42e0444SCameron Smith phastaio_addReadBytes(nUnits*sizeof(double));
157259599516SKenneth E. Jansen }
157359599516SKenneth E. Jansen //MR CHANGE
157459599516SKenneth E. Jansen // else if ( cscompare ( datatype, "integer"))
157559599516SKenneth E. Jansen else if ( cscompare ( "integer" , ts2))
157659599516SKenneth E. Jansen //MR CHANGE END
157759599516SKenneth E. Jansen {
15788f9016f6SCameron Smith phastaioTime t0,t1;
1579389f1f6aSCameron Smith phastaio_time(&t0);
158059599516SKenneth E. Jansen MPI_File_read_at_all_begin(PhastaIOActiveFiles[i]->file_handle,
158159599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset + DB_HEADER_SIZE,
158259599516SKenneth E. Jansen valueArray,
158359599516SKenneth E. Jansen nUnits,
158459599516SKenneth E. Jansen MPI_INT );
158559599516SKenneth E. Jansen MPI_File_read_at_all_end( PhastaIOActiveFiles[i]->file_handle,
158659599516SKenneth E. Jansen valueArray,
158759599516SKenneth E. Jansen &read_data_status );
158859599516SKenneth E. Jansen data_size=4*nUnits;
1589389f1f6aSCameron Smith phastaio_time(&t1);
1590389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1591389f1f6aSCameron Smith phastaio_addReadTime(elapsed);
1592f42e0444SCameron Smith phastaio_addReadBytes(nUnits*sizeof(int));
159359599516SKenneth E. Jansen }
159459599516SKenneth E. Jansen else
159559599516SKenneth E. Jansen {
159659599516SKenneth E. Jansen *fileDescriptor = DATA_TYPE_ILLEGAL;
159759599516SKenneth E. Jansen printf("readdatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
159859599516SKenneth E. Jansen endTimer(&timer_end);
159959599516SKenneth E. Jansen printPerf("readdatablock", timer_start, timer_end, 0, 0, "");
160059599516SKenneth E. Jansen return;
160159599516SKenneth E. Jansen }
160259599516SKenneth E. Jansen free(ts2);
160359599516SKenneth E. Jansen
160459599516SKenneth E. Jansen
160559599516SKenneth E. Jansen // printf("%d Read finishe\n",PhastaIOActiveFiles[i]->myrank);
160659599516SKenneth E. Jansen
160759599516SKenneth E. Jansen // Swap data byte order if endianess is different ...
160859599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->Wrong_Endian )
160959599516SKenneth E. Jansen {
161059599516SKenneth E. Jansen SwapArrayByteOrder( valueArray, type_size, nUnits );
161159599516SKenneth E. Jansen }
161259599516SKenneth E. Jansen }
161359599516SKenneth E. Jansen
161459599516SKenneth E. Jansen endTimer(&timer_end);
161559599516SKenneth E. Jansen char extra_msg[1024];
161659599516SKenneth E. Jansen memset(extra_msg, '\0', 1024);
161759599516SKenneth E. Jansen char* key = StringStripper(keyphrase);
161859599516SKenneth E. Jansen sprintf(extra_msg, " field is %s ", key);
161959599516SKenneth E. Jansen printPerf("readdatablock", timer_start, timer_end, data_size, 1, extra_msg);
162059599516SKenneth E. Jansen free(key);
162159599516SKenneth E. Jansen
162259599516SKenneth E. Jansen }
162359599516SKenneth E. Jansen
writeHeader(FILE * f,const char keyphrase[],const void * valueArray,const int nItems,const int ndataItems,const char datatype[],const char iotype[])1624103be424SCameron Smith void writeHeader(
1625103be424SCameron Smith FILE* f,
1626ea868eb1SCameron Smith const char keyphrase[],
1627ea868eb1SCameron Smith const void* valueArray,
1628ea868eb1SCameron Smith const int nItems,
1629ea868eb1SCameron Smith const int ndataItems,
1630ea868eb1SCameron Smith const char datatype[],
1631ea868eb1SCameron Smith const char iotype[])
1632ea868eb1SCameron Smith {
1633ea868eb1SCameron Smith isBinary( iotype );
1634ea868eb1SCameron Smith
1635ea868eb1SCameron Smith const int _newline =
1636ea868eb1SCameron Smith ( ndataItems > 0 ) ? sizeof( char ) : 0;
1637ea868eb1SCameron Smith int size_of_nextblock =
1638ea868eb1SCameron Smith ( binary_format ) ? typeSize(datatype) * ndataItems + _newline : ndataItems;
1639ea868eb1SCameron Smith
1640ea868eb1SCameron Smith fprintf( f, "%s : < %d > ", keyphrase, size_of_nextblock );
1641ea868eb1SCameron Smith for( int i = 0; i < nItems; i++ )
1642ea868eb1SCameron Smith fprintf( f, "%d ", *((int*)((int*)valueArray+i)));
1643ea868eb1SCameron Smith fprintf( f, "\n");
1644ea868eb1SCameron Smith }
1645ea868eb1SCameron Smith
writeheader(const int * fileDescriptor,const char keyphrase[],const void * valueArray,const int * nItems,const int * ndataItems,const char datatype[],const char iotype[])1646103be424SCameron Smith void writeheader(
1647103be424SCameron Smith const int* fileDescriptor,
164859599516SKenneth E. Jansen const char keyphrase[],
164959599516SKenneth E. Jansen const void* valueArray,
165059599516SKenneth E. Jansen const int* nItems,
165159599516SKenneth E. Jansen const int* ndataItems,
165259599516SKenneth E. Jansen const char datatype[],
165359599516SKenneth E. Jansen const char iotype[])
165459599516SKenneth E. Jansen {
165559599516SKenneth E. Jansen
165659599516SKenneth E. Jansen //if(irank == 0) printf("entering writeheader()\n");
165759599516SKenneth E. Jansen
165859599516SKenneth E. Jansen double timer_start, timer_end;
165959599516SKenneth E. Jansen startTimer(&timer_start);
166059599516SKenneth E. Jansen
166159599516SKenneth E. Jansen int i = *fileDescriptor;
166259599516SKenneth E. Jansen checkFileDescriptor("writeheader",&i);
166359599516SKenneth E. Jansen
166459599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) {
166559599516SKenneth E. Jansen int filePtr = *fileDescriptor - 1;
166659599516SKenneth E. Jansen if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
166759599516SKenneth E. Jansen fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
166859599516SKenneth E. Jansen fprintf(stderr,"openfile_ function has to be called before \n") ;
166959599516SKenneth E. Jansen fprintf(stderr,"acessing the file\n ") ;
167059599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
167159599516SKenneth E. Jansen endTimer(&timer_end);
167259599516SKenneth E. Jansen printPerf("writeheader", timer_start, timer_end, 0, 0, "");
167359599516SKenneth E. Jansen return;
167459599516SKenneth E. Jansen }
167559599516SKenneth E. Jansen
1676961a4ff6SCameron Smith LastHeaderKey[filePtr] = std::string(keyphrase);
167759599516SKenneth E. Jansen DataSize = *ndataItems;
1678ea868eb1SCameron Smith FILE* fileObject = fileArray[ filePtr ] ;
1679ea868eb1SCameron Smith header_type[ filePtr ] = typeSize( datatype );
1680ea868eb1SCameron Smith writeHeader(fileObject,keyphrase,valueArray,*nItems,
1681ea868eb1SCameron Smith *ndataItems,datatype,iotype);
168259599516SKenneth E. Jansen }
168359599516SKenneth E. Jansen else { // else it's parallel I/O
168459599516SKenneth E. Jansen DataSize = *ndataItems;
168559599516SKenneth E. Jansen size_t type_size = typeSize( datatype );
168659599516SKenneth E. Jansen isBinary( iotype );
168759599516SKenneth E. Jansen char mpi_tag[MAX_FIELDS_NAME_LENGTH];
168859599516SKenneth E. Jansen
168959599516SKenneth E. Jansen int string_length = strlen( keyphrase );
169059599516SKenneth E. Jansen char* buffer = (char*) malloc ( string_length+1 );
169159599516SKenneth E. Jansen
169259599516SKenneth E. Jansen strcpy ( buffer, keyphrase);
169359599516SKenneth E. Jansen buffer[ string_length ] = '\0';
169459599516SKenneth E. Jansen
169559599516SKenneth E. Jansen char* st2 = strtok ( buffer, "@" );
169659599516SKenneth E. Jansen st2 = strtok (NULL, "@");
169759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->GPid = atoi(st2);
169859599516SKenneth E. Jansen
169959599516SKenneth E. Jansen if ( char* p = strpbrk(buffer, "@") )
170059599516SKenneth E. Jansen *p = '\0';
1701*8aafd871SPranav Subramanian assert(PhastaIOActiveFiles[i]->field_count < MAX_FIELDS_NUMBER);
1702*8aafd871SPranav Subramanian assert(PhastaIOActiveFiles[i]->part_count < PhastaIOActiveFiles[i]->nppp);
170359599516SKenneth E. Jansen bzero((void*)mpi_tag,MAX_FIELDS_NAME_LENGTH);
170459599516SKenneth E. Jansen sprintf(mpi_tag, "\n%s : %d\n", buffer, PhastaIOActiveFiles[i]->field_count);
17052dd307a1SCameron Smith unsigned long offset_value;
170659599516SKenneth E. Jansen
170759599516SKenneth E. Jansen int temp = *ndataItems;
17082dd307a1SCameron Smith unsigned long number_of_items = (unsigned long)temp;
170959599516SKenneth E. Jansen MPI_Barrier(PhastaIOActiveFiles[i]->local_comm);
171059599516SKenneth E. Jansen
171159599516SKenneth E. Jansen MPI_Scan( &number_of_items,
171259599516SKenneth E. Jansen &offset_value,
171359599516SKenneth E. Jansen 1,
171459599516SKenneth E. Jansen MPI_LONG_LONG_INT,
171559599516SKenneth E. Jansen MPI_SUM,
171659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm );
171759599516SKenneth E. Jansen
171859599516SKenneth E. Jansen offset_value = (offset_value - number_of_items) * type_size;
171959599516SKenneth E. Jansen
172059599516SKenneth E. Jansen offset_value += PhastaIOActiveFiles[i]->local_myrank *
172159599516SKenneth E. Jansen DB_HEADER_SIZE +
172259599516SKenneth E. Jansen PhastaIOActiveFiles[i]->next_start_address;
172359599516SKenneth E. Jansen // This offset is the starting address of each datablock header...
172459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset = offset_value;
172559599516SKenneth E. Jansen
172659599516SKenneth E. Jansen // Write in my offset table ...
1727*8aafd871SPranav Subramanian PhastaIOActiveFiles[i]->my_offset_table[PhastaIOActiveFiles[i]->field_count][PhastaIOActiveFiles[i]->part_count] = PhastaIOActiveFiles[i]->my_offset;
172859599516SKenneth E. Jansen
172959599516SKenneth E. Jansen // Update the next-start-address ...
173059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->next_start_address = offset_value +
173159599516SKenneth E. Jansen number_of_items * type_size +
173259599516SKenneth E. Jansen DB_HEADER_SIZE;
173359599516SKenneth E. Jansen MPI_Bcast( &(PhastaIOActiveFiles[i]->next_start_address),
173459599516SKenneth E. Jansen 1,
173559599516SKenneth E. Jansen MPI_LONG_LONG_INT,
173659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_numprocs-1,
173759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->local_comm );
173859599516SKenneth E. Jansen
173959599516SKenneth E. Jansen // Prepare datablock header ...
174059599516SKenneth E. Jansen int _newline = (*ndataItems>0)?sizeof(char):0;
174159599516SKenneth E. Jansen unsigned int size_of_nextblock = type_size * (*ndataItems) + _newline;
174259599516SKenneth E. Jansen
174359599516SKenneth E. Jansen //char datablock_header[255];
174459599516SKenneth E. Jansen //bzero((void*)datablock_header,255);
174559599516SKenneth E. Jansen char datablock_header[DB_HEADER_SIZE];
174659599516SKenneth E. Jansen bzero((void*)datablock_header,DB_HEADER_SIZE);
174759599516SKenneth E. Jansen
174859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->GPid = PhastaIOActiveFiles[i]->nppp*PhastaIOActiveFiles[i]->myrank+PhastaIOActiveFiles[i]->part_count;
174959599516SKenneth E. Jansen sprintf( datablock_header,
175059599516SKenneth E. Jansen "\n%s : < %u >",
175159599516SKenneth E. Jansen keyphrase,
175259599516SKenneth E. Jansen size_of_nextblock );
175359599516SKenneth E. Jansen
175459599516SKenneth E. Jansen for ( int j = 0; j < *nItems; j++ )
175559599516SKenneth E. Jansen {
175659599516SKenneth E. Jansen sprintf( datablock_header,
175759599516SKenneth E. Jansen "%s %d ",
175859599516SKenneth E. Jansen datablock_header,
175959599516SKenneth E. Jansen *((int*)((int*)valueArray+j)));
176059599516SKenneth E. Jansen }
176159599516SKenneth E. Jansen sprintf( datablock_header,
176259599516SKenneth E. Jansen "%s\n ",
176359599516SKenneth E. Jansen datablock_header );
176459599516SKenneth E. Jansen
176559599516SKenneth E. Jansen // Write datablock header ...
176659599516SKenneth E. Jansen //MR CHANGE
176759599516SKenneth E. Jansen // if ( cscompare(datatype,"double") )
176859599516SKenneth E. Jansen char* ts1 = StringStripper( datatype );
176959599516SKenneth E. Jansen if ( cscompare("double",ts1) )
177059599516SKenneth E. Jansen //MR CHANGE END
177159599516SKenneth E. Jansen {
177259599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->double_chunk );
177359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->double_chunk = ( double * )malloc( (sizeof( double )*number_of_items+ DB_HEADER_SIZE));
177459599516SKenneth E. Jansen
177559599516SKenneth E. Jansen double * aa = ( double * )datablock_header;
177659599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->double_chunk, aa, DB_HEADER_SIZE);
177759599516SKenneth E. Jansen }
177859599516SKenneth E. Jansen //MR CHANGE
177959599516SKenneth E. Jansen // if ( cscompare(datatype,"integer") )
178059599516SKenneth E. Jansen else if ( cscompare("integer",ts1) )
178159599516SKenneth E. Jansen //MR CHANGE END
178259599516SKenneth E. Jansen {
178359599516SKenneth E. Jansen free ( PhastaIOActiveFiles[i]->int_chunk );
178459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->int_chunk = ( int * )malloc( (sizeof( int )*number_of_items+ DB_HEADER_SIZE));
178559599516SKenneth E. Jansen
178659599516SKenneth E. Jansen int * aa = ( int * )datablock_header;
178759599516SKenneth E. Jansen memcpy(PhastaIOActiveFiles[i]->int_chunk, aa, DB_HEADER_SIZE);
178859599516SKenneth E. Jansen }
178959599516SKenneth E. Jansen else {
179059599516SKenneth E. Jansen // *fileDescriptor = DATA_TYPE_ILLEGAL;
179159599516SKenneth E. Jansen printf("writeheader - DATA_TYPE_ILLEGAL - %s\n",datatype);
179259599516SKenneth E. Jansen endTimer(&timer_end);
179359599516SKenneth E. Jansen printPerf("writeheader", timer_start, timer_end, 0, 0, "");
179459599516SKenneth E. Jansen return;
179559599516SKenneth E. Jansen }
179659599516SKenneth E. Jansen free(ts1);
179759599516SKenneth E. Jansen
179859599516SKenneth E. Jansen PhastaIOActiveFiles[i]->part_count++;
179959599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->part_count == PhastaIOActiveFiles[i]->nppp ) {
180059599516SKenneth E. Jansen //A new field will be written
180159599516SKenneth E. Jansen if ( PhastaIOActiveFiles[i]->local_myrank == 0 ) {
180259599516SKenneth E. Jansen memcpy( PhastaIOActiveFiles[i]->master_header +
180359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->field_count *
180459599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH +
180559599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH * 2,
180659599516SKenneth E. Jansen mpi_tag,
180759599516SKenneth E. Jansen MAX_FIELDS_NAME_LENGTH-1);
180859599516SKenneth E. Jansen }
180959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->field_count++;
181059599516SKenneth E. Jansen PhastaIOActiveFiles[i]->part_count=0;
181159599516SKenneth E. Jansen }
181259599516SKenneth E. Jansen free(buffer);
181359599516SKenneth E. Jansen }
181459599516SKenneth E. Jansen
181559599516SKenneth E. Jansen endTimer(&timer_end);
181659599516SKenneth E. Jansen printPerf("writeheader", timer_start, timer_end, 0, 0, "");
181759599516SKenneth E. Jansen }
181859599516SKenneth E. Jansen
writeDataBlock(FILE * f,const void * valueArray,const int nItems,const char datatype[],const char iotype[])1819103be424SCameron Smith void writeDataBlock(
1820103be424SCameron Smith FILE* f,
1821ea868eb1SCameron Smith const void* valueArray,
1822ea868eb1SCameron Smith const int nItems,
1823ea868eb1SCameron Smith const char datatype[],
1824ea868eb1SCameron Smith const char iotype[] )
1825ea868eb1SCameron Smith {
1826ea868eb1SCameron Smith isBinary( iotype );
1827ea868eb1SCameron Smith size_t type_size = typeSize( datatype );
18288f9016f6SCameron Smith phastaioTime t0,t1;
1829389f1f6aSCameron Smith phastaio_time(&t0);
1830ea868eb1SCameron Smith if ( binary_format ) {
1831ea868eb1SCameron Smith fwrite( valueArray, type_size, nItems, f );
1832ea868eb1SCameron Smith fprintf( f,"\n");
1833ea868eb1SCameron Smith } else {
1834ea868eb1SCameron Smith char* ts1 = StringStripper( datatype );
1835ea868eb1SCameron Smith if ( cscompare( "integer", ts1 ) ) {
1836be3da47bSCameron Smith const int* vals = (int*) valueArray;
1837ea868eb1SCameron Smith for( int n=0; n < nItems ; n++ )
1838be3da47bSCameron Smith fprintf(f,"%d\n",vals[n]);
1839ea868eb1SCameron Smith } else if ( cscompare( "double", ts1 ) ) {
1840be3da47bSCameron Smith const double* vals = (double*) valueArray;
1841ea868eb1SCameron Smith for( int n=0; n < nItems ; n++ )
1842be3da47bSCameron Smith fprintf(f,"%f\n",vals[n]);
1843ea868eb1SCameron Smith }
1844ea868eb1SCameron Smith free (ts1);
1845ea868eb1SCameron Smith }
1846389f1f6aSCameron Smith phastaio_time(&t1);
1847389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1848389f1f6aSCameron Smith phastaio_addWriteTime(elapsed);
1849f42e0444SCameron Smith phastaio_addWriteBytes(nItems*type_size);
1850ea868eb1SCameron Smith }
1851ea868eb1SCameron Smith
writedatablock(const int * fileDescriptor,const char keyphrase[],const void * valueArray,const int * nItems,const char datatype[],const char iotype[])1852103be424SCameron Smith void writedatablock(
1853103be424SCameron Smith const int* fileDescriptor,
185459599516SKenneth E. Jansen const char keyphrase[],
185559599516SKenneth E. Jansen const void* valueArray,
185659599516SKenneth E. Jansen const int* nItems,
185759599516SKenneth E. Jansen const char datatype[],
185859599516SKenneth E. Jansen const char iotype[] )
185959599516SKenneth E. Jansen {
186059599516SKenneth E. Jansen //if(irank == 0) printf("entering writedatablock()\n");
186159599516SKenneth E. Jansen
18622dd307a1SCameron Smith unsigned long data_size = 0;
186359599516SKenneth E. Jansen double timer_start, timer_end;
186459599516SKenneth E. Jansen startTimer(&timer_start);
186559599516SKenneth E. Jansen
186659599516SKenneth E. Jansen int i = *fileDescriptor;
186759599516SKenneth E. Jansen checkFileDescriptor("writedatablock",&i);
186859599516SKenneth E. Jansen
186959599516SKenneth E. Jansen if ( PhastaIONextActiveIndex == 0 ) {
187059599516SKenneth E. Jansen int filePtr = *fileDescriptor - 1;
187159599516SKenneth E. Jansen
187259599516SKenneth E. Jansen if ( *fileDescriptor < 1 || *fileDescriptor > (int)fileArray.size() ) {
187359599516SKenneth E. Jansen fprintf(stderr,"No file associated with Descriptor %d\n",*fileDescriptor);
187459599516SKenneth E. Jansen fprintf(stderr,"openfile_ function has to be called before \n") ;
187559599516SKenneth E. Jansen fprintf(stderr,"acessing the file\n ") ;
187659599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n");
187759599516SKenneth E. Jansen endTimer(&timer_end);
187859599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
187959599516SKenneth E. Jansen return;
188059599516SKenneth E. Jansen }
188159599516SKenneth E. Jansen // since we require that a consistant header always preceed the data block
188259599516SKenneth E. Jansen // let us check to see that it is actually the case.
188359599516SKenneth E. Jansen
1884961a4ff6SCameron Smith if ( ! cscompare( LastHeaderKey[ filePtr ].c_str(), keyphrase ) ) {
188559599516SKenneth E. Jansen fprintf(stderr, "Header not consistant with data block\n");
1886961a4ff6SCameron Smith fprintf(stderr, "Header: %s\n", LastHeaderKey[ filePtr ].c_str() );
188759599516SKenneth E. Jansen fprintf(stderr, "DataBlock: %s\n ", keyphrase );
188859599516SKenneth E. Jansen fprintf(stderr, "Please recheck write sequence \n");
188959599516SKenneth E. Jansen if( Strict_Error ) {
189059599516SKenneth E. Jansen fprintf(stderr, "fatal error: cannot continue, returning out of call\n");
189159599516SKenneth E. Jansen endTimer(&timer_end);
189259599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
189359599516SKenneth E. Jansen return;
189459599516SKenneth E. Jansen }
189559599516SKenneth E. Jansen }
189659599516SKenneth E. Jansen
189759599516SKenneth E. Jansen FILE* fileObject = fileArray[ filePtr ] ;
189859599516SKenneth E. Jansen size_t type_size=typeSize( datatype );
189959599516SKenneth E. Jansen isBinary( iotype );
190059599516SKenneth E. Jansen
1901bae968b9SCameron Smith LastHeaderKey.erase(filePtr);
1902bae968b9SCameron Smith
190359599516SKenneth E. Jansen if ( header_type[filePtr] != (int)type_size ) {
190459599516SKenneth E. Jansen fprintf(stderr,"header and datablock differ on typeof data in the block for\n");
190559599516SKenneth E. Jansen fprintf(stderr,"keyphrase : %s\n", keyphrase);
190659599516SKenneth E. Jansen if( Strict_Error ) {
190759599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
190859599516SKenneth E. Jansen endTimer(&timer_end);
190959599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
191059599516SKenneth E. Jansen return;
191159599516SKenneth E. Jansen }
191259599516SKenneth E. Jansen }
191359599516SKenneth E. Jansen
191459599516SKenneth E. Jansen int nUnits = *nItems;
191559599516SKenneth E. Jansen
191659599516SKenneth E. Jansen if ( nUnits != DataSize ) {
191759599516SKenneth E. Jansen fprintf(stderr,"header and datablock differ on number of data items for\n");
191859599516SKenneth E. Jansen fprintf(stderr,"keyphrase : %s\n", keyphrase);
191959599516SKenneth E. Jansen if( Strict_Error ) {
192059599516SKenneth E. Jansen fprintf(stderr,"fatal error: cannot continue, returning out of call\n" );
192159599516SKenneth E. Jansen endTimer(&timer_end);
192259599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
192359599516SKenneth E. Jansen return;
192459599516SKenneth E. Jansen }
192559599516SKenneth E. Jansen }
1926ea868eb1SCameron Smith writeDataBlock(fileObject,valueArray,*nItems,datatype,iotype);
192759599516SKenneth E. Jansen }
192859599516SKenneth E. Jansen else { // syncIO case
192959599516SKenneth E. Jansen MPI_Status write_data_status;
193059599516SKenneth E. Jansen isBinary( iotype );
193159599516SKenneth E. Jansen int nUnits = *nItems;
193259599516SKenneth E. Jansen
193359599516SKenneth E. Jansen //MR CHANGE
193459599516SKenneth E. Jansen // if ( cscompare(datatype,"double") )
193559599516SKenneth E. Jansen char* ts1 = StringStripper( datatype );
193659599516SKenneth E. Jansen if ( cscompare("double",ts1) )
193759599516SKenneth E. Jansen //MR CHANGE END
193859599516SKenneth E. Jansen {
193959599516SKenneth E. Jansen memcpy((PhastaIOActiveFiles[i]->double_chunk+DB_HEADER_SIZE/sizeof(double)), valueArray, nUnits*sizeof(double));
19408f9016f6SCameron Smith phastaioTime t0,t1;
1941389f1f6aSCameron Smith phastaio_time(&t0);
194259599516SKenneth E. Jansen MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
194359599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset,
194459599516SKenneth E. Jansen PhastaIOActiveFiles[i]->double_chunk,
194559599516SKenneth E. Jansen //BLOCK_SIZE/sizeof(double),
194659599516SKenneth E. Jansen nUnits+DB_HEADER_SIZE/sizeof(double),
194759599516SKenneth E. Jansen MPI_DOUBLE );
194859599516SKenneth E. Jansen MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
194959599516SKenneth E. Jansen PhastaIOActiveFiles[i]->double_chunk,
195059599516SKenneth E. Jansen &write_data_status );
195159599516SKenneth E. Jansen data_size=8*nUnits;
1952389f1f6aSCameron Smith phastaio_time(&t1);
1953389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1954389f1f6aSCameron Smith phastaio_addWriteTime(elapsed);
1955f42e0444SCameron Smith phastaio_addWriteBytes((nUnits*sizeof(double))+DB_HEADER_SIZE);
195659599516SKenneth E. Jansen }
195759599516SKenneth E. Jansen //MR CHANGE
195859599516SKenneth E. Jansen // else if ( cscompare ( datatype, "integer"))
195959599516SKenneth E. Jansen else if ( cscompare("integer",ts1) )
196059599516SKenneth E. Jansen //MR CHANGE END
196159599516SKenneth E. Jansen {
196259599516SKenneth E. Jansen memcpy((PhastaIOActiveFiles[i]->int_chunk+DB_HEADER_SIZE/sizeof(int)), valueArray, nUnits*sizeof(int));
19638f9016f6SCameron Smith phastaioTime t0,t1;
1964389f1f6aSCameron Smith phastaio_time(&t0);
196559599516SKenneth E. Jansen MPI_File_write_at_all_begin( PhastaIOActiveFiles[i]->file_handle,
196659599516SKenneth E. Jansen PhastaIOActiveFiles[i]->my_offset,
196759599516SKenneth E. Jansen PhastaIOActiveFiles[i]->int_chunk,
196859599516SKenneth E. Jansen nUnits+DB_HEADER_SIZE/sizeof(int),
196959599516SKenneth E. Jansen MPI_INT );
197059599516SKenneth E. Jansen MPI_File_write_at_all_end( PhastaIOActiveFiles[i]->file_handle,
197159599516SKenneth E. Jansen PhastaIOActiveFiles[i]->int_chunk,
197259599516SKenneth E. Jansen &write_data_status );
197359599516SKenneth E. Jansen data_size=4*nUnits;
1974389f1f6aSCameron Smith phastaio_time(&t1);
1975389f1f6aSCameron Smith const size_t elapsed = phastaio_time_diff(&t0,&t1);
1976389f1f6aSCameron Smith phastaio_addWriteTime(elapsed);
1977f42e0444SCameron Smith phastaio_addWriteBytes((nUnits*sizeof(int))+DB_HEADER_SIZE);
197859599516SKenneth E. Jansen }
197959599516SKenneth E. Jansen else {
198059599516SKenneth E. Jansen printf("Error: writedatablock - DATA_TYPE_ILLEGAL - %s\n",datatype);
198159599516SKenneth E. Jansen endTimer(&timer_end);
198259599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, 0, 0, "");
198359599516SKenneth E. Jansen return;
198459599516SKenneth E. Jansen }
198559599516SKenneth E. Jansen free(ts1);
198659599516SKenneth E. Jansen }
198759599516SKenneth E. Jansen
198859599516SKenneth E. Jansen endTimer(&timer_end);
198959599516SKenneth E. Jansen char extra_msg[1024];
199059599516SKenneth E. Jansen memset(extra_msg, '\0', 1024);
199159599516SKenneth E. Jansen char* key = StringStripper(keyphrase);
199259599516SKenneth E. Jansen sprintf(extra_msg, " field is %s ", key);
199359599516SKenneth E. Jansen printPerf("writedatablock", timer_start, timer_end, data_size, 1, extra_msg);
199459599516SKenneth E. Jansen free(key);
199559599516SKenneth E. Jansen
199659599516SKenneth E. Jansen }
199759599516SKenneth E. Jansen
SwapArrayByteOrder(void * array,int nbytes,int nItems)1998103be424SCameron Smith void SwapArrayByteOrder( void* array, int nbytes, int nItems )
199959599516SKenneth E. Jansen {
200059599516SKenneth E. Jansen /* This swaps the byte order for the array of nItems each
200159599516SKenneth E. Jansen of size nbytes , This will be called only locally */
200259599516SKenneth E. Jansen int i,j;
200359599516SKenneth E. Jansen unsigned char* ucDst = (unsigned char*)array;
200459599516SKenneth E. Jansen for(i=0; i < nItems; i++) {
200559599516SKenneth E. Jansen for(j=0; j < (nbytes/2); j++)
200659599516SKenneth E. Jansen std::swap( ucDst[j] , ucDst[(nbytes - 1) - j] );
200759599516SKenneth E. Jansen ucDst += nbytes;
200859599516SKenneth E. Jansen }
200959599516SKenneth E. Jansen }
201059599516SKenneth E. Jansen
writestring(int * fileDescriptor,const char inString[])2011103be424SCameron Smith void writestring( int* fileDescriptor, const char inString[] )
201259599516SKenneth E. Jansen {
201359599516SKenneth E. Jansen int filePtr = *fileDescriptor - 1;
201459599516SKenneth E. Jansen FILE* fileObject = fileArray[filePtr] ;
201559599516SKenneth E. Jansen fprintf(fileObject,"%s",inString );
201659599516SKenneth E. Jansen return;
201759599516SKenneth E. Jansen }
201859599516SKenneth E. Jansen
Gather_Headers(int * fileDescriptor,vector<string> & headers)2019103be424SCameron Smith void Gather_Headers( int* fileDescriptor, vector< string >& headers )
202059599516SKenneth E. Jansen {
202159599516SKenneth E. Jansen FILE* fileObject;
202259599516SKenneth E. Jansen char Line[1024];
202359599516SKenneth E. Jansen
202459599516SKenneth E. Jansen fileObject = fileArray[ (*fileDescriptor)-1 ];
202559599516SKenneth E. Jansen
202659599516SKenneth E. Jansen while( !feof(fileObject) ) {
202759599516SKenneth E. Jansen fgets( Line, 1024, fileObject);
202859599516SKenneth E. Jansen if ( Line[0] == '#' ) {
202959599516SKenneth E. Jansen headers.push_back( Line );
203059599516SKenneth E. Jansen } else {
203159599516SKenneth E. Jansen break;
203259599516SKenneth E. Jansen }
203359599516SKenneth E. Jansen }
203459599516SKenneth E. Jansen rewind( fileObject );
203559599516SKenneth E. Jansen clearerr( fileObject );
203659599516SKenneth E. Jansen }
203759599516SKenneth E. Jansen
isWrong(void)2038103be424SCameron Smith void isWrong( void ) {
2039103be424SCameron Smith (Wrong_Endian) ? fprintf(stdout,"YES\n") : fprintf(stdout,"NO\n");
2040103be424SCameron Smith }
204159599516SKenneth E. Jansen
togglestrictmode(void)2042103be424SCameron Smith void togglestrictmode( void ) { Strict_Error = !Strict_Error; }
2043103be424SCameron Smith
isLittleEndian(void)2044103be424SCameron Smith int isLittleEndian( void )
204559599516SKenneth E. Jansen {
204659599516SKenneth E. Jansen // this function returns a 1 if the current running architecture is
204759599516SKenneth E. Jansen // LittleEndian Byte Ordered, else it returns a zero
204859599516SKenneth E. Jansen union {
204959599516SKenneth E. Jansen long a;
205059599516SKenneth E. Jansen char c[sizeof( long )];
205159599516SKenneth E. Jansen } endianUnion;
205259599516SKenneth E. Jansen endianUnion.a = 1 ;
205359599516SKenneth E. Jansen if ( endianUnion.c[sizeof(long)-1] != 1 ) return 1 ;
205459599516SKenneth E. Jansen else return 0;
205559599516SKenneth E. Jansen }
205659599516SKenneth E. Jansen
205759599516SKenneth E. Jansen namespace PHASTA {
205859599516SKenneth E. Jansen const char* const PhastaIO_traits<int>::type_string = "integer";
205959599516SKenneth E. Jansen const char* const PhastaIO_traits<double>::type_string = "double";
206059599516SKenneth E. Jansen }
2061