src/hdf5/H5ecloud/Bench.c

Go to the documentation of this file.
00001 #include <stdio.h>
00002 #include <stdlib.h>
00003 #include <string.h>
00004 #include <mpi.h>
00005 /* #include <mpio.h> */
00006 #include <unistd.h>
00007 #include <sys/types.h>
00008 #ifndef PARALLEL_IO
00009 #define PARALLEL_IO
00010 #endif
00011 #include "H5Part.h"
00012 
00013 #define FILENAME "testio"
00014 /* normally 64 steps for real benchmark */
00015 #define NSTEPS 64
00016 /* normally 51e6 for real benchmark */
00017 #define NPARTICLES 51e6
00018 #define NTRIALS 3
00019 
00020 int main(int argc,char *argv[]){
00021   MPI_Info info;
00022   int npdims=1;
00023   int nprocs,rank;
00024   int trial;
00025   int i,j,n; /* iteration variables */
00026   double starttime,curtime, endtime;
00027   int nparticles = NPARTICLES;
00028   double *x,*y,*z,*px,*py,*pz;
00029   typedef double *ddouble;
00030   ddouble data[6];
00031   int64_t *id;
00032   MPI_Datatype chunktype;
00033   int offset;
00034   int localnp;
00035   int err1,err2;
00036   char filename[128]; /*= FILENAME; */
00037   H5PartFile *f;
00038   char newfilename[128];
00039   FILE *fd;
00040   MPI_File file;
00041   MPI_Info bogusinfo;
00042   MPI_Offset foffset;
00043 
00044   MPI_Comm comm,dcomm = MPI_COMM_WORLD;
00045   
00046   MPI_Init(&argc,&argv);
00047   MPI_Comm_rank(dcomm,&rank);
00048   MPI_Comm_size(dcomm,&nprocs);
00049 
00050   localnp=nparticles/(int64_t)nprocs;
00051   for(offset=0,i=0;i<rank;i++){
00052     offset+=localnp;
00053   }
00054   
00055   data[0]=x=(double*)malloc(sizeof(double)*(size_t)localnp);
00056   data[1]=y=(double*)malloc(sizeof(double)*(size_t)localnp);
00057   data[2]=z=(double*)malloc(sizeof(double)*(size_t)localnp);
00058   data[3]=px=(double*)malloc(sizeof(double)*(size_t)localnp);
00059   data[4]=py=(double*)malloc(sizeof(double)*(size_t)localnp);
00060   data[5]=pz=(double*)malloc(sizeof(double)*(size_t)localnp);
00061   
00062   /* printf("about to call create subarray with nparticles=%u localnp=%u offset=%u\n",
00063         nparticles,localnp,offset); */
00064   MPI_Type_create_subarray(1, /* rank */
00065                            &nparticles, /* size of the global array */
00066                            &localnp, /* size of my local chunk */
00067                            &offset, /* offset of this chunk in global */
00068                            MPI_ORDER_FORTRAN, /* fortran storage order */
00069                            MPI_DOUBLE,
00070                            &chunktype);
00071   MPI_Type_commit(&chunktype);
00072   MPI_Info_create(&info);
00073   if(rank==0) printf("Nprocs=%u Particles=%u*6attribs*sizeof(double) Particles/proc=%u Nsteps=%u Ntrials=%u\n",
00074         nprocs,nparticles,localnp,NSTEPS,NTRIALS);
00075         
00076 
00077   for(trial=0;trial<NTRIALS;trial++){
00078   if(rank==0) printf("---------------------- Trial %u of %u ---------------------\n",trial+1,NTRIALS);
00079   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00080   sprintf(filename,"%s.%u.mpio.dat",FILENAME,nprocs);
00081   if(rank==0) unlink(filename);
00082   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00083   MPI_File_open(MPI_COMM_WORLD,filename,
00084         MPI_MODE_CREATE | MPI_MODE_RDWR,
00085         info,&file);
00086         
00087   MPI_File_set_view(file,0,MPI_DOUBLE,chunktype,"native",info);
00088   /* now a barrier to get the start timers roughly synced*/
00089   MPI_Barrier(MPI_COMM_WORLD);
00090   curtime = starttime = MPI_Wtime();
00091   endtime = starttime+5.0*60.0; /* end in 5 minutes */
00092   MPI_Bcast(&endtime,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
00093   /* must touch the entire array after each write */
00094   /* ensures cache-invalidation */
00095   foffset=0;
00096   i=0;
00097   curtime=starttime;
00098   for(i=0;i<NSTEPS;i++){
00099     int n;
00100     MPI_Status status;
00101     for(j=0;j<6;j++){
00102       /* touch data */
00103       for(n=0;n<localnp;n++)
00104         (data[j])[n]=(double)rank;
00105       /* write to that file */
00106       /*  MPI_File_set_view(file,foffset,MPI_DOUBLE,chunktype,"native",info);*/
00107       MPI_File_write_at_all(file,
00108                             foffset,
00109                             data[j],
00110                             localnp,
00111                             MPI_DOUBLE,&status);
00112       foffset+=nparticles;
00113     }
00114     curtime=MPI_Wtime(); /* ensure no race condition by broadcasting time */
00115     MPI_Bcast(&curtime,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
00116   }
00117   MPI_File_close(&file);
00118   MPI_Barrier(MPI_COMM_WORLD);
00119   endtime=MPI_Wtime();
00120   sprintf(filename,"%s.%u.h5.dat",FILENAME,nprocs);
00121   if(rank==0){
00122     puts("*");
00123     unlink(filename);
00124     puts("======================================================");
00125     printf("Raw MPI-IO Total Duration %lf seconds, iterations=%u %lf Megabytes written\n",
00126            (endtime-starttime),i,((double)foffset)/(1024.0*1024.0));
00127     printf("Raw MPI-IO Effective Data Rate = %lf Megabytes/sec global and %lf Megabytes/sec per task\n",
00128            (double)(nprocs*localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0),
00129            (double)(localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0));
00130     puts("======================================================");
00131   }
00132 
00133   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00134   /* OK, now we do this using H5Part */
00135   sprintf(newfilename,"testio%u.%u.dat",rank,nprocs);
00136   unlink(newfilename);
00137   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00138   fd = fopen(newfilename,"w");
00139   /* start the timer */
00140   starttime=endtime=MPI_Wtime();
00141   for(i=0;i<NSTEPS;i++){
00142     for(j=0;j<6;j++){
00143       /* touch data */
00144       for(n=0;n<localnp;n++)
00145         (data[j])[n]=(double)rank;
00146       fwrite(data[j],sizeof(double),localnp,fd);
00147     }
00148     curtime=MPI_Wtime(); /* ensure no race condition by broadcasting time */
00149     MPI_Bcast(&curtime,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
00150   }
00151   fclose(fd);
00152   MPI_Barrier(MPI_COMM_WORLD);
00153   endtime=MPI_Wtime();
00154   if(rank==0) puts("*");
00155   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00156   unlink(newfilename);
00157   MPI_Barrier(MPI_COMM_WORLD);
00158   if(rank==0){
00159     puts("======================================================");
00160     printf("Raw 1-file-per-proc Total Duration %lf seconds, iterations=%u %lf Megabytes written\n",
00161            (endtime-starttime),NSTEPS,((double)foffset)/(1024.0*1024.0));
00162     printf("Raw 1-file-per-proc Effective Data Rate = %lf Megabytes/sec global and %lf Megabytes/sec per task\n",
00163            (double)(nprocs*localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0),
00164            (double)(localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0));
00165     puts("======================================================");
00166   }
00167 
00168   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00169   /* OK, now we do this using H5Part */
00170   f = H5PartOpenFileParallel(filename,H5PART_WRITE,MPI_COMM_WORLD);
00171   MPI_Barrier(MPI_COMM_WORLD); /* to prevent unlink from interfering with file open */
00172   /* start the timer */
00173   starttime=endtime=MPI_Wtime();
00174   H5PartSetNumParticles(f,localnp);
00175   for(i=0;i<NSTEPS;i++){
00176     for(j=0;j<6;j++){
00177       /* touch data */
00178       for(n=0;n<localnp;n++)
00179         (data[j])[n]=(double)rank;
00180     }
00181     H5PartSetStep(f,i);
00182     H5PartWriteDataFloat64(f,"x",x);
00183     H5PartWriteDataFloat64(f,"y",y);
00184     H5PartWriteDataFloat64(f,"z",z);
00185     H5PartWriteDataFloat64(f,"px",px);
00186     H5PartWriteDataFloat64(f,"py",py);
00187     H5PartWriteDataFloat64(f,"pz",pz);
00188 
00189     curtime=MPI_Wtime(); /* ensure no race condition by broadcasting time */
00190     MPI_Bcast(&curtime,1,MPI_DOUBLE,0,MPI_COMM_WORLD);
00191   }
00192   H5PartCloseFile(f);
00193   MPI_Barrier(MPI_COMM_WORLD);
00194   endtime=MPI_Wtime();
00195   if(rank==0){
00196     puts("*");
00197     unlink(filename);
00198     puts("======================================================");
00199     printf("H5Part Total Duration %lf seconds, iterations=%u %lf Megabytes written\n",
00200            (endtime-starttime),NSTEPS,((double)foffset)/(1024.0*1024.0));
00201     printf("H5Part Effective Data Rate = %lf Megabytes/sec global and %lf Megabytes/sec per task\n",
00202            (double)(nprocs*localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0),
00203            (double)(localnp*sizeof(double))*((double)NSTEPS)*6.0/((endtime-starttime)*1024.0*1024.0));
00204     puts("======================================================");
00205   }
00206   MPI_Barrier(MPI_COMM_WORLD);
00207   } /* trials */
00208   MPI_Finalize();
00209 }

Generated on Mon Jan 16 13:23:47 2006 for IPPL by  doxygen 1.4.6