Actual source code: viewerhdf5impl.h

  1: #ifndef __VIEWERHDF5IMPL_H

  4: #if defined(H5_VERSION)
  5:   #error "viewerhdf5impl.h must be included *before* any other HDF5 headers"
  6: #else
  7:   #define H5_USE_18_API
  8: #endif
  9: #include <petscviewerhdf5.h>
 10: #include <petsc/private/viewerimpl.h>

 12: #if defined(PETSC_HAVE_HDF5)

 14:   /*
 15:   HDF5 function specifications usually read:
 16:   Returns a non-negative value if successful; otherwise returns a negative value.
 17:   (see e.g. https://support.hdfgroup.org/HDF5/doc/RM/RM_H5O.html#Object-Close)
 18: */
 19:   #define PetscCallHDF5(func, args) \
 20:     do { \
 21:       herr_t _status; \
 22:       PetscStackPushExternal(#func); \
 23:       _status = func args; \
 24:       PetscStackPop; \
 25:       PetscCheck(_status >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HDF5 call %s() Status %d", #func, (int)_status); \
 26:     } while (0)

 28:   #define PetscCallHDF5ReturnNoCheck(ret, func, args) \
 29:     do { \
 30:       PetscStackPushExternal(#func); \
 31:       ret = func args; \
 32:       PetscStackPop; \
 33:     } while (0)

 35:   #define PetscCallHDF5Return(ret, func, args) \
 36:     do { \
 37:       PetscCallHDF5ReturnNoCheck(ret, func, args); \
 38:       PetscCheck(ret >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HDF5 call %s() Status %d", #func, (int)ret); \
 39:     } while (0)

 41: typedef struct PetscViewerHDF5GroupList {
 42:   const char                      *name;
 43:   struct PetscViewerHDF5GroupList *next;
 44: } PetscViewerHDF5GroupList;

 46: typedef struct {
 47:   char                     *filename;
 48:   PetscFileMode             btype;
 49:   hid_t                     file_id;
 50:   hid_t                     dxpl_id;         /* H5P_DATASET_XFER property list controlling raw data transfer (read/write). Properties are modified using H5Pset_dxpl_* functions. */
 51:   PetscBool                 timestepping;    /* Flag to indicate whether objects are stored by time index */
 52:   PetscInt                  timestep;        /* The time index to look for an object at */
 53:   PetscBool                 defTimestepping; /* Default if timestepping attribute is not found. Support for legacy files with no timestepping attribute */
 54:   PetscViewerHDF5GroupList *groups;
 55:   PetscBool                 basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
 56:   PetscBool                 spoutput;       /* write data in single precision even if PETSc is compiled with double precision PetscReal */
 57:   PetscBool                 horizontal;     /* store column vectors as blocks (needed for MATDENSE I/O) */
 58: } PetscViewer_HDF5;

 60: PETSC_EXTERN PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer, const char[]); /* currently used in src/dm/impls/da/gr2.c so needs to be extern */
 61: PETSC_INTERN PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer, const char *[]);

 63:   /* DMPlex-specific support */
 64:   #define DMPLEX_STORAGE_VERSION_READING_KEY "_dm_plex_storage_version_reading"
 65:   #define DMPLEX_STORAGE_VERSION_WRITING_KEY "_dm_plex_storage_version_writing"

 67: static inline PetscErrorCode PetscViewerHDF5ResetAttachedDMPlexStorageVersion(PetscViewer v)
 68: {
 69:   PetscFunctionBegin;
 70:   PetscCall(PetscObjectCompose((PetscObject)v, DMPLEX_STORAGE_VERSION_READING_KEY, NULL));
 71:   PetscCall(PetscObjectCompose((PetscObject)v, DMPLEX_STORAGE_VERSION_WRITING_KEY, NULL));
 72:   PetscFunctionReturn(PETSC_SUCCESS);
 73: }
 74: #endif
 75: #endif