xref: /petsc/include/petscviewerhdf5.h (revision 577e0e045ba1fcf91fa9450468fd1cbfa08e7ca8)
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)
82c1a2d08SJed Brown 
9d70abbfaSBarry Smith #include <hdf5.h>
10d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer,hid_t*);
11d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *);
127adab82aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer, const char[], PetscInt *, PetscInt *);
13d70abbfaSBarry Smith 
14d70abbfaSBarry Smith /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
15d70abbfaSBarry Smith #define PETSC_HDF5_INT_MAX  2147483647
16d70abbfaSBarry Smith #define PETSC_HDF5_INT_MIN -2147483647
17d70abbfaSBarry Smith 
18d70abbfaSBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
19d70abbfaSBarry Smith {
20d70abbfaSBarry Smith   PetscFunctionBegin;
21d70abbfaSBarry Smith #if defined(PETSC_USE_64BIT_INDICES) && (PETSC_SIZEOF_SIZE_T == 4)
22d70abbfaSBarry Smith   if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
23d70abbfaSBarry Smith #endif
24d70abbfaSBarry Smith   *b =  (hsize_t)(a);
25d70abbfaSBarry Smith   PetscFunctionReturn(0);
26d70abbfaSBarry Smith }
2736bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
2836bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);
292c1a2d08SJed Brown 
3070efba33SBarry Smith #define PetscStackCallHDF5(func,args) do {                        \
31729ab6d9SBarry Smith     herr_t _status;                                               \
32729ab6d9SBarry 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); \
33729ab6d9SBarry Smith   } while (0)
34729ab6d9SBarry Smith 
35729ab6d9SBarry Smith #define PetscStackCallHDF5Return(ret,func,args) do {              \
36729ab6d9SBarry 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); \
3770efba33SBarry Smith   } while (0)
3870efba33SBarry Smith 
392c1a2d08SJed Brown #endif  /* defined(PETSC_HAVE_HDF5) */
40d70abbfaSBarry Smith 
41ecce1506SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
42ad0c61feSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
4336bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
44e7bdbf8eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *);
45*577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*);
46*577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
47*577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);
48d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float *,int,int *,int);
49d70abbfaSBarry Smith 
50d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
51d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char *);
52d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
53d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char **);
540a9f38efSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,PetscBool*);
55d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
56d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
57d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);
58d70abbfaSBarry Smith 
5982ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
6082ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);
6182ea9b62SBarry Smith 
629a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
639a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);
649a0502c6SHåkon Strandenes 
65058bd781SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer,const char[],const char[],const char[],const char[]);
66058bd781SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer,const char*[],const char*[],const char*[],const char*[]);
67d70abbfaSBarry Smith #endif
68