xref: /petsc/include/petscviewerhdf5.h (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1d70abbfaSBarry Smith 
226bd1501SBarry Smith #if !defined(PETSCVIEWERHDF5_H)
326bd1501SBarry Smith #define PETSCVIEWERHDF5_H
4d70abbfaSBarry Smith 
5d70abbfaSBarry Smith #include <petscviewer.h>
6d70abbfaSBarry Smith 
72c1a2d08SJed Brown #if defined(PETSC_HAVE_HDF5)
8d70abbfaSBarry Smith #include <hdf5.h>
9ee8b9732SVaclav Hapla #if !defined(H5_VERSION_GE)
10ee8b9732SVaclav Hapla /* H5_VERSION_GE was introduced in HDF5 1.8.7, we support >= 1.8.0 */
11ee8b9732SVaclav Hapla /* So beware this will automatically 0 for HDF5 1.8.0 - 1.8.6 */
12ee8b9732SVaclav Hapla #define H5_VERSION_GE(a, b, c) 0
13ee8b9732SVaclav Hapla #endif
14ee8b9732SVaclav Hapla 
15d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer, hid_t *);
16d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *);
17d70abbfaSBarry Smith 
18d70abbfaSBarry Smith /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
1984ccb19eSJed Brown #define PETSC_HDF5_INT_MAX (~(hsize_t)0)
20d70abbfaSBarry Smith 
2187de8cd9SSajid Ali /* As per https://portal.hdfgroup.org/display/HDF5/Chunking+in+HDF5, max. chunk size is 4GB */
2287de8cd9SSajid Ali #define PETSC_HDF5_MAX_CHUNKSIZE 2147483647
2387de8cd9SSajid Ali 
24*9371c9d4SSatish Balay static inline PetscErrorCode PetscViewerHDF5PathIsRelative(const char path[], PetscBool emptyIsRelative, PetscBool *has) {
254302210dSVaclav Hapla   size_t len;
264302210dSVaclav Hapla 
274302210dSVaclav Hapla   PetscFunctionBegin;
284302210dSVaclav Hapla   *has = emptyIsRelative;
299566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(path, &len));
305f80ce2aSJacob Faibussowitsch   if (len) *has = (PetscBool)(path[0] != '/');
314302210dSVaclav Hapla   PetscFunctionReturn(0);
324302210dSVaclav Hapla }
334302210dSVaclav Hapla 
34*9371c9d4SSatish Balay static inline PetscErrorCode PetscHDF5IntCast(PetscInt a, hsize_t *b) {
35d70abbfaSBarry Smith   PetscFunctionBegin;
362c71b3e2SJacob Faibussowitsch   PetscCheck(a >= 0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Cannot convert negative size");
3784ccb19eSJed Brown #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4)
382c71b3e2SJacob Faibussowitsch   PetscCheck(a >= PETSC_HDF5_INT_MAX, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Array too long for HDF5");
39d70abbfaSBarry Smith #endif
40d70abbfaSBarry Smith   *b = (hsize_t)(a);
41d70abbfaSBarry Smith   PetscFunctionReturn(0);
42d70abbfaSBarry Smith }
4336bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType, hid_t *);
4436bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t, PetscDataType *);
452c1a2d08SJed Brown 
4689e0ef10SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer, const char[], PetscBool *);
47ecce1506SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer, PetscObject, PetscBool *);
48d70ec8faSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer, const char[], const char[], PetscDataType, const void *, void *);
4936bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer, const char[], const char[], PetscDataType, const void *);
50e7bdbf8eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *);
51a2d6be1bSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer, PetscObject, const char[], PetscDataType, void *, void *);
52577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer, PetscObject, const char[], PetscDataType, const void *);
53577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer, PetscObject, const char[], PetscBool *);
54d70abbfaSBarry Smith 
55d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm, const char[], PetscFileMode, PetscViewer *);
56be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer, const char[]);
571ad4d0e4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroupRelative(PetscViewer, const char[]);
58d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
59be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char *[]);
60535716cbSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer, const char[], PetscBool *);
61d7dd068bSVaclav Hapla 
6219a20e4cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer, PetscBool);
6319a20e4cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer, PetscBool *);
64d7dd068bSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer);
65d7dd068bSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer);
66d7dd068bSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer, PetscBool *);
67d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
68d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer, PetscInt);
69d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer, PetscInt *);
70d70abbfaSBarry Smith 
7182ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer, PetscBool);
7282ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer, PetscBool *);
7382ea9b62SBarry Smith 
749a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer, PetscBool);
759a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer, PetscBool *);
76ee8b9732SVaclav Hapla 
77ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer, PetscBool);
78ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer, PetscBool *);
797847244bSVaclav Hapla #endif /* defined(PETSC_HAVE_HDF5) */
80d70abbfaSBarry Smith #endif
81