src/hdf5/parallel/H5Part.c

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.h"
00006 
00007 static unsigned h5part_debug=0;
00008 
00009 /********* Private Function Declarations *************/
00010 herr_t H5PartIOcounter(hid_t group_id,
00011                    const char *member_name,
00012                    void *operator_data);
00013 
00014 /********** Definitions of Functions/Subroutines ******/
00015 
00016 /*========== File Opening/Closing ===============*/
00017 /*
00018   H5PartOpenFile:  Opens file with specified filename. 
00019     If you open with flag H5PART_WRITE, it will truncate any
00020     file with the specified filename and start writing to it.
00021     If you open with H5PART_READ, then it will open the file
00022     readonly.  I could probably do an APPEND mode as well, but
00023     just haven't done so yet.  APPEND would need to check the
00024     file to determine its current state and modify the 
00025     H5PartFile datastructure accordingly.
00026 
00027     H5PartFile should be treated as an essentially opaque
00028     datastructure.  It acts as the file handle, but internally
00029     it maintains several key state variables associated with 
00030     the file.
00031 */
00032 
00033 H5PartFile *H5PartOpenFileParallel(const char *filename,unsigned flags
00034 #ifdef PARALLEL_IO
00035                 ,                  MPI_Comm comm
00036 #endif
00037 ){
00038   H5PartFile *f=(H5PartFile*)malloc(sizeof(H5PartFile));
00039 #ifdef PARALLEL_IO
00040   MPI_Info info = MPI_INFO_NULL; /* for the SP2... perhaps different for linux */
00041 #endif
00042   f->xfer_prop = f->create_prop = f->access_prop=H5P_DEFAULT;
00043 #ifdef PARALLEL_IO
00044   MPI_Comm_size(comm,&f->nprocs);
00045   MPI_Comm_rank(comm,&f->myproc);
00046   f->pnparticles=(long long*)malloc(sizeof(long long)*f->nprocs);
00047   /* access_prop: properties like chunking or parallel I/O */
00048   f->access_prop = H5Pcreate(H5P_FILE_ACCESS);
00049   if(H5Pset_fapl_mpio(f->access_prop,comm,info)<0){
00050         fprintf(stderr,"Total failure trying to setup mpi-io for access\n");
00051         exit(0);
00052   }
00053   /* create_prop: tunable parameters like blocksize and btree sizes */
00054   /* f->create_prop = H5Pcreate(H5P_FILE_CREATE); */
00055   f->create_prop = H5P_DEFAULT;
00056   /* currently create_prop is empty */
00057   /* xfer_prop:  also used for parallel I/O, during actual writes
00058           rather than the access_prop which is for file creation. */
00059   f->xfer_prop = H5Pcreate(H5P_DATASET_XFER);
00060   H5Pset_dxpl_mpio(f->xfer_prop,H5FD_MPIO_COLLECTIVE);
00061   f->comm = comm;
00062 #endif
00063   if(flags==H5PART_READ){
00064     f->file=H5Fopen(filename,H5F_ACC_RDONLY,f->access_prop);
00065     /* just do this serially 
00066        f->file = H5Fopen(filename,H5F_ACC_RDONLY,H5P_DEFAULT); */
00067   }
00068   else if(flags==H5PART_WRITE){
00069     f->file=H5Fcreate(filename,H5F_ACC_TRUNC,f->create_prop,f->access_prop);
00070   }
00071   else {
00072     f->file=-1;
00073     fprintf(stderr,"H5Open: Invalid file access type %d\n",flags);
00074   }
00075   if(f->file<=0) {
00076     free(f);
00077     fprintf(stderr,"H5Open: File open failed for file [%s]\n",
00078             filename);
00079         exit(0); 
00080         return 0;
00081   }
00082 #ifdef PARALLEL_IO
00083   else {
00084         if(h5part_debug) fprintf(stderr,"Proc[%u] Opened file %s val=%d\n",
00085                 f->myproc,filename,f->file);
00086   }
00087 #endif
00088   f->mode=flags;
00089   f->maxstep=0;
00090   f->timegroup=0;
00091   f->timestep=0;
00092   f->shape=0;
00093   f->diskshape=H5S_ALL;
00094   f->memshape=H5S_ALL;
00095   f->viewstart=-1;
00096   f->viewend=-1;
00097   return f;
00098 }
00099 
00100 H5PartFile *H5PartOpenFile(const char *filename,unsigned flags){
00101   H5PartFile *f=(H5PartFile*)malloc(sizeof(H5PartFile));
00102   /* printf("filename is [%s] flags=%0X\n",filename,flags); */
00103   
00104   f->xfer_prop = f->create_prop = f->access_prop=H5P_DEFAULT;
00105 #ifdef PARALLEL_IO
00106   f->pnparticles=0;
00107   f->comm = MPI_COMM_WORLD;
00108   f->nprocs=1;
00109   f->myproc=0;
00110 #endif
00111   if(flags==H5PART_READ){
00112     /*  f->file=H5Fopen(filename,H5F_ACC_RDONLY,f->access_prop); */
00113     /* just do this serially */
00114     f->file = H5Fopen(filename,H5F_ACC_RDONLY,H5P_DEFAULT);
00115     if(f->file==-1){
00116       fprintf(stderr,"H5Open: File %s not found\n",filename);
00117     }
00118   }
00119   else if(flags==H5PART_WRITE){
00120     f->file=H5Fcreate(filename,H5F_ACC_TRUNC,H5P_DEFAULT,H5P_DEFAULT);
00121   }
00122   else {
00123     f->file=-1;
00124     fprintf(stderr,"H5Open: Invalid file access type %d\n",flags);
00125   }
00126   if(f->file<0) {
00127     free(f);
00128     fprintf(stderr,"H5Open: File open failed for file [%s]\n",
00129             filename);
00130     return 0;
00131   }
00132   f->mode=flags;
00133   f->maxstep=0;
00134   f->timegroup=0;
00135   f->timestep=0;
00136   f->shape=0;
00137   f->diskshape=H5S_ALL;
00138   f->memshape=H5S_ALL;
00139   f->viewstart=-1;
00140   f->viewend=-1;
00141   return f;
00142 }
00143 int H5PartFileIsValid(H5PartFile *f){
00144   if(f->file > 0) return 1;
00145   else return 0;
00146 }
00147 /*
00148   H5PartCloseFile:  Does what it says it does.
00149 */
00150 void H5PartCloseFile(H5PartFile *f){
00151   if(f->shape>0) H5Sclose(f->shape); f->shape=0;
00152   if(f->timegroup>0) H5Gclose(f->timegroup); f->timegroup=0;
00153   if(f->diskshape!=H5S_ALL) H5Sclose(f->diskshape);
00154 #ifdef PARALLEL_IO
00155   if(f->pnparticles) free(f->pnparticles);
00156 #endif
00157   if(f->xfer_prop!=H5P_DEFAULT) {
00158     H5Pclose(f->xfer_prop); f->xfer_prop=H5P_DEFAULT;
00159   }
00160   if(f->access_prop!=H5P_DEFAULT) {
00161     H5Pclose(f->access_prop); f->access_prop=H5P_DEFAULT;
00162   }  
00163   if(f->create_prop!=H5P_DEFAULT) {
00164     H5Pclose(f->create_prop); f->create_prop=H5P_DEFAULT;
00165   }
00166   if(f->file) H5Fclose(f->file); f->file=0;
00167   free(f);
00168 }
00169 
00170 
00171 /*============== File Writing Functions ==================== */
00172 /*
00173   H5PartSetNumParticles:  This is a private function used by the 
00174   other functions that write arrays.  Its sole purpose is to 
00175   prevent needless creation of new HDF5 DataSpace handles if
00176   the number of particles is invariant throughout the sim.
00177   That's its only reason for existence.
00178 */
00179 void H5PartSetNumParticles(H5PartFile *f,long long nparticles){
00180   hssize_t start[1];
00181   hsize_t stride[1],count[1],total;
00182   register int i,r;
00183   /*  printf("Set Number of Particles to %lld\n",nparticles); */
00184 #ifndef PARALLEL_IO
00185   /* if we are not using parallel-IO, there is enough information
00186      to know that we can short circuit this routine.  However,
00187      for parallel IO, this is going to cause problems because
00188      we don't know if things have changed globally */
00189   if(f->nparticles==nparticles) return;  /* no change  */
00190 #endif
00191   if(f->diskshape!=H5S_ALL) H5Sclose(f->diskshape);
00192   if(f->memshape!=H5S_ALL) H5Sclose(f->memshape);
00193   f->memshape=f->diskshape = H5S_ALL;
00194   if(f->shape) H5Sclose(f->shape);
00195   f->nparticles=(hsize_t)nparticles;
00196 #ifndef PARALLEL_IO
00197   f->shape=H5Screate_simple(1, /* 1 dimensional shape */
00198                             &(f->nparticles), /* with nparticles elements */
00199                             NULL); /* ignore this :-) */
00200 #else /* PARALLEL_IO */
00201   /* The Gameplan here is to declare the overall size of the on-disk
00202      data structure the same way we do for the serial case.  But
00203      then we must have additional "DataSpace" structures to define
00204      our in-memory layout of our domain-decomposed portion of the particle
00205      list as well as a "selection" of a subset of the on-disk 
00206      data layout that will be written in parallel to mutually exclusive
00207      regions by all of the processors during a parallel I/O operation.
00208      These are f->shape, f->memshape and f->diskshape respectively. */
00209   /* acquire the number of particles to be written from each MPI process */
00210   MPI_Allgather(&nparticles,1,MPI_LONG_LONG,f->pnparticles,1,MPI_LONG_LONG,f->comm);
00211   if(h5part_debug) {
00212     if(f->myproc==0){
00213       puts("AllGather:  Particle offsets:\n");
00214       for(i=0;i<f->nprocs;i++) 
00215         printf("\tnp=%u\n",f->pnparticles[i]); /* compute total nparticles */
00216     }
00217   }
00218   /* should I create a selection here? */
00219   stride[0]=1;
00220   start[0]=0;
00221   for(i=0;i<f->myproc;i++) start[0]+=f->pnparticles[i]; /* compute start offsets */
00222   total=0;
00223   for(i=0;i<f->nprocs;i++) total+=f->pnparticles[i]; /* compute total nparticles */
00224   f->shape = H5Screate_simple(1,&total,NULL); /* declare overall datasize */
00225   f->diskshape = H5Screate_simple(1,&total,NULL); /* declare overall data size  but then will select a subset */
00226   f->memshape = H5Screate_simple(1,&(f->nparticles),NULL); /* declare local memory datasize */
00227 
00228   if(f->shape<0 || f->memshape<0 || f->diskshape<0){
00229         fprintf(stderr,"Abort:  shape construction failed\n");
00230         exit(0);
00231 }
00232   count[0]=nparticles; /* based on local nparticles (for the selection */
00233   /* and then set the subset of the data you will write to */
00234   r=H5Sselect_hyperslab(f->diskshape,H5S_SELECT_SET,start,stride,count,NULL);
00235   if(r<0){
00236         fprintf(stderr,"Abort: Selection Failed!\n");
00237         exit(0);
00238   }
00239   if(f->timegroup<0){
00240     H5PartSetStep(f,0);
00241   }
00242 #endif
00243 }
00244 
00245 /*
00246   H5PartWriteArray:  This writes an individual array of data
00247     for a given timestep.  It probably would be convenient if this
00248     were public, but its private because it there are some order/
00249     state dependencies that get kind of nasty to manage (and therefore
00250     are not managed).  It is otherwise assuming that the timegroup
00251     has already been created and other such details.
00252     It is called by the mongo H5PartWriteStep() routine.
00253  */
00254 int H5PartWriteDataFloat64(H5PartFile *f,char *name,double *array){
00255   register int r;
00256   hid_t dataset;
00257  if(f->mode==H5PART_READ){
00258         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00259         return 0;
00260   }
00261   /* fprintf(stderr,"Create a dataset[%s]  mounted on the timegroup %u\n",
00262          name,f->timestep); */
00263   if(f->timegroup<=0) {fprintf(stderr,"Eeeerrrrroooorrrr!!!!\n");}
00264   dataset = H5Dcreate(f->timegroup,name,H5T_NATIVE_DOUBLE,f->shape,H5P_DEFAULT);
00265   if(dataset<0)
00266     fprintf(stderr,"Dataset open failed for name=[%s] step=%u!\n",
00267             name,f->timestep);
00268   r=H5Dwrite(dataset,H5T_NATIVE_DOUBLE,f->memshape,f->diskshape,H5P_DEFAULT,array);
00269   H5Dclose(dataset);
00270   return r;
00271 }
00272 int H5PartWriteDataInt64(H5PartFile *f,char *name,long long *array){
00273   register int r;
00274   hid_t dataset;
00275   /*fprintf(stderr,"Create a dataset[%s]  mounted on the timegroup %u\n",
00276          name,f->timestep); */
00277  if(f->mode==H5PART_READ){
00278         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00279         return 0;
00280   }
00281   if(f->timegroup<=0) {fprintf(stderr,"Eeeerrrrroooorrrr!!!!\n");}
00282   dataset = H5Dcreate(f->timegroup,name,H5T_NATIVE_INT64,f->shape,H5P_DEFAULT);
00283   if(dataset<0){
00284     fprintf(stderr,"Dataset open failed for name=[%s] step=%u!\n",
00285             name,f->timestep);
00286     exit(0);
00287   }
00288   r=H5Dwrite(dataset,H5T_NATIVE_INT64,f->memshape,f->diskshape,H5P_DEFAULT,array);
00289   if(r<0) {
00290         fprintf(stderr,"Attempt to write dataset failed!\n");
00291         exit(0);
00292   }
00293   H5Dclose(dataset);
00294   return r;
00295 }
00296 int H5PartWriteFileAttribString(H5PartFile *f,char *name,
00297                                 char *attrib){
00298  if(f->mode==H5PART_READ){
00299         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00300         return 0;
00301   }
00302   return H5PartWriteFileAttrib(f,name,H5T_NATIVE_CHAR,attrib,strlen(attrib)+1);
00303 }
00304 int H5PartWriteStepAttribString(H5PartFile *f,char *name,
00305                                 char *attrib){
00306  if(f->mode==H5PART_READ){
00307         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00308         return 0;
00309   }
00310   return H5PartWriteStepAttrib(f,name,H5T_NATIVE_CHAR,attrib,strlen(attrib)+1);
00311 }
00312 int H5PartWriteStepAttrib(H5PartFile *f,char *name,
00313                           hid_t type,void *value,int nelem){
00314   register int r;
00315   hid_t attrib;
00316   hid_t space;
00317   hsize_t len;
00318  if(f->mode==H5PART_READ){
00319         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00320         return 0;
00321   }
00322   /*fprintf(stderr,"Create a attribute[%s]  mounted on the timegroup %u\n",
00323          name,f->timestep); */
00324   if(f->timegroup<=0) {fprintf(stderr,"Eeeerrrrroooorrrr!!!!\n");}
00325   len = nelem;
00326   space = H5Screate_simple(1,&len,NULL);
00327   attrib = H5Acreate(f->timegroup,name,type,space,H5P_DEFAULT);
00328   r=H5Awrite(attrib,type,value);
00329   H5Aclose(attrib);
00330   H5Sclose(space);
00331   return r;
00332 }
00333 
00334 int H5PartWriteAttrib(H5PartFile *f,char *name,
00335     hid_t type,void *value,int nelem){
00336  if(f->mode==H5PART_READ){
00337         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00338         return 0;
00339   }
00340     return H5PartWriteStepAttrib(f,name,type,value,nelem);
00341 }
00342 
00343 
00344 int H5PartWriteFileAttrib(H5PartFile *f,char *name,
00345                           hid_t type,void *value,int nelem){
00346   register int r;
00347   hid_t attrib;
00348   hid_t space;
00349   hsize_t len;
00350  if(f->mode==H5PART_READ){
00351         fprintf(stderr,"H5PartWriteDataFloat64:  Error!  Attempting to write to read-only file\n");
00352         return 0;
00353   }
00354   /*fprintf(stderr,"Create a attribute[%s]  mounted on the timegroup %u\n",
00355          name,f->timestep); */
00356   if(f->file<=0) {fprintf(stderr,"Eeeerrrrroooorrrr!!!! File is not open!\n");}
00357   len = nelem;
00358   space = H5Screate_simple(1,&len,NULL);
00359   attrib = H5Acreate(f->file,name,type,space,H5P_DEFAULT);
00360   r=H5Awrite(attrib,type,value);
00361   H5Aclose(attrib);
00362   H5Sclose(space);
00363   return r;
00364 }
00365 herr_t H5PartAttribcounter(hid_t group_id,
00366                    const char *member_name,
00367                    void *operator_data){
00368   int *count = (int*)operator_data;
00369   (*count)++;
00370   return 0;
00371 }
00372 int H5PartGetNumStepAttribs(H5PartFile *f){ /* for current filestep */
00373   return H5Aget_num_attrs(f->timegroup);
00374 }
00375 int H5PartGetNumFileAttribs(H5PartFile *f){
00376   /* must walk the attributes */
00377   /* iterate to get numsteps */
00378   int nattribs;
00379   hid_t rootgroup=H5Gopen(f->file,"/");;
00380   nattribs=H5Aget_num_attrs(rootgroup); /* open / for the file?
00381                                        or is there a problem with 
00382                                        file attributes? */
00383   H5Gclose(rootgroup);
00384   return nattribs;
00385 }
00386 
00387 hid_t H5PartNormType(hid_t type){
00388   H5T_class_t tclass = H5Tget_class(type);
00389   int size = H5Tget_size(type);
00390   switch(tclass){
00391   case H5T_INTEGER:
00392     if(size==8) return H5T_NATIVE_INT64;
00393     else if(size==1) return H5T_NATIVE_CHAR;
00394     else fprintf(stderr,"Unknown type %u.  Cannot infer normalized type\n",type);
00395     break;
00396   case H5T_FLOAT:
00397     return H5T_NATIVE_DOUBLE;
00398     break;
00399   }
00400   return -1;
00401 }
00402 
00403 void H5PartGetStepAttribInfo(H5PartFile *f,int idx,
00404                          char *name,size_t maxname,hid_t *type,int *nelem){
00405   hid_t attrib=H5Aopen_idx(f->timegroup,idx);
00406   hid_t mytype=H5Aget_type(attrib);
00407   hid_t space = H5Aget_space(attrib);
00408   if(nelem)
00409     *nelem=H5Sget_simple_extent_npoints(space);
00410   if(name)
00411     H5Aget_name(attrib,maxname,name);
00412   if(type)
00413     *type=H5PartNormType(mytype); /* normalize the type to be machine-native
00414                                  this means one of 
00415                                 H5T_NATIVE_DOUBLE
00416                                 H5T_NATIVE_INT64
00417                                 H5T_NATIVE_CHAR */
00418   H5Sclose(space);
00419   H5Tclose(mytype);
00420   H5Aclose(attrib);
00421 }
00422 void H5PartGetFileAttribInfo(H5PartFile *f,int idx,
00423                          char *name,size_t maxname,hid_t *type,int *nelem){
00424   hid_t attrib=H5Aopen_idx(f->file,idx);
00425   hid_t mytype=H5Aget_type(attrib);
00426   hid_t space = H5Aget_space(attrib);
00427   if(nelem)
00428     *nelem=H5Sget_simple_extent_npoints(space);
00429   if(name)
00430     H5Aget_name(attrib,maxname,name);
00431   if(type)
00432     *type=H5PartNormType(mytype); /* normalize the type to be machine-native
00433                                  this means one of 
00434                                 H5T_NATIVE_DOUBLE
00435                                 H5T_NATIVE_INT64
00436                                 H5T_NATIVE_CHAR */
00437   H5Sclose(space);
00438   H5Tclose(mytype);
00439   H5Aclose(attrib);
00440 }
00441 
00442 int H5PartReadStepAttrib(H5PartFile *f,char *name,void *data){
00443   /* use the open attribute by name mode of operation */
00444   hid_t attrib=H5Aopen_name(f->timegroup,name);
00445   if(attrib<=0) return -1;
00446   hid_t mytype=H5Aget_type(attrib);
00447   hid_t space = H5Aget_space(attrib);
00448   hsize_t nelem=H5Sget_simple_extent_npoints(space);
00449   hid_t type=H5PartNormType(mytype); /* normalize the type to be machine-native
00450                                  this means one of 
00451                                 H5T_NATIVE_DOUBLE
00452                                 H5T_NATIVE_INT64
00453                                 H5T_NATIVE_CHAR */
00454   H5Aread(attrib,type,data);
00455   H5Sclose(space);
00456   H5Tclose(mytype);
00457   H5Aclose(attrib);
00458   return 1;
00459 }
00460 
00461 void H5PartReadAttrib(H5PartFile *f,char *name,void *data){
00462     /* use the open attribute by name mode of operation */
00463     H5PartReadStepAttrib(f,name,data);
00464 }
00465 
00466 
00467 
00468 
00469 int H5PartReadFileAttrib(H5PartFile *f,char *name,void *data){
00470   hid_t attrib=H5Aopen_name(f->file,name);
00471   if(attrib<0) return -1;
00472   hid_t mytype=H5Aget_type(attrib);
00473   hid_t space = H5Aget_space(attrib);
00474   hsize_t nelem=H5Sget_simple_extent_npoints(space);
00475   hid_t type=H5PartNormType(mytype); /* normalize the type to be machine-native
00476                                  this means one of 
00477                                 H5T_NATIVE_DOUBLE
00478                                 H5T_NATIVE_INT64
00479                                 H5T_NATIVE_CHAR */
00480   H5Aread(attrib,type,data);
00481   H5Sclose(space);
00482   H5Tclose(mytype);
00483   H5Aclose(attrib);
00484   return 1;
00485 }
00486 
00487 
00488 /*================== File Reading Routines =================*/
00489 /*
00490   H5PartSetStep:  This routine is *only* useful when you are reading
00491   data.  Using it while you are writing will generate nonsense results!
00492   (This file format is only half-baked... robustness is not std equipment!)
00493   So you use this to random-access the file for a particular timestep.
00494   Failure to explicitly set the timestep on each read will leave you
00495   stuck on the same timestep for *all* of your reads.  That is to say
00496   the writes auto-advance the file pointer, but the reads do not
00497   (they require explicit advancing by selecting a particular timestep).
00498 */
00499 void H5PartSetStep(H5PartFile *f,int step){
00500   char name[128];
00501 #ifdef PARALLEL_IO
00502   if(h5part_debug) printf("Proc[%u] SetStep to %u for file %u\n",f->myproc,step,f->file);
00503 #else
00504   if(h5part_debug) fprintf(stderr,"SetStep to %d for file\n",step);
00505 #endif
00506   f->timestep=step;   /* because we start at 0 */
00507   if(f->timegroup>0) {
00508         H5Gclose(f->timegroup); 
00509   }
00510   f->timegroup=-1;
00511   sprintf(name,"Particles#%u",f->timestep);
00512   /*if(h5part_debug) fprintf(stderr,"P[%u] Open timegroup [%s] for file %d\n",
00513         f->myproc,name,f->file);*/
00514   if(f->mode==H5PART_READ) { /* kludge to prevent error messages */
00515 /*      if(h5part_debug) fprintf(stderr,"P[%u] open group\n"); */
00516     f->timegroup=H5Gopen(f->file,name);
00517   }
00518   if(f->timegroup<=0){
00519     /* timegroup doesn't exist, so we need to create one */
00520         /* if(h5part_debug) fprintf(stderr,"P[%u] create group\n"); */
00521     if(f->mode == H5PART_READ){
00522       fprintf(stderr,"Error in H5PartSetStep() Step #%d does not exist!\n",step);
00523     }
00524     f->timegroup = H5Gcreate(f->file,name,0);
00525     if(f->timegroup<0) {
00526 #ifdef PARALLEL_IO
00527         fprintf(stderr,"p[%u] timegroup creation failed!\n",
00528                 f->myproc);
00529 #endif
00530         return;
00531     }
00532   }
00533   /*
00534     #ifdef PARALLEL_IO
00535     H5PartWriteStepAttrib(f,"pnparticles",H5T_NATIVE_INT64,f->pnparticles,f->nprocs);
00536     #endif
00537   */
00538 }
00539 
00540 /*
00541   H5PartIOcounter:  This is an entirely internal callback function 
00542    which is used in conjunction with HDF5 iterators.  The HDF5 Group
00543    iterator will call this repeatedly in order to count how many
00544    timesteps of data have been stored in a particular file.
00545    This is used by H5PartGetNumSteps().
00546 */
00547 herr_t H5PartIOcounter(hid_t group_id,
00548                    const char *member_name,
00549                    void *operator_data){
00550   int *count = (int*)operator_data;
00551   H5G_stat_t objinfo;
00552   /* only count the particle groups... ignore all others */
00553   if(!strncmp(member_name,"Particles",9)) (*count)++;
00554   return 0;
00555 }
00556 herr_t H5PartDScounter(hid_t group_id,
00557                    const char *member_name,
00558                    void *operator_data){
00559   int *count = (int*)operator_data;
00560   H5G_stat_t objinfo;
00561   /* only count the datasets... ignore all others */
00562   if(H5Gget_objinfo(group_id, 
00563                  member_name,
00564                  1 /* follow links */, 
00565                     &objinfo)<0) {
00566     return 0; /* error (probably bad symlink) */
00567   }
00568   if(objinfo.type==H5G_DATASET)    (*count)++;
00569   return 0;
00570 }
00571 
00572 typedef struct H5IO_getname_t {
00573   int index,count;
00574   char *name;
00575 }H5IO_getname_t;
00576 
00577 herr_t H5IOgetname(hid_t group_id,
00578                           const char *member_name, 
00579                           void *operator_data){
00580   H5IO_getname_t *getn = (H5IO_getname_t*)operator_data;
00581   /* check type first (only respond if it is a dataset) */
00582   H5G_stat_t objinfo;
00583   /* request info about the type of objects in root group */
00584   if(H5Gget_objinfo(group_id, 
00585                     member_name,
00586                     1 /* follow links */, 
00587                     &objinfo)<0) return 0; /* error (probably bad symlink) */
00588   /* only count objects that are datasets (not subgroups) */
00589   if(objinfo.type!=H5G_DATASET) 
00590     return 0; /* do not increment count if it isn't a dataset. */
00591   if(getn->index==getn->count){
00592     strcpy(getn->name,member_name);
00593     return 1; /* success */
00594   }
00595   getn->count++;
00596   return 0;
00597 }
00598 
00599 /* 
00600    H5PartGetNumSteps:  This reads the number of datasteps that are 
00601      currently stored in the datafile.  (simple return of an int).
00602      It works for both reading and writing of files, but is probably
00603      only typically used when you are reading.
00604 */
00605 int H5PartGetNumSteps(H5PartFile *f){
00606   /* iterate to get numsteps */
00607   int count=0;
00608   int idx=0;
00609   while(H5Giterate(f->file, /* hid_t loc_id, */
00610                    "/", /*const char *name, */
00611                    &idx, /* int *idx, */
00612                    H5PartIOcounter,
00613                    &count)<0){}  
00614   return count;
00615 }
00616 int H5PartGetNumDatasets(H5PartFile *f){
00617   int count=0;
00618   int idx=0;
00619   char name[128];
00620   sprintf(name,"Particles#%u",f->timestep);
00621   while(H5Giterate(f->file, /* hid_t loc_id */
00622                    name,
00623                    &idx,
00624                    H5PartDScounter,
00625                    &count)<0){}
00626   return count;   
00627 }
00628 int H5PartGetDatasetName(H5PartFile *f,
00629                          int _index, /* index of the dataset (0 to N-1) */
00630                          char *name, /* buffer to store name of dataset */
00631                          int maxlen){ /* max size of the "name" buffer */
00632   H5IO_getname_t getn;
00633   int idx=_index;
00634   char gname[128];
00635   sprintf(gname,"Particles#%u",f->timestep);
00636   getn.index=_index; getn.name=name; getn.count=_index;
00637   while(H5Giterate(f->file, /* hid_t loc_id, */
00638                    gname, /*const char *name, */
00639                    &idx, /* int *idx, */
00640                    H5IOgetname,
00641                    &getn)<0){}
00642   return 1;
00643 }
00644 
00645 hid_t H5PartGetDiskShape(H5PartFile *f,hid_t dataset){
00646   hid_t space = H5Dget_space(dataset);
00647   if(h5part_debug) fprintf(stderr,"H5PartGetDiskShape\n");
00648   if(H5PartHasView(f)){ 
00649     int r;
00650     hsize_t total,  stride, count;
00651     hssize_t range[2];
00652     if(h5part_debug) fprintf(stderr,"\tSelection is available\n");
00653     /* so, is this selection inclusive or exclusive? */
00654     range[0]=f->viewstart;
00655     range[1]=f->viewend;
00656     count = range[1]-range[0]+1; /* to be inclusive */
00657     stride=1;
00658     /* now we select a subset */
00659     if(f->diskshape>0)
00660       r=H5Sselect_hyperslab(f->diskshape,H5S_SELECT_SET,
00661                             range/* only using first element */,
00662                             &stride,&count,NULL);
00663     /* now we select a subset */
00664     r=H5Sselect_hyperslab(space,H5S_SELECT_SET,
00665                           range,&stride,&count,NULL);
00666     if(h5part_debug) fprintf(stderr,"\tSelection: range=%u:%u, npoints=%u s=%u\n",
00667             (int)range[0],(int)range[1],
00668             (int)H5Sget_simple_extent_npoints(space),
00669             (int)H5Sget_select_npoints(space));
00670     if(r<0){
00671       fprintf(stderr,"Abort: Selection Failed!\n");
00672       return space;
00673     }
00674   }
00675   else {
00676     if(h5part_debug) fprintf(stderr,"\tNo Selection\n");
00677   }
00678   return space;
00679 }
00680 
00681 hid_t H5PartGetMemShape(H5PartFile *f,hid_t dataset){  
00682   if(h5part_debug) fprintf(stderr,"H5PartGetMemShape\n");
00683   if(H5PartHasView(f)) {
00684     hsize_t len = f->viewend - f->viewstart+1;
00685     return H5Screate_simple(1,&len,NULL);
00686   }
00687   else {
00688     return H5S_ALL;
00689   }
00690 }
00691 /*
00692   H5PartGetNumParticles:  This gets the number of particles 
00693     that are stored in the current timestep.  It will arbitrarily
00694     select a timestep if you haven't already set the timestep
00695     with H5PartSetStep().
00696  */
00697 long long H5PartGetNumParticles(H5PartFile *f){
00698   hid_t space,dataset;
00699   hsize_t nparticles;
00700   if(f->timegroup<0) {
00701     if(f->timestep<0)
00702       H5PartSetStep(f,0);
00703     else
00704       H5PartSetStep(f,f->timestep); /* choose a step */
00705   }
00706   if(f->timegroup<=0){
00707     fprintf(stderr,"Damnit!!! tried to select a timestep %d\n", f->timestep);
00708     exit(0);
00709   }
00710   /* lets open up a dataset in this group... we know "x" is there */
00711   dataset=H5Dopen(f->timegroup,"x");/* choice of X is a kludge of sorts */
00712   /*  space = H5Dget_space(dataset); */
00713   /* need to use H5PartGetDiskShape for any changes to f->diskshape */
00714   space = H5PartGetDiskShape(f,dataset); /* kludge (or utility depending on point of view) */
00715   if(H5PartHasView(f)){
00716     if(h5part_debug) fprintf(stderr,"\tGetNumParticles::get_selection\n");
00717     nparticles=H5Sget_select_npoints(space);
00718   }
00719   else {
00720     if(h5part_debug) fprintf(stderr,"\tGetNumParticles::get_simple_extent from space\n");
00721     nparticles= H5Sget_simple_extent_npoints(space);
00722   }
00723   if(space!=H5S_ALL) H5Sclose(space); /* release data shape handle */
00724   H5Dclose(dataset); /* release the dataset handle */
00725   return (long long)nparticles;
00726 }
00727 
00728 
00729 /*
00730   H5SetView: For parallel I/O or for subsetting operations on
00731   the datafile, the H5SetView subroutine allows you to define
00732   a subset of the total particle dataset to read.  The concept
00733   of "view" works for both serial and for parallel I/O.  The
00734   "view" will remain in effect until a new view is set, or the
00735   number of particles in a dataset changes, or the view is
00736     "unset" by calling H5SetView(file,-1,-1);
00737 
00738     Before you set a view, the H5PartGetNumParticles will
00739     return the total number of particles in a file (even for
00740     the parallel reads).  However, after you set a view, it
00741     will return the number of particles contained in the view.
00742 */
00743 void H5PartSetView(H5PartFile *f,long long start,long long end){
00744   hsize_t total,  stride, count;
00745   hssize_t range[2];
00746   int r;
00747   if(h5part_debug) fprintf(stderr,"SetView(%d,%d)\n",
00748           (int)start,(int)end);
00749   if(f->mode==H5PART_WRITE){
00750     fprintf(stderr,"H5PartSetView(): SetView does not make sense for write-only files.  It is meant to be used for read-only files. (maybe later this will change)\n");
00751     return;
00752   }
00753   /* if there is already a view selected, lets destroy it */ 
00754   f->viewstart=-1;
00755   f->viewend=-1;
00756   if(f->shape!=0){
00757     H5Sclose(f->shape);
00758     f->shape=0;
00759   }
00760   if(f->diskshape!=0 && f->diskshape!=H5S_ALL){
00761     H5Sclose(f->diskshape);
00762     f->diskshape=H5S_ALL;
00763   }
00764   f->diskshape = H5S_ALL;
00765   if(f->memshape!=0 && f->memshape!=H5S_ALL){
00766     H5Sclose(f->memshape);
00767     f->memshape=H5S_ALL;
00768   }
00769   if(start==-1 && end==-1) {
00770     if(h5part_debug) fprintf(stderr,"\tEarly Termination: Unsetting View\n");
00771     return; /* all done */
00772   }
00773   /* for now, we interpret start=-1 to mean 0 and 
00774      end==-1 to mean end of file */
00775   total=H5PartGetNumParticles(f);
00776   if(h5part_debug) fprintf(stderr,"\tTotal nparticles=%u\n",(int)total);
00777   if(start==-1) start=0;
00778   if(end==-1) end=total-1; /* can we trust nparticles (no)? 
00779                                                fortunately, view has been reset
00780                                               so H5PartGetNumParticles will tell
00781                                               us the total number of particles. */
00782 
00783   /* so, is this selection inclusive or exclusive? */
00784   range[0]=start;
00785   range[1]=end;
00786   /* setting up the new view */
00787   f->viewstart=range[0];
00788   f->viewend=range[1];
00789   f->nparticles=range[1]-range[0]+1;
00790   if(h5part_debug) fprintf(stderr,"\tRange is now %u:%u\n",(int)range[0],(int)range[1]);
00791   /* OK, now we must create a selection from this */
00792 
00793   /* how to check shape... is default H5S_ALL? */
00794   f->shape = H5Screate_simple(1,&total,NULL); /* declare overall datasize */
00795   /*  f->diskshape = H5S_ALL; */
00796   f->diskshape= H5Screate_simple(1,&total,NULL); /* declare overall data size  but then will select a subset */
00797   f->memshape = H5Screate_simple(1,&(f->nparticles),NULL); /* declare local memory datasize */
00798   
00799   if(f->shape<0 || f->memshape<0 || f->diskshape<0){
00800     fprintf(stderr,"Abort: shape construction failed\n");
00801     exit(0);
00802   }
00803   count = range[1]-range[0];
00804   if(h5part_debug) fprintf(stderr,"\tcount=%d\n",(int)count);
00805   stride=1;
00806   /* now we select a subset */
00807   r=H5Sselect_hyperslab(f->diskshape,H5S_SELECT_SET,range,&stride,&count,NULL);
00808   if(r<0){
00809     fprintf(stderr,"Abort: Selection Failed!\n");
00810     exit(0);
00811   }
00812      /* OK, now we have selected a reasonable hyperslab (all done) */
00813 }
00814 
00815 /*
00816   H5PartGetView: Allows you to query the current view. Start and End
00817   will be -1 if there is no current view established.
00818  */
00819 int H5PartGetView(H5PartFile *f,long long *start,long long *end){
00820   long long range[2];
00821   range[0]=(f->viewstart>=0)?f->viewstart:0;
00822   range[1]=(f->viewend>=0)?f->viewend:H5PartGetNumParticles(f);
00823   if(start) {
00824     *start=range[0];
00825   }
00826   if(end) {
00827     *end=range[1];
00828   }
00829   /* we could return the number of elements in the View as a convenience */
00830   return range[1]-range[0];
00831 }
00832 
00833 /* 
00834    H5SetCanonicalView: If it is too tedious to manually set the
00835    start and end coordinates for a view, the H5SetCanonicalView()
00836    will automatically select an appropriate domain decomposition of
00837    the data arrays for the degree of parallelism and set the "view" 
00838    accordingly.
00839 */
00840 void H5PartSetCanonicalView(H5PartFile *f){
00841   /* if a read_only file, search for on-disk canonical view */
00842   /* if this view does not exist, then if MPI, subdivide by numprocs */
00843   /* else, "unset" any existing View */
00844   if(f->mode != H5PART_READ){
00845     fprintf(stderr,"H5PartSetCanonicalView(): Views do not make sense for write-only files.  It is meant to be used for read-only files. (maybe later this will change)\n");
00846     return;
00847   }
00848   H5PartSetView(f,-1,-1); /* unset the view */
00849 #ifdef PARALLEL_IO
00850   if(f->timegroup<0){
00851     H5PartSetStep(f,0); /* set to first step in file */
00852   }
00853   /* 
00854      now lets query the attributes for this group to see if there is
00855      a 'pnparticles' group that contains the offsets for the processors.
00856   */
00857   
00858   if(H5PartReadStepAttrib(f,"pnparticles",f->pnparticles)<0){ /* try to read pnparticles right off of the disk */
00859     /* automatically subdivide the view into NP mostly equal pieces */
00860     int i;
00861     long long total=0, n = H5PartGetNumParticles(f);
00862     n/=f->nprocs;
00863     f->pnparticles[0]=n;
00864     total=n;
00865     for(i=1;i<f->nprocs;i++){
00866       f->pnparticles[i]=n;
00867       total+=n;
00868     }
00869   }
00870   {
00871     int i;
00872     long long total = 0, n = H5PartGetNumParticles(f);
00873     /* now we set the view for this processor */
00874     for(i=0;i<f->myproc;i++){
00875       total+=f->pnparticles[i];
00876     }
00877     H5PartSetView(f,total,total+f->pnparticles[f->myproc]-1);
00878   }
00879 #endif
00880   /* the canonical view is to see everything if this is serial
00881      so there is nothing left to do */
00882 }
00883 
00884 /*
00885   H5PartReadArray:  This reads in an individual array from a
00886     particlar timestep.  Unlike its "Write" mode sibling, this routine
00887     is completely public and will not have any wacky state dependencies.
00888     If you haven't selected a particular timestep, it will pick
00889     the current one.
00890 */
00891 int H5PartReadDataFloat64(H5PartFile *f,char *name,double *array){
00892   hid_t space,memspace,dataset,datatype;
00893   if(!f->timegroup) H5PartSetStep(f,f->timestep); /* choose current step */
00894   dataset=H5Dopen(f->timegroup,name);
00895   space = H5PartGetDiskShape(f,dataset); /* gets space with selection if view is set */
00896   memspace = H5PartGetMemShape(f,dataset);
00897   /*  datatype=H5Dget_type(dataset); */
00898   H5Dread(dataset, /* handle for the dataset */
00899           H5T_NATIVE_DOUBLE, /* the datatype we use in memory 
00900                               you can change it to FLOAT if you want */
00901           memspace, /* shape/size of data in memory (the complement to disk hyperslab) */
00902           space, /* shape/size of data on disk  (get hyperslab if needed) */
00903           H5P_DEFAULT,/* ignore... its for parallel reads */
00904           array); /* the data array we are reading into */
00905   if(space!=H5S_ALL) H5Sclose(space); /* release data shape handle */
00906   if(memspace!=H5S_ALL) H5Sclose(memspace);
00907   H5Dclose(dataset); /* release the dataset handle */
00908   return 1;
00909 }
00910 
00911 int H5PartReadDataInt64(H5PartFile *f,char *name,long long *array){
00912   hid_t space,memspace,dataset,datatype;
00913   if(!f->timegroup) H5PartSetStep(f,f->timestep); /* choose a step */
00914   dataset=H5Dopen(f->timegroup,name);
00915   space = H5PartGetDiskShape(f,dataset); /* H5Dget_space(dataset); */
00916   memspace = H5PartGetMemShape(f,dataset);
00917   /*  datatype=H5Dget_type(dataset); */
00918   H5Dread(dataset, /* handle for the dataset */
00919           H5T_NATIVE_INT64, /* the datatype we use in memory 
00920                               you can change it to FLOAT if you want */
00921           memspace, /* shape/size of data in memory (complement to disk hyperslab) */
00922           space, /* shape/size of data on disk  (currently get all) */
00923           H5P_DEFAULT,/* ignore... its for parallel reads */
00924           array); /* the data array we are reading into */
00925   if(space!=H5S_ALL) H5Sclose(space); /* release data shape handle */
00926   if(memspace!=H5S_ALL) H5Sclose(memspace);
00927   H5Dclose(dataset); /* release the dataset handle */
00928   return 1; /* totally bogus */
00929 }
00930 
00931 /*
00932   H5PartReadParticleStep:  This is the mongo read function that
00933     pulls in all of the data for a given timestep in one shot.
00934     It also takes the timestep as an argument and will call
00935     H5PartSetStep() internally so that you don't have to 
00936     make that call separately.
00937     See also: H5PartReadArray() if you want to just
00938     read in one of the many arrays.
00939 */
00940 int H5PartReadParticleStep(H5PartFile *f,
00941                            int step,
00942                            double *x,double *y,double *z,
00943                            double *px,double *py,double *pz,
00944                            long long *id){
00945   H5PartSetStep(f,step);
00946   /* or smuggle it into the array names */
00947   H5PartReadDataFloat64(f,"x",x);
00948   H5PartReadDataFloat64(f,"y",y);
00949   H5PartReadDataFloat64(f,"z",z);
00950   H5PartReadDataFloat64(f,"px",px);
00951   H5PartReadDataFloat64(f,"py",py);
00952   H5PartReadDataFloat64(f,"pz",pz);
00953   H5PartReadDataInt64(f,"id",id);
00954   return 1;
00955 }
00956 

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