xref: /petsc/include/petscviewerhdf5.h (revision fc8a9adeb7fcdc98711d755fa2dc544ddccf0f3e)
1 
2 #if !defined(__PETSCVIEWERHDF5_H)
3 #define __PETSCVIEWERHDF5_H
4 
5 #include <petscviewer.h>
6 
7 #if defined(PETSC_HAVE_HDF5)
8 
9 #include <hdf5.h>
10 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer,hid_t*);
11 PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *);
12 
13 /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
14 #define PETSC_HDF5_INT_MAX  2147483647
15 #define PETSC_HDF5_INT_MIN -2147483647
16 
17 PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
18 {
19   PetscFunctionBegin;
20 #if defined(PETSC_USE_64BIT_INDICES) && (PETSC_SIZEOF_SIZE_T == 4)
21   if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
22 #endif
23   *b =  (hsize_t)(a);
24   PetscFunctionReturn(0);
25 }
26 PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
27 PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);
28 
29 #define PetscStackCallHDF5(func,args) do {                        \
30     herr_t _status;                                               \
31     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); \
32   } while (0)
33 
34 #define PetscStackCallHDF5Return(ret,func,args) do {              \
35     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); \
36   } while (0)
37 
38 #endif  /* defined(PETSC_HAVE_HDF5) */
39 
40 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
41 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
42 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
43 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *);
44 PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*);
45 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
46 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);
47 PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float *,int,int *,int);
48 
49 PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
50 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char *);
51 PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
52 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char **);
53 PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,PetscBool*);
54 PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
55 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
56 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);
57 
58 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
59 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);
60 
61 PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
62 PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);
63 #endif
64