src/hdf5/H5ecloud/H5PartTest.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 /*
00008   A simple regression test that shows how you use this API
00009   to write and read multi-timestep files of particle data.
00010 */
00011 #ifdef PARALLEL_IO
00012 
00013 int main(int argc,char *argv[]){
00014   int sz=5;
00015   double *x,*y,*z;
00016   long long *id;
00017   char name[64];
00018   H5PartFile *file;
00019   int i,t,nt,nds;
00020   int nprocs,myproc;
00021   hid_t gid;
00022   MPI_Comm comm=MPI_COMM_WORLD;
00023 
00024   MPI_Init(&argc,&argv);
00025   MPI_Comm_size(comm,&nprocs);
00026   MPI_Comm_rank(comm,&myproc);
00027 
00028   x=(double*)malloc(sz*nprocs*sizeof(double));
00029   y=(double*)malloc(sz*nprocs*sizeof(double));
00030   z=(double*)malloc(sz*nprocs*sizeof(double));
00031   id=(long long*)malloc(sz*nprocs*sizeof(long long));
00032   /* parallel file creation */
00033   file=H5PartOpenFileParallel("parttest.h5",H5PART_WRITE,comm);
00034   if(!file) {
00035     perror("File open failed:  exiting!");
00036     exit(0);
00037   }
00038 
00039   for(t=0;t<5;t++){
00040     MPI_Barrier(comm);
00041     for(i=0;i<sz;i++) {
00042       x[i]=(double)(i+t)+10.0*(double)myproc;
00043       y[i]=0.1 + (double)(i+t);
00044       z[i]=0.2 + (double)(i+t*10);
00045       id[i]=i+sz*myproc;
00046     }
00047     printf("Proc[%u] Writing timestep %u file=%u\n",myproc,t,file->file);
00048     H5PartSetStep(file,t); /* must set the current timestep in file */
00049     H5PartSetNumParticles(file,sz); /* then set number of particles to store */
00050     /* now write different tuples of data into this timestep of the file */
00051     H5PartWriteDataFloat64(file,"x",x); 
00052     H5PartWriteDataFloat64(file,"y",y);
00053     H5PartWriteDataFloat64(file,"z",z);
00054 
00055     H5PartWriteDataFloat64(file,"px",x); 
00056     H5PartWriteDataFloat64(file,"py",y);
00057     H5PartWriteDataFloat64(file,"pz",z);
00058 
00059     H5PartWriteDataInt64(file,"id",id);
00060   }
00061 
00062   unsigned int idStart = 0+sz*myproc;
00063   unsigned int idEnd = (sz-1)+sz*myproc;
00064 
00065   printf("AllDone p[%u]\n",myproc);
00066   H5PartCloseFile(file);
00067   fprintf(stderr,"Closed files p[%u]\n",myproc);
00068   MPI_Barrier(comm);
00069   
00070   fprintf(stderr,"p[%u:%u] : OK, close file and reopen for reading idStart %u  idEnd %u \n",myproc,nprocs,idStart,idEnd);
00071 
00072   file=H5PartOpenFileParallel("parttest.h5",H5PART_READ,comm);
00073   H5PartSetStep(file,0);
00074  // unsigned int np = 0;
00075   unsigned int np = (int)H5PartGetNumParticles(file);
00076   nt=H5PartGetNumSteps(file); /* get number of steps in file */
00077   nds=H5PartGetNumDatasets(file); /* get number of datasets in timestep 0 */
00078 
00079   MPI_Barrier(comm);
00080   
00081   H5PartSetView(file,idStart,idEnd);
00082 
00083   np = (int)H5PartGetNumParticles(file);
00084   printf("After SetView(%d,%d): steps= %u  datasets= %u particles= %u\n",
00085         (int)idStart,(int)idEnd,
00086         nt,nds,np);
00087 
00088   if(x) 
00089     free(x); 
00090   if(y) 
00091     free(y);
00092   if(z) 
00093     free(z);
00094   if(id) 
00095     free(id);
00096 
00097   H5PartCloseFile(file);
00098   MPI_Barrier(comm);
00099   fprintf(stderr,"proc[%u]:  done\n",myproc);
00100   return MPI_Finalize();
00101 }
00102 
00103 #else
00104 int main(int argc,char *argv[]){
00105   int sz=10;
00106   double *x,*y,*z;
00107   long long *id;
00108   char name[64];
00109   H5PartFile *file;
00110   int i,t,nt,nds,np;
00111   hid_t gid;
00112   long long idStart = 0;
00113   long long idEnd = 0;
00114 
00115 
00116   x=(double*)malloc(sz*sizeof(double));
00117   y=(double*)malloc(sz*sizeof(double));
00118   z=(double*)malloc(sz*sizeof(double));
00119   id=(long long*)malloc(sz*sizeof(long long));
00120   /* parallel file creation */
00121   file=H5PartOpenFile("parttest.h5",H5PART_WRITE);
00122   if(!file) {
00123     perror("File open failed:  exiting!");
00124     exit(0);
00125   }
00126 
00127   for(t=0;t<5;t++){
00128     fprintf(stdout,"Writing timestep %u\n",t);
00129     for(i=0;i<sz;i++) {
00130       x[i]=(double)(i+t);
00131       y[i]=0.1 + (double)(i+t);
00132       z[i]=0.2 + (double)(i+t*10);
00133       id[i]=i;
00134       fprintf(stdout,"\tp[%u] x=%f y=%f z=%f id=%d\n",
00135               i,x[i],y[i],z[i],(int)id[i]);
00136     }
00137     H5PartSetStep(file,t); /* must set the current timestep in file */
00138     H5PartSetNumParticles(file,sz); /* then set number of particles to store */
00139     /* now write different tuples of data into this timestep of the file */
00140     H5PartWriteDataFloat64(file,"x",x); 
00141     H5PartWriteDataFloat64(file,"y",y);
00142     H5PartWriteDataFloat64(file,"z",z);
00143 
00144     H5PartWriteDataFloat64(file,"px",x); 
00145     H5PartWriteDataFloat64(file,"py",y);
00146     H5PartWriteDataFloat64(file,"pz",z);
00147 
00148     H5PartWriteDataInt64(file,"id",id);
00149   }
00150 
00151   printf("AllDone writing\n");
00152   H5PartCloseFile(file);
00153   
00154  
00155   /*+++++++++++++ Reopen File for Reading +++++++++++*/
00156   file=H5PartOpenFile("parttest.h5",H5PART_READ);
00157 
00158   
00159   /********************************************/
00160   H5PartSetStep(file,0);
00161   nt=H5PartGetNumSteps(file); /* get number of steps in file */
00162   nds=H5PartGetNumDatasets(file); /* get number of datasets in timestep 0 */
00163   np=H5PartGetNumParticles(file);
00164   
00165 
00166   fprintf(stdout,"OK, close file and reopen for reading\n");
00167   fprintf(stdout,"steps= %u\tdatasets=%u\tparticles= %u\n",
00168          nt,nds,np);
00169   
00170   // clear the particles
00171   for(i=0;i<np;i++){
00172     x[i]=y[i]=z[i]=0.0;
00173     id[i]=0;
00174   }
00175 
00176   H5PartReadDataFloat64(file,"x",x);
00177   H5PartReadDataFloat64(file,"y",y);
00178   H5PartReadDataFloat64(file,"z",z);
00179   H5PartReadDataInt64(file,"id",id);
00180 
00181   for(i=0;i<np;i++){
00182     fprintf(stdout,"\tp[%3u] x=%lf y=%lf z=%lf id=%lu\n",
00183             i,x[i],y[i],z[i],(unsigned)(id[i]));
00184   }
00185   /********************************************/
00186   printf("Set to last step and reload data\n");
00187   H5PartSetStep(file,nt-1);
00188   H5PartReadDataFloat64(file,"x",x);
00189   H5PartReadDataFloat64(file,"y",y);
00190   H5PartReadDataFloat64(file,"z",z);
00191   H5PartReadDataInt64(file,"id",id);
00192   for(i=0;i<np;i++){
00193     fprintf(stdout,"\tp[%3u] x=%lf y=%lf z=%lf id=%lu\n",
00194             i,x[i],y[i],z[i],(unsigned)(id[i]));
00195   }
00196   
00197   /********************************************/
00198   idEnd=np;
00199   printf("Old View is %d:%d\n",(int)idStart,(int)idEnd);
00200   H5PartSetView(file,idStart,idEnd>>1);
00201   printf("Set new view = %d:%d\n",(int)idStart,(int)(idEnd>>1));
00202   H5PartGetView(file,&idStart,&idEnd);
00203   np=H5PartGetNumParticles(file);
00204   printf("steps= %u  datasets= %u particles= %d with view %d:%d\n",
00205          nt,nds,(int)np,(int)idStart,(int)idEnd);
00206   H5PartSetStep(file,nt-1); // set to last step
00207   printf("Setting to last step = %u\n",nt-1);
00208   for(i=0;i<10;i++){ x[i]=y[i]=z[i]=0.0; id[i]=0; } /* clear the arrays */
00209   H5PartReadDataFloat64(file,"x",x);
00210   H5PartReadDataFloat64(file,"y",y);
00211   H5PartReadDataFloat64(file,"z",z);
00212   H5PartReadDataInt64(file,"id",id);
00213 
00214   for(i=0;i<np;i++){
00215     fprintf(stdout,"\tp[%3u] x=%lf y=%lf z=%lf id=%lu\n",i,x[i],y[i],z[i],(unsigned)id[i]);
00216   }
00217 
00218   /********************************************/
00219   printf("Now set the view to the latter half of the data in step #%u\n",nt-1);
00220   H5PartResetView(file);
00221   H5PartGetView(file,&idStart,&idEnd);
00222   printf("Reset view = %d:%d\nSetting to %u:%u\n",
00223          (int)idStart,(int)idEnd,
00224          (int)idEnd>>1,(int)idEnd);
00225   H5PartSetView(file,(idEnd>>1),idEnd);
00226   np=H5PartGetNumParticles(file);
00227   printf("Now particles in selection are %d\n",np);
00228   printf("doubleCheck=%d\n",H5PartGetView(file,0,0));
00229   
00230   for(i=0;i<10;i++){ x[i]=y[i]=z[i]=0.0; id[i]=0; } /* clear the arrays */
00231  
00232   H5PartReadDataFloat64(file,"x",x);
00233   H5PartReadDataFloat64(file,"y",y);
00234   H5PartReadDataFloat64(file,"z",z);
00235   H5PartReadDataInt64(file,"id",id);
00236   for(i=0;i<np;i++){
00237     fprintf(stdout,"\tp[%3u] x=%lf y=%lf z=%lf id=%lu\n",i,x[i],y[i],z[i],(unsigned)id[i]);
00238   }
00239   if(x) 
00240     free(x); 
00241   if(y) 
00242     free(y);
00243   if(z) 
00244     free(z);
00245   if(id) 
00246     free(id);
00247 
00248   H5PartCloseFile(file);
00249   fprintf(stderr,"done\n");
00250 }
00251 
00252 #endif
00253 
00254 
00255 

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