xref: /phasta/phSolver/common/test/phIOreaddatablock.cc (revision 4a32e26d9a154ddd316ffc54cbda2fac2e6f6ed1)
1*4a32e26dSCameron Smith #include <mpi.h>
2*4a32e26dSCameron Smith #include <stdio.h>
3*4a32e26dSCameron Smith #include <stdlib.h>
4*4a32e26dSCameron Smith #include <unistd.h>
5*4a32e26dSCameron Smith #include "phIO.h"
6*4a32e26dSCameron Smith 
7*4a32e26dSCameron Smith int main(int argc, char* argv[]) {
8*4a32e26dSCameron Smith   MPI_Init(&argc,&argv);
9*4a32e26dSCameron Smith   int rank;
10*4a32e26dSCameron Smith   MPI_Comm_rank(MPI_COMM_WORLD, &rank);
11*4a32e26dSCameron Smith   if( argc != 2 ) {
12*4a32e26dSCameron Smith     fprintf(stderr, "Usage: %s <numSyncFiles>\n",argv[0]);
13*4a32e26dSCameron Smith     MPI_Finalize();
14*4a32e26dSCameron Smith     return 1;
15*4a32e26dSCameron Smith   }
16*4a32e26dSCameron Smith   int nfiles[2] = {atoi(argv[1]), 1};
17*4a32e26dSCameron Smith   const char* phrase = "co-ordinates";
18*4a32e26dSCameron Smith   const char* type = "double";
19*4a32e26dSCameron Smith   const char* iotype = "binary";
20*4a32e26dSCameron Smith   const char* dir[2] = {"4-procs_case-SyncIO-2", "4-procs_case-Posix"};
21*4a32e26dSCameron Smith   const char* filename[2] = {"geombc-dat.", "geombc.dat."};
22*4a32e26dSCameron Smith   double* coords[2] = {NULL, NULL};
23*4a32e26dSCameron Smith   int len[2] = {0, 0};
24*4a32e26dSCameron Smith   int numpts[2];
25*4a32e26dSCameron Smith   phio_fp file;
26*4a32e26dSCameron Smith   int two = 2;
27*4a32e26dSCameron Smith   for(int i=0; i<2; i++) {
28*4a32e26dSCameron Smith     chdir(dir[i]);
29*4a32e26dSCameron Smith     MPI_Barrier(MPI_COMM_WORLD);
30*4a32e26dSCameron Smith     phio_openfile_read(filename[i], &(nfiles[i]), &file);
31*4a32e26dSCameron Smith     phio_readheader(file, phrase, numpts, &two, type, iotype);
32*4a32e26dSCameron Smith     len[i] = numpts[0]*3; //numPts * 3 dimensions
33*4a32e26dSCameron Smith     fprintf(stderr, "%d %s len %d\n", rank, __func__, len[i]);
34*4a32e26dSCameron Smith     coords[i] = (double*) malloc(len[i]*sizeof(double));
35*4a32e26dSCameron Smith     phio_readdatablock(file, phrase, coords[i], &(len[i]), type, iotype);
36*4a32e26dSCameron Smith     phio_closefile_read(file);
37*4a32e26dSCameron Smith     chdir("..");
38*4a32e26dSCameron Smith     MPI_Barrier(MPI_COMM_WORLD);
39*4a32e26dSCameron Smith     if(!rank)
40*4a32e26dSCameron Smith       fprintf(stderr, "-------------------\n");
41*4a32e26dSCameron Smith   }
42*4a32e26dSCameron Smith   int match = (len[0] == len[1]);
43*4a32e26dSCameron Smith   if(!rank && match)
44*4a32e26dSCameron Smith     fprintf(stderr, "number of points match!\n");
45*4a32e26dSCameron Smith   if(!rank && !match) {
46*4a32e26dSCameron Smith     fprintf(stderr, "number of points don't match... :(\n");
47*4a32e26dSCameron Smith     return 1;
48*4a32e26dSCameron Smith   }
49*4a32e26dSCameron Smith 
50*4a32e26dSCameron Smith   match = true;
51*4a32e26dSCameron Smith   for(int i=0; i<len[0]; i++)
52*4a32e26dSCameron Smith     match = match && ( coords[0][i] == coords[1][i] );
53*4a32e26dSCameron Smith 
54*4a32e26dSCameron Smith   if(!rank && match)
55*4a32e26dSCameron Smith     fprintf(stderr, "points match!\n");
56*4a32e26dSCameron Smith   if(!rank && !match) {
57*4a32e26dSCameron Smith     fprintf(stderr, "points don't match... :(\n");
58*4a32e26dSCameron Smith     return 1;
59*4a32e26dSCameron Smith   }
60*4a32e26dSCameron Smith 
61*4a32e26dSCameron Smith   for(int i=0; i<2; i++)
62*4a32e26dSCameron Smith     free(coords[i]);
63*4a32e26dSCameron Smith 
64*4a32e26dSCameron Smith   MPI_Finalize();
65*4a32e26dSCameron Smith   return !match;
66*4a32e26dSCameron Smith }
67