xref: /petsc/include/petscviewerhdf5.h (revision 89e0ef10bc51934e7d3e7d7ba7751b56bc779fa0)
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 
244302210dSVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5PathIsRelative(const char path[], PetscBool emptyIsRelative, PetscBool *has)
254302210dSVaclav Hapla {
264302210dSVaclav Hapla   PetscErrorCode ierr;
274302210dSVaclav Hapla   size_t         len;
284302210dSVaclav Hapla 
294302210dSVaclav Hapla   PetscFunctionBegin;
304302210dSVaclav Hapla   *has = emptyIsRelative;
314302210dSVaclav Hapla   ierr = PetscStrlen(path,&len);CHKERRQ(ierr);
324302210dSVaclav Hapla   if (len) {
334302210dSVaclav Hapla     *has = (PetscBool) (path[0] != '/');
344302210dSVaclav Hapla   }
354302210dSVaclav Hapla   PetscFunctionReturn(0);
364302210dSVaclav Hapla }
374302210dSVaclav Hapla 
38d70abbfaSBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
39d70abbfaSBarry Smith {
40d70abbfaSBarry Smith   PetscFunctionBegin;
4184ccb19eSJed Brown   if (a < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Cannot convert negative size");
4284ccb19eSJed Brown #if defined(PETSC_USE_64BIT_INDICES) && (H5_SIZEOF_HSIZE_T == 4)
4384ccb19eSJed Brown   if (a > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
44d70abbfaSBarry Smith #endif
45d70abbfaSBarry Smith   *b =  (hsize_t)(a);
46d70abbfaSBarry Smith   PetscFunctionReturn(0);
47d70abbfaSBarry Smith }
4836bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
4936bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);
502c1a2d08SJed Brown 
51*89e0ef10SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer,const char[],PetscBool*);
52ecce1506SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObject(PetscViewer,PetscObject,PetscBool*);
53ad0c61feSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
5436bce6e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
55e7bdbf8eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer,const char[],const char[],PetscBool*);
56577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,void*);
57577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer,PetscObject,const char[],PetscDataType,const void*);
58577e0e04SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer,PetscObject,const char[],PetscBool*);
59d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float*,int,int*,int);
60d70abbfaSBarry Smith 
61d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
62be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char[]);
63d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
64be7aa430SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer,const char*[]);
65535716cbSVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer,const char[],PetscBool*);
66d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
67d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
68d70abbfaSBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);
69d70abbfaSBarry Smith 
7082ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
7182ea9b62SBarry Smith PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);
7282ea9b62SBarry Smith 
739a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
749a0502c6SHåkon Strandenes PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);
75ee8b9732SVaclav Hapla 
76ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer,PetscBool);
77ee8b9732SVaclav Hapla PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer,PetscBool*);
787847244bSVaclav Hapla #endif  /* defined(PETSC_HAVE_HDF5) */
79d70abbfaSBarry Smith #endif
80