1 2 #if !defined(PETSCVIEWERHDF5_H) 3 #define PETSCVIEWERHDF5_H 4 5 #include <petscviewer.h> 6 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 14 15 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer,hid_t*); 16 PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *); 17 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) 20 21 /* As per https://portal.hdfgroup.org/display/HDF5/Chunking+in+HDF5, max. chunk size is 4GB */ 22 #define PETSC_HDF5_MAX_CHUNKSIZE 2147483647 23 24 static inline PetscErrorCode PetscViewerHDF5PathIsRelative(const char path[], PetscBool emptyIsRelative, PetscBool *has) 25 { 26 size_t len; 27 28 PetscFunctionBegin; 29 *has = emptyIsRelative; 30 PetscCall(PetscStrlen(path,&len)); 31 if (len) *has = (PetscBool) (path[0] != '/'); 32 PetscFunctionReturn(0); 33 } 34 35 static inline PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b) 36 { 37 PetscFunctionBegin; 38 PetscCheck(a >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot convert negative size"); 39 #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4) 40 PetscCheck(a >= PETSC_HDF5_INT_MAX,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5"); 41 #endif 42 *b = (hsize_t)(a); 43 PetscFunctionReturn(0); 44 } 45 PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*); 46 PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*); 47 48 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer,const char[],PetscBool*); 49 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*); 50 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*,void*); 51 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*); 52 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer,const char[],const char[],PetscBool*); 53 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*,void*); 54 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*); 55 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*); 56 57 PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*); 58 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char[]); 59 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer); 60 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer,const char*[]); 61 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,const char[],PetscBool*); 62 63 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer,PetscBool); 64 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer,PetscBool*); 65 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer); 66 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer); 67 PETSC_EXTERN PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer,PetscBool*); 68 PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer); 69 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt); 70 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*); 71 72 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool); 73 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*); 74 75 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool); 76 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*); 77 78 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer,PetscBool); 79 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer,PetscBool*); 80 #endif /* defined(PETSC_HAVE_HDF5) */ 81 #endif 82