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