1 2 #ifndef 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 // HDF5 1.13.0 switched from hsize_t being typedef'd to unsigned long long to being uint64_t and introduced the 16 // PRIuHSIZE macro for printing. Definition of PRIuHSIZE actually precedes the change of typedef in the Git history, 17 // though there was never a release with the old definitions. Nonetheless, the logic here will work any commit. 18 #if !defined(PRIuHSIZE) 19 #define PRIuHSIZE "llu" 20 #endif 21 22 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer, hid_t *); 23 24 /* On 32-bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */ 25 #define PETSC_HDF5_INT_MAX (~(hsize_t)0) 26 27 /* As per https://portal.hdfgroup.org/display/HDF5/Chunking+in+HDF5, max. chunk size is 4GB */ 28 #define PETSC_HDF5_MAX_CHUNKSIZE 2147483647 29 30 static inline PetscErrorCode PetscViewerHDF5PathIsRelative(const char path[], PetscBool emptyIsRelative, PetscBool *has) 31 { 32 size_t len; 33 34 PetscFunctionBegin; 35 *has = emptyIsRelative; 36 PetscCall(PetscStrlen(path, &len)); 37 if (len) *has = (PetscBool)(path[0] != '/'); 38 PetscFunctionReturn(PETSC_SUCCESS); 39 } 40 41 static inline PetscErrorCode PetscHDF5IntCast(PetscInt a, hsize_t *b) 42 { 43 PetscFunctionBegin; 44 PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot convert negative size"); 45 #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4) 46 PetscCheck(a >= PETSC_HDF5_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Array too long for HDF5"); 47 #endif 48 *b = (hsize_t)(a); 49 PetscFunctionReturn(PETSC_SUCCESS); 50 } 51 PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType, hid_t *); 52 PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t, PetscDataType *); 53 54 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer, const char[], PetscBool *); 55 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer, PetscObject, PetscBool *); 56 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer, const char[], const char[], PetscDataType, const void *, void *); 57 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer, const char[], const char[], PetscDataType, const void *); 58 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *); 59 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer, PetscObject, const char[], PetscDataType, void *, void *); 60 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer, PetscObject, const char[], PetscDataType, const void *); 61 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer, PetscObject, const char[], PetscBool *); 62 63 PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm, const char[], PetscFileMode, PetscViewer *); 64 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer, const char[]); 65 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer); 66 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char[], char *[]); 67 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer, const char[], PetscBool *); 68 PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, const char[], hid_t *, hid_t *); 69 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer, const char[]); 70 71 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer, PetscBool); 72 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer, PetscBool *); 73 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer); 74 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer); 75 PETSC_EXTERN PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer, PetscBool *); 76 PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer); 77 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer, PetscInt); 78 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer, PetscInt *); 79 80 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer, PetscBool); 81 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer, PetscBool *); 82 83 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer, PetscBool); 84 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer, PetscBool *); 85 86 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer, PetscBool); 87 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer, PetscBool *); 88 #endif /* defined(PETSC_HAVE_HDF5) */ 89 #endif 90