Actual source code: petscviewerhdf5.h

petsc-3.6.1 2015-08-06
Report Typos and Errors

  5: #include <petscviewer.h>

  7: #if defined(PETSC_HAVE_HDF5)

  9: #include <hdf5.h>
 10: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer,hid_t*);
 11: PETSC_EXTERN PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer, hid_t *, hid_t *);
 12: PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer, const char[], PetscInt *, PetscInt *);

 14: /* On 32 bit systems HDF5 is limited by size of integer, because hsize_t is defined as size_t */
 15: #define PETSC_HDF5_INT_MAX  2147483647
 16: #define PETSC_HDF5_INT_MIN -2147483647

 20: PETSC_STATIC_INLINE PetscErrorCode PetscHDF5IntCast(PetscInt a,hsize_t *b)
 21: {
 23: #if defined(PETSC_USE_64BIT_INDICES) && (PETSC_SIZEOF_SIZE_T == 4)
 24:   if ((a) > PETSC_HDF5_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for HDF5");
 25: #endif
 26:   *b =  (hsize_t)(a);
 27:   return(0);
 28: }
 29: PETSC_EXTERN PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType,hid_t*);
 30: PETSC_EXTERN PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t,PetscDataType*);

 32: #define PetscStackCallHDF5(func,args) do {                        \
 33:     herr_t _status;                                               \
 34:     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); \
 35:   } while (0)

 37: #define PetscStackCallHDF5Return(ret,func,args) do {              \
 38:     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); \
 39:   } while (0)

 41: #endif  /* defined(PETSC_HAVE_HDF5) */

 43: PETSC_EXTERN PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer,const char[],const char[],PetscDataType,void*);
 44: PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer,const char[],const char[],PetscDataType,const void*);
 45: PETSC_EXTERN PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer, const char[], const char[], PetscBool *);
 46: PETSC_EXTERN PetscErrorCode PetscViewerHDF5WriteSDS(PetscViewer,float *,int,int *,int);

 48: PETSC_EXTERN PetscErrorCode PetscViewerHDF5Open(MPI_Comm,const char[],PetscFileMode,PetscViewer*);
 49: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer,const char *);
 50: PETSC_EXTERN PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer);
 51: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer, const char **);
 52: PETSC_EXTERN PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer);
 53: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer,PetscInt);
 54: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer,PetscInt*);

 56: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer,PetscBool);
 57: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer,PetscBool*);

 59: PETSC_EXTERN PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer,PetscBool);
 60: PETSC_EXTERN PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer,PetscBool*);

 62: /* Reset __FUNCT__ in case the user does not define it themselves */

 66: #endif