xref: /petsc/include/petscviewerhdf5.h (revision ee8b9732ea65a94c6a91384e50cdc92f602fa401)
1d70abbfaSBarry Smith 
2d70abbfaSBarry Smith #if !defined(__PETSCVIEWERHDF5_H)
3d70abbfaSBarry 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>
9*ee8b9732SVaclav Hapla #if !defined(H5_VERSION_GE)
10*ee8b9732SVaclav Hapla /* H5_VERSION_GE was introduced in HDF5 1.8.7, we support >= 1.8.0 */
11*ee8b9732SVaclav Hapla /* So beware this will automatically 0 for HDF5 1.8.0 - 1.8.6 */
12*ee8b9732SVaclav Hapla #define H5_VERSION_GE(a,b,c) 0
13*ee8b9732SVaclav Hapla #endif
14*ee8b9732SVaclav 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 
21d70abbfaSBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
22d70abbfaSBarry Smith {
23d70abbfaSBarry Smith   PetscFunctionBegin;
2484ccb19eSJed Brown   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot convert negative size");
2584ccb19eSJed Brown #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4)
2684ccb19eSJed Brown   if (a > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
27d70abbfaSBarry Smith #endif
28d70abbfaSBarry Smith   *b =  (hsize_t)(a);
29d70abbfaSBarry Smith   PetscFunctionReturn(0);
30d70abbfaSBarry Smith }
3136bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
3236bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);
332c1a2d08SJed Brown 
3470efba33SBarry Smith #define PetscStackCallHDF5(func,args) do {                        \
35729ab6d9SBarry Smith     herr_t _status;                                               \
36729ab6d9SBarry Smith     PetscStackPush(#func);_status = func args;PetscStackPop; if (_status) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HDF5 call %s() Status %d",#func,(int)_status); \
37729ab6d9SBarry Smith   } while (0)
38729ab6d9SBarry Smith 
39729ab6d9SBarry Smith #define PetscStackCallHDF5Return(ret,func,args) do {              \
40729ab6d9SBarry Smith     PetscStackPush(#func);ret = func args;PetscStackPop; if (ret < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Error in HDF5 call %s() Status %d",#func,(int)ret); \
4170efba33SBarry Smith   } while (0)
4270efba33SBarry Smith 
432c1a2d08SJed Brown #endif  /* defined(PETSC_HAVE_HDF5) */
44d70abbfaSBarry Smith 
45ecce1506SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
46ad0c61feSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
4736bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
48e7bdbf8eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *);
49577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*);
50577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
51577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);
52d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float *,int,int *,int);
53d70abbfaSBarry Smith 
54d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
55d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char *);
56d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
57d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char **);
580a9f38efSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,PetscBool*);
59d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
60d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
61d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);
62d70abbfaSBarry Smith 
6382ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
6482ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);
6582ea9b62SBarry Smith 
669a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
679a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);
68*ee8b9732SVaclav Hapla 
69*ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer,PetscBool);
70*ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer,PetscBool*);
71d70abbfaSBarry Smith #endif
72