src/hdf5/H5ecloud/ORG/H5Part.h

Go to the documentation of this file.
00001 #ifndef _H5Part_H_
00002 #define _H5Part_H_
00003 
00004 #include <hdf5.h>
00005 #ifdef PARALLEL_IO
00006 #include <mpi.h>
00007 #endif
00008 
00009 /* 
00010    H5PartFile:  This is an essentially opaque datastructure that
00011    acts as the filehandle for all practical purposes.
00012    It is created by H5PartOpenFile<xx>() and destroyed by
00013    H5PartCloseFile().  
00014 */
00015 typedef struct H5PartFile {
00016   hid_t file;
00017   int timestep;
00018 
00019   hid_t timegroup;
00020   hid_t properties;
00021   int myproc;
00022   hsize_t nparticles;
00023   hid_t shape;
00024   unsigned mode;
00025   int maxstep;
00026   hid_t xfer_prop,create_prop,access_prop;
00027   hid_t diskshape,memshape; /* for parallel I/O (this is on-disk) H5S_ALL 
00028                     if serial I/O */
00029 #ifdef PARALLEL_IO
00030   long long *pnparticles;
00031   MPI_Comm comm;
00032   int nprocs,myproc;
00033 #endif
00034 }H5PartFile;
00035 
00036 #define H5PART_READ 0x01
00037 #define H5PART_WRITE 0x02
00038 
00039 /*========== File Opening/Closing ===============*/
00040 /*
00041   H5PartOpenFile:  Opens file with specified filename. 
00042     If you open with flag H5PART_WRITE, it will truncate any
00043     file with the specified filename and start writing to it.
00044     If you open with H5PART_READ, then it will open the file
00045     readonly.  I could probably do an APPEND mode as well, but
00046     just haven't done so yet.  APPEND would need to check the
00047     file to determine its current state and modify the 
00048     H5PartFile datastructure accordingly.
00049 
00050     H5PartFile should be treated as an essentially opaque
00051     datastructure.  It acts as the file handle, but internally
00052     it maintains several key state variables associated with 
00053     the file.
00054 
00055     returns: A new filehandle with an open file or NULL if error.
00056 */
00057 #ifdef PARALLEL_IO
00058 #include <mpi.h>
00059 H5PartFile *H5PartOpenFileParallel(const char *filename,
00060                                    unsigned flags,
00061                                    MPI_Comm communicator);
00062 #endif
00063 /* the pure serial version */
00064 #define H5PartOpenFileSerial(x,y) H5PartOpenFile(x,y)
00065 H5PartFile *H5PartOpenFile(const char *filename, /* name of datafile */
00066                                  unsigned flags); /* H5PART_READ | H5PART_WRITE */
00067 
00068 int H5PartFileIsValid(H5PartFile *f);
00069 
00070 /*
00071   H5PartCloseFile:  Does what it says it does.
00072 */
00073 void H5PartCloseFile(H5PartFile *f);
00074 
00075 
00076 /*============== File Writing Functions ==================== */
00077 /*
00078   H5PartSetNumParticles:  This function's sole purpose is to 
00079     prevent needless creation of new HDF5 DataSpace handles if
00080     the number of particles is invariant throughout the sim.
00081     That's its only reason for existence.  After you call this
00082     subroutine, all subsequent operations will assume this 
00083     number of particles will be written.
00084  */
00085 void H5PartSetNumParticles(H5PartFile *f,long long nparticles);
00086 /*
00087   H5PartWriteData:  This writes a dataset with the indicated
00088   name to the currently selected Timestep in the datafile.
00089   It will be an error if you don't call
00090      H5PartSetStep()
00091   before you use this function (but once those things are set, they
00092   are persistent values... no need to call them before *every* write.)
00093    Also, you need to have set the
00094   number of particles that you are storing in this timestep using
00095      H5PartSetNumParticles()
00096   It returns a negative value if the write fails.  This either
00097   indicates that the Timestep or number of particles were not 
00098   set, *or* you are storing two datasets in the same timestep 
00099   that have the same name (not allowed).
00100  */
00101 int H5PartWriteDataFloat64(H5PartFile *f,char *name,double *array);
00102 /* same deal as above, but for 64bit integers */
00103 int H5PartWriteDataInt64(H5PartFile *f,char *name,long long *array);
00104 
00105 /*================== File Reading Routines =================*/
00106 /*
00107   H5PartSetStep:  This routine is *only* useful when you are reading
00108   data.  Using it while you are writing will generate nonsense results!
00109   (This file format is only half-baked... robustness is not std equipment!)
00110   So you use this to random-access the file for a particular timestep.
00111   Failure to explicitly set the timestep on each read will leave you
00112   stuck on the same timestep for *all* of your reads.  That is to say
00113   the writes auto-advance the file pointer, but the reads do not
00114   (they require explicit advancing by selecting a particular timestep).
00115 */
00116 void H5PartSetStep(H5PartFile *f, /* file handle */
00117                    int step); /* current timestep to select (0 to n-1) */
00118 
00119 
00120 /* 
00121    H5PartGetNumSteps:  This reads the number of datasteps that are 
00122      currently stored in the datafile.  (simple return of an int).
00123      It works for both reading and writing of files, but is probably
00124      only typically used when you are reading.
00125 
00126      returns: number of timesteps currently stored in the file.
00127 */
00128 int H5PartGetNumSteps(H5PartFile *f);
00129 
00130 /*
00131   H5PartGetNumDatasets/H5PartGetDatasetName() are both used together
00132   to find out how many datasets are stored at a particular timestep
00133   and what their names are if you don't know what they are a-priori.
00134  */
00135 int H5PartGetNumDatasets(H5PartFile *f);
00136 int H5PartGetDatasetName(H5PartFile *f,int index,char *name,int maxlen);
00137 /*
00138   H5PartGetNumParticles:  This gets the number of particles 
00139     that are stored in the current timestep.  It will arbitrarily
00140     select a timestep if you haven't already set the timestep
00141     with H5PartSetStep().
00142 
00143     returns: number of particles in current timestep
00144  */
00145 long long H5PartGetNumParticles(H5PartFile *f);
00146 
00147 
00148 /*
00149   H5PartReadArray:  This reads in an individual array from a
00150     particlar timestep.  Unlike its "Write" mode sibling, this routine
00151     is completely public and will not have any wacky state dependencies.
00152     If you haven't selected a particular timestep, it will pick
00153     the current one.  Also, it assumes you allocated enough space to
00154     read in all of the particles in this timestep.
00155 
00156     returns: ignore returnval for now (should be used for error handling)
00157 */
00158 int H5PartReadDataFloat64(H5PartFile *f,
00159                           char *name, /* name of the array to read
00160                                          "x"=position in x direction (y,z)
00161                                          "vx"=velocity in x directio (y,z)
00162                                          "px"=position in x dir (y,z) */
00163                           double *array); /* array to read data into */
00164 int H5PartReadDataInt64(H5PartFile *f,
00165                           char *name, /* name of the array to read
00166                                          "x"=position in x direction (y,z)
00167                                          "vx"=velocity in x directio (y,z)
00168                                          "px"=position in x dir (y,z) */
00169                           long long *array); /* array to read data into */
00170 
00171 /* the following is a back-door for extensions to the data writing */
00172 #if 0
00173 int H5PartReadData(H5PartFile *f,char *name,void *array,hid_t type);
00174 int H5PartWriteData(H5PartFile *f,char *name,void *array,hid_t type);
00175 #endif
00176 /*
00177   H5PartReadParticleStep:  This is the mongo read function that
00178     pulls in all of the data for a given timestep in one shot.
00179     It also takes the timestep as an argument and will call
00180     H5PartSetStep() internally so that you don't have to 
00181     make that call separately.
00182     See also: H5PartReadArray() if you want to just
00183     read in one of the many arrays.
00184 
00185     At this point, I don't think it will be used (too rigid), 
00186     but it is for example purposes.
00187 */
00188 int H5PartReadParticleStep(H5PartFile *f, /* filehandle */
00189                            int step, /* selects timestep to read from*/
00190                            double *x,double *y,double *z, /* particle positions */
00191                            double *px,double *py,double *pz, /* particle momenta */
00192                            long long *id); /* and phase */
00193 /**********==============Attributes Interface============***************/
00194 /* currently there is file attributes:  Attributes bound to the file
00195    and step attributes which are bound to the current timestep.  You 
00196    must set the timestep explicitly before writing the attributes (just
00197    as you must do when you write a new dataset.  Currently there are no
00198    attributes that are bound to a particular data array, but this could
00199    easily be done if requred.
00200 */
00201 int H5PartWriteStepAttrib(H5PartFile *f,char *name,
00202     hid_t type,void *attrib,int nelem);
00203 int H5PartWriteFileAttrib(H5PartFile *f,char *name,
00204     hid_t type,void *attrib,int nelem);
00205 int H5PartWriteAttrib(H5PartFile *f,char *name,
00206                       hid_t type,void *attrib,int nelem); // this should be deprecated
00207 int H5PartWriteFileAttribString(H5PartFile *f,char *name,
00208     char *attrib);
00209 int H5PartWriteStepAttribString(H5PartFile *f,char *name,
00210     char *attrib);
00211 int H5PartGetNumStepAttribs(H5PartFile *f); /* for current filestep */
00212 int H5PartGetNumFileAttribs(H5PartFile *f);
00213 void H5PartGetStepAttribInfo(H5PartFile *f,int idx,
00214     char *name,size_t maxnamelen,
00215                          hid_t *type,int *nelem);
00216 void H5PartGetFileAttribInfo(H5PartFile *f,int idx,
00217     char *name,size_t maxnamelen,
00218     hid_t *type,int *nelem);
00219 void H5PartReadStepAttrib(H5PartFile *f,char *name,void *data);
00220 void H5PartReadAttrib(H5PartFile *f,char *name,void *data);
00221 void H5PartReadFileAttrib(H5PartFile *f,char *name,void *data);
00222 
00223 
00224 #endif
00225 

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