Actual source code: petscviewerhdf5.h


  2: #if !defined(PETSCVIEWERHDF5_H)
  3: #define PETSCVIEWERHDF5_H

  5: #include <petscviewer.h>

  7: #if defined(PETSC_HAVE_HDF5)
  8: #include <hdf5.h>
  9: #if !defined(H5_VERSION_GE)
 10: /* H5_VERSION_GE was introduced in HDF5 1.8.7, we support >= 1.8.0 */
 11: /* So beware this will automatically 0 for HDF5 1.8.0 - 1.8.6 */
 12: #define H5_VERSION_GE(a,b,c) 0
 13: #endif

 15: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer,hid_t*);
 16: PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *);

 18: /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
 19: #define PETSC_HDF5_INT_MAX  (~(hsize_t)0)

 21: /* As per https://portal.hdfgroup.org/display/HDF5/Chunking+in+HDF5, max. chunk size is 4GB */
 22: #define PETSC_HDF5_MAX_CHUNKSIZE 2147483647

 24: PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5PathIsRelative(const char path[], PetscBool emptyIsRelative, PetscBool *has)
 25: {
 27:   size_t         len;

 30:   *has = emptyIsRelative;
 31:   PetscStrlen(path,&len);
 32:   if (len) {
 33:     *has = (PetscBool) (path[0] != '/');
 34:   }
 35:   return(0);
 36: }

 38: PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
 39: {
 41:   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot convert negative size");
 42: #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4)
 43:   if (a > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
 44: #endif
 45:   *b =  (hsize_t)(a);
 46:   return(0);
 47: }
 48: PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
 49: PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);

 51: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer,const char[],PetscBool*);
 52: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
 53: PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*,void*);
 54: PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
 55: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer,const char[],const char[],PetscBool*);
 56: PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*,void*);
 57: PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
 58: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);

 60: PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
 61: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char[]);
 62: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
 63: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer,const char*[]);
 64: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,const char[],PetscBool*);

 66: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer);
 67: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer);
 68: PETSC_EXTERN PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer,PetscBool*);
 69: PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
 70: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
 71: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);

 73: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
 74: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);

 76: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
 77: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);

 79: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer,PetscBool);
 80: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer,PetscBool*);
 81: #endif  /* defined(PETSC_HAVE_HDF5) */
 82: #endif