src/hdf5/H5ecloud/H5PartAndreasTest.cc

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

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