Actual source code: viewerhdf5impl.h
1: #pragma once
3: #if defined(H5_VERSION)
4: #error "viewerhdf5impl.h must be included *before* any other HDF5 headers"
5: #else
6: #define H5_USE_18_API
7: #endif
8: #include <petscviewerhdf5.h>
9: #include <petsc/private/viewerimpl.h>
11: #if defined(PETSC_HAVE_HDF5)
13: /*
14: HDF5 function specifications usually read:
15: Returns a non-negative value if successful; otherwise returns a negative value.
16: (see e.g. https://support.hdfgroup.org/HDF5/doc/RM/RM_H5O.html#Object-Close)
17: */
18: #define PetscCallHDF5(func, args) \
19: do { \
20: herr_t _status; \
21: PetscStackPushExternal(#func); \
22: _status = func args; \
23: PetscStackPop; \
24: PetscCheck(_status >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HDF5 call %s() Status %d", #func, (int)_status); \
25: } while (0)
27: #define PetscCallHDF5ReturnNoCheck(ret, func, args) \
28: do { \
29: PetscStackPushExternal(#func); \
30: ret = func args; \
31: PetscStackPop; \
32: } while (0)
34: #define PetscCallHDF5Return(ret, func, args) \
35: do { \
36: PetscCallHDF5ReturnNoCheck(ret, func, args); \
37: PetscCheck(ret >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in HDF5 call %s() Status %d", #func, (int)ret); \
38: } while (0)
40: typedef struct PetscViewerHDF5GroupList {
41: const char *name;
42: struct PetscViewerHDF5GroupList *next;
43: } PetscViewerHDF5GroupList;
45: typedef struct {
46: char *filename;
47: PetscFileMode btype;
48: hid_t file_id;
49: hid_t dxpl_id; /* H5P_DATASET_XFER property list controlling raw data transfer (read/write). Properties are modified using H5Pset_dxpl_* functions. */
50: PetscBool timestepping; /* Flag to indicate whether objects are stored by time index */
51: PetscInt timestep; /* The time index to look for an object at */
52: PetscBool defTimestepping; /* Default if timestepping attribute is not found. Support for legacy files with no timestepping attribute */
53: PetscViewerHDF5GroupList *groups;
54: PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
55: PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */
56: PetscBool horizontal; /* store column vectors as blocks (needed for MATDENSE I/O) */
57: } PetscViewer_HDF5;
59: PETSC_SINGLE_LIBRARY_INTERN PetscErrorCode PetscViewerHDF5CheckTimestepping_Internal(PetscViewer, const char[]); /* currently used in src/dm/impls/da/gr2.c so needs to be extern */
61: /* DMPlex-specific support */
62: #define DMPLEX_STORAGE_VERSION_READING_KEY "_dm_plex_storage_version_reading"
63: #define DMPLEX_STORAGE_VERSION_WRITING_KEY "_dm_plex_storage_version_writing"
65: static inline PetscErrorCode PetscViewerHDF5ResetAttachedDMPlexStorageVersion(PetscViewer v)
66: {
67: PetscFunctionBegin;
68: PetscCall(PetscObjectCompose((PetscObject)v, DMPLEX_STORAGE_VERSION_READING_KEY, NULL));
69: PetscCall(PetscObjectCompose((PetscObject)v, DMPLEX_STORAGE_VERSION_WRITING_KEY, NULL));
70: PetscFunctionReturn(PETSC_SUCCESS);
71: }
72: #endif