src/hdf5/H5ecloud/H5PartTestParallel.cc

Go to the documentation of this file.
00001 #include <stdio.h>
00002 #include <stdlib.h>
00003 #include <string.h>
00004 #include <hdf5.h>
00005 #include "H5Part.hh"
00006 
00007 #ifdef PARALLEL_IO
00008 
00009 /*
00010   This regression test is used to ensure parallel I/O is 
00011   working correctly and that Views are working for 
00012   parallel reads.
00013  */
00014 int main(int argc,char *argv[]){
00015   const int sz=5000;
00016   double *x,*y,*z;
00017   long long *id;
00018   char name[64];
00019   H5PartFile *file;
00020   int i,t,nt,nds;
00021   int nprocs,myproc;
00022   hid_t gid;
00023   MPI_Comm comm=MPI_COMM_WORLD;
00024 
00025   MPI_Init(&argc,&argv);
00026   MPI_Comm_size(comm,&nprocs);
00027   MPI_Comm_rank(comm,&myproc);
00028 
00029   x=(double*)malloc(sz*nprocs*sizeof(double));
00030   y=(double*)malloc(sz*nprocs*sizeof(double));
00031   z=(double*)malloc(sz*nprocs*sizeof(double));
00032   id=(long long*)malloc(sz*nprocs*sizeof(long long));
00033   /* parallel file creation */
00034   file=H5PartOpenFileParallel("parttest.h5",H5PART_WRITE,comm);
00035   if(!file) {
00036     perror("File open failed:  exiting!");
00037     exit(0);
00038   }
00039 
00040   for(t=0;t<5;t++){
00041     MPI_Barrier(comm);
00042     for(i=0;i<sz;i++) {
00043       x[i]=(double)(i+t)+10.0*(double)myproc;
00044       y[i]=0.1 + (double)(i+t);
00045       z[i]=0.2 + (double)(i+t*10);
00046       id[i]=i+sz*myproc;
00047     }
00048     printf("Proc[%u] Writing timestep %u\n",myproc,t);
00049     if(t==0){
00050         printf("Proc[%u]: data values x[first,last]=%f:%f y[%u:%u]=%f:%f z[:]=%f:%f id[:]=%f:%f\n",
00051                 myproc,x[0],x[sz-1],0,sz-1,y[0],y[sz-1],z[0],z[sz-1],(int)id[0],(int)id[sz-1]);
00052     }
00053     H5PartSetStep(file,t); /* must set the current timestep in file */
00054     H5PartSetNumParticles(file,sz); /* then set number of particles to store */
00055     /* now write different tuples of data into this timestep of the file */
00056     H5PartWriteDataFloat64(file,"x",x); 
00057     H5PartWriteDataFloat64(file,"y",y);
00058     H5PartWriteDataFloat64(file,"z",z);
00059 
00060     H5PartWriteDataFloat64(file,"px",x); 
00061     H5PartWriteDataFloat64(file,"py",y);
00062     H5PartWriteDataFloat64(file,"pz",z);
00063 
00064     H5PartWriteDataInt64(file,"id",id);
00065   }
00066 
00067   printf("AllDone p[%u]\n",myproc);
00068   H5PartCloseFile(file);
00069   MPI_Barrier(comm);
00070   
00071   printf("p[%u:%u] : OK, close file and reopen for reading\n",myproc,nprocs);
00072 
00073   file=H5PartOpenFileParallel("parttest.h5",H5PART_READ,comm);
00074   H5PartSetStep(file,0);
00075   unsigned int np,total_np = (int)H5PartGetNumParticles(file);
00076   nt=H5PartGetNumSteps(file); /* get number of steps in file */
00077   nds = H5PartGetNumDatasets(file);
00078   if(myproc==0){
00079     fprintf(stdout,"steps= %u\tdatasets=%u\tparticles= %u\n",
00080             nt,nds,total_np);
00081   }
00082   MPI_Barrier(comm);
00083 
00084   /* now lets compute the appropriate idStart and idEnd 
00085      for this particular processor */
00086 
00087   unsigned long long idStart = sz*myproc;
00088   unsigned long long idEnd   = (sz-1)+sz*myproc;
00089   H5PartSetView(file,idStart,idEnd);
00090   np=H5PartGetNumParticles(file);
00091   printf("Proc[%u]: View=%u:%u : particles= %u\n",
00092          myproc,(int)idStart,(int)idEnd,H5PartGetNumParticles(file));
00093   /* now lets read them and print some out */
00094   H5PartReadDataFloat64(file,"x",x); 
00095   H5PartReadDataFloat64(file,"y",y);
00096   H5PartReadDataFloat64(file,"z",z);
00097   H5PartReadDataInt64(file,"id",id);
00098         printf("Proc[%u]: data values x[first,last]=%f:%f y[%u:%u]=%f:%f z[:]=%f:%f id[:]=%f:%f\n",
00099                 myproc,x[0],x[sz-1],(int)idStart,(int)idEnd,y[0],y[sz-1],z[0],z[sz-1],(int)id[0],(int)id[sz-1]);
00100   
00101   /* H5PartCloseFile(file); MPI_Finalize(); exit(0); */
00102 
00103   if(x) 
00104     free(x); 
00105   if(y) 
00106     free(y);
00107   if(z) 
00108     free(z);
00109   if(id) 
00110     free(id);
00111 
00112   H5PartCloseFile(file);
00113   MPI_Barrier(comm);
00114   fprintf(stderr,"proc[%u]:  done\n",myproc);
00115   return MPI_Finalize();
00116 }
00117 
00118 #else
00119 #error This file only works when PARALLEL_IO is enabled.
00120 #endif

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