xref: /petsc/include/petscviewerhdf5.h (revision 535716cbde47d9d3b0faeb019505407dc3bf6fc7)
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 
24d70abbfaSBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
25d70abbfaSBarry Smith {
26d70abbfaSBarry Smith   PetscFunctionBegin;
2784ccb19eSJed Brown   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot convert negative size");
2884ccb19eSJed Brown #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4)
2984ccb19eSJed Brown   if (a > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
30d70abbfaSBarry Smith #endif
31d70abbfaSBarry Smith   *b =  (hsize_t)(a);
32d70abbfaSBarry Smith   PetscFunctionReturn(0);
33d70abbfaSBarry Smith }
3436bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
3536bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);
362c1a2d08SJed Brown 
37ecce1506SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
38ad0c61feSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
3936bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
40e7bdbf8eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer,const char[],const char[],PetscBool*);
41577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*);
42577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
43577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);
44d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float*,int,int*,int);
45d70abbfaSBarry Smith 
46d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
47be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char[]);
48d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
49be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer,const char*[]);
50*535716cbSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,const char[],PetscBool*);
51d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
52d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
53d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);
54d70abbfaSBarry Smith 
5582ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
5682ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);
5782ea9b62SBarry Smith 
589a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
599a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);
60ee8b9732SVaclav Hapla 
61ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer,PetscBool);
62ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer,PetscBool*);
637847244bSVaclav Hapla #endif  /* defined(PETSC_HAVE_HDF5) */
64d70abbfaSBarry Smith #endif
65