Actual source code: ihdf5v.c
1: #include <petsc/private/viewerhdf5impl.h>
3: static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
4: {
5: htri_t exists;
6: hid_t group;
8: PetscFunctionBegin;
9: PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
10: if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
11: if (!exists && createGroup) {
12: hid_t plist_id;
13: PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_GROUP_CREATE));
14: PetscCallHDF5(H5Pset_link_creation_order, (plist_id, H5P_CRT_ORDER_TRACKED | H5P_CRT_ORDER_INDEXED));
15: PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, plist_id, H5P_DEFAULT));
16: PetscCallHDF5(H5Pclose, (plist_id));
17: PetscCallHDF5(H5Gclose, (group));
18: exists = PETSC_TRUE;
19: }
20: *exists_ = (PetscBool)exists;
21: PetscFunctionReturn(PETSC_SUCCESS);
22: }
24: static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
25: {
26: const char rootGroupName[] = "/";
27: hid_t h5;
28: PetscBool exists = PETSC_FALSE;
29: PetscInt i;
30: int n;
31: char **hierarchy;
32: char buf[PETSC_MAX_PATH_LEN] = "";
34: PetscFunctionBegin;
36: if (name) PetscAssertPointer(name, 2);
37: else name = rootGroupName;
38: if (has) {
39: PetscAssertPointer(has, 4);
40: *has = PETSC_FALSE;
41: }
42: if (otype) {
43: PetscAssertPointer(otype, 5);
44: *otype = H5O_TYPE_UNKNOWN;
45: }
46: PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
48: /*
49: Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
50: Hence, each of them needs to be tested separately:
51: 1) whether it's a valid link
52: 2) whether this link resolves to an object
53: See H5Oexists_by_name() documentation.
54: */
55: PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
56: if (!n) {
57: /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
58: if (has) *has = PETSC_TRUE;
59: if (otype) *otype = H5O_TYPE_GROUP;
60: PetscCall(PetscStrToArrayDestroy(n, hierarchy));
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
63: for (i = 0; i < n; i++) {
64: PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
65: PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf)));
66: PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
67: if (!exists) break;
68: }
69: PetscCall(PetscStrToArrayDestroy(n, hierarchy));
71: /* If the object exists, get its type */
72: if (exists && otype) {
73: H5O_info_t info;
75: /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
76: PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
77: *otype = info.type;
78: }
79: if (has) *has = exists;
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
83: static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
84: {
85: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
87: PetscFunctionBegin;
89: PetscAssertPointer(name, 2);
90: if (hdf5->groups) *name = hdf5->groups->name;
91: else *name = NULL;
92: PetscFunctionReturn(PETSC_SUCCESS);
93: }
95: static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
96: {
97: hid_t h5;
98: htri_t hhas;
100: PetscFunctionBegin;
101: PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
102: PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
103: *has = hhas ? PETSC_TRUE : PETSC_FALSE;
104: PetscFunctionReturn(PETSC_SUCCESS);
105: }
107: static PetscErrorCode PetscViewerHDF5GetGroup_HDF5(PetscViewer viewer, const char path[], const char *abspath[])
108: {
109: size_t len;
110: PetscBool relative = PETSC_FALSE;
111: const char *group;
112: char buf[PETSC_MAX_PATH_LEN] = "";
114: PetscFunctionBegin;
115: PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
116: PetscCall(PetscStrlen(path, &len));
117: relative = (PetscBool)(!len || path[0] != '/');
118: if (relative) {
119: PetscCall(PetscStrncpy(buf, group, sizeof(buf)));
120: if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
121: PetscCall(PetscStrlcat(buf, path, sizeof(buf)));
122: PetscCall(PetscStrallocpy(buf, (char **)abspath));
123: } else {
124: PetscCall(PetscStrallocpy(path, (char **)abspath));
125: }
126: PetscFunctionReturn(PETSC_SUCCESS);
127: }
129: static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems PetscOptionsObject)
130: {
131: PetscBool flg = PETSC_FALSE, set;
132: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
134: PetscFunctionBegin;
135: PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
136: PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
137: PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
138: PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
139: if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
140: flg = PETSC_FALSE;
141: PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
142: if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
143: PetscCall(PetscOptionsBool("-viewer_hdf5_compress", "Enable compression", "PetscViewerHDF5SetCompress", hdf5->compress, &hdf5->compress, NULL));
144: PetscOptionsHeadEnd();
145: PetscFunctionReturn(PETSC_SUCCESS);
146: }
148: static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
149: {
150: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
151: PetscBool flg;
153: PetscFunctionBegin;
154: if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
155: PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
156: PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
157: PetscCall(PetscViewerHDF5GetCollective(v, &flg));
158: PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
159: PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
160: PetscCall(PetscViewerASCIIPrintf(viewer, "Compression: %s\n", PetscBools[hdf5->compress]));
161: PetscFunctionReturn(PETSC_SUCCESS);
162: }
164: static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
165: {
166: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
168: PetscFunctionBegin;
169: PetscCall(PetscFree(hdf5->filename));
170: if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
171: PetscFunctionReturn(PETSC_SUCCESS);
172: }
174: static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
175: {
176: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
178: PetscFunctionBegin;
179: if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
180: PetscFunctionReturn(PETSC_SUCCESS);
181: }
183: static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
184: {
185: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
187: PetscFunctionBegin;
188: PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
189: PetscCall(PetscViewerFileClose_HDF5(viewer));
190: while (hdf5->groups) {
191: PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
193: PetscCall(PetscFree(hdf5->groups->name));
194: PetscCall(PetscFree(hdf5->groups));
195: hdf5->groups = tmp;
196: }
197: PetscCall(PetscFree(hdf5));
198: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
199: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
200: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
201: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
202: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
203: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetBaseDimension2_C", NULL));
204: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
205: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
206: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
207: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
208: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
209: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCompress_C", NULL));
210: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCompress_C", NULL));
211: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5WriteGroup_C", NULL));
212: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PushGroup_C", NULL));
213: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PopGroup_C", NULL));
214: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetGroup_C", NULL));
215: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PushTimestepping_C", NULL));
216: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5PopTimestepping_C", NULL));
217: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5IsTimestepping_C", NULL));
218: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5IncrementTimestep_C", NULL));
219: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetTimestep_C", NULL));
220: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetTimestep_C", NULL));
221: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5WriteAttribute_C", NULL));
222: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5ReadAttribute_C", NULL));
223: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasGroup_C", NULL));
224: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasDataset_C", NULL));
225: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5HasAttribute_C", NULL));
226: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetSPOutput_C", NULL));
227: PetscFunctionReturn(PETSC_SUCCESS);
228: }
230: static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
231: {
232: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
234: PetscFunctionBegin;
235: hdf5->btype = type;
236: PetscFunctionReturn(PETSC_SUCCESS);
237: }
239: static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
240: {
241: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
243: PetscFunctionBegin;
244: *type = hdf5->btype;
245: PetscFunctionReturn(PETSC_SUCCESS);
246: }
248: static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
249: {
250: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
252: PetscFunctionBegin;
253: hdf5->basedimension2 = flg;
254: PetscFunctionReturn(PETSC_SUCCESS);
255: }
257: static PetscErrorCode PetscViewerHDF5GetBaseDimension2_HDF5(PetscViewer viewer, PetscBool *flg)
258: {
259: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
261: PetscFunctionBegin;
262: *flg = hdf5->basedimension2;
263: PetscFunctionReturn(PETSC_SUCCESS);
264: }
266: static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
267: {
268: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
270: PetscFunctionBegin;
271: hdf5->spoutput = flg;
272: PetscFunctionReturn(PETSC_SUCCESS);
273: }
275: static PetscErrorCode PetscViewerHDF5GetSPOutput_HDF5(PetscViewer viewer, PetscBool *flg)
276: {
277: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
279: PetscFunctionBegin;
280: *flg = hdf5->spoutput;
281: PetscFunctionReturn(PETSC_SUCCESS);
282: }
284: static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
285: {
286: PetscFunctionBegin;
287: /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
288: - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
289: #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
290: {
291: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
292: PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
293: }
294: #else
295: if (flg) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n"));
296: #endif
297: PetscFunctionReturn(PETSC_SUCCESS);
298: }
300: static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
301: {
302: #if defined(H5_HAVE_PARALLEL)
303: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
304: H5FD_mpio_xfer_t mode;
305: #endif
307: PetscFunctionBegin;
308: #if !defined(H5_HAVE_PARALLEL)
309: *flg = PETSC_FALSE;
310: #else
311: PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
312: *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
313: #endif
314: PetscFunctionReturn(PETSC_SUCCESS);
315: }
317: static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
318: {
319: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
320: hid_t plist_access_id;
321: hid_t plist_create_id;
323: PetscFunctionBegin;
324: if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
325: if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
326: PetscCall(PetscStrallocpy(name, &hdf5->filename));
327: PetscCallHDF5Return(plist_create_id, H5Pcreate, (H5P_FILE_CREATE));
328: PetscCallHDF5(H5Pset_link_creation_order, (plist_create_id, H5P_CRT_ORDER_TRACKED | H5P_CRT_ORDER_INDEXED));
329: PetscCallHDF5Return(plist_access_id, H5Pcreate, (H5P_FILE_ACCESS));
330: #if defined(H5_HAVE_PARALLEL)
331: PetscCallHDF5(H5Pset_fapl_mpio, (plist_access_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
332: #endif
333: /* Create or open the file collectively */
334: switch (hdf5->btype) {
335: case FILE_MODE_READ:
336: if (PetscDefined(USE_DEBUG)) {
337: PetscMPIInt rank;
338: PetscBool flg;
340: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
341: if (rank == 0) {
342: PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
343: PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
344: }
345: PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
346: }
347: PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_access_id));
348: break;
349: case FILE_MODE_APPEND:
350: case FILE_MODE_UPDATE: {
351: PetscBool flg;
352: PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
353: if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_access_id));
354: else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, plist_create_id, plist_access_id));
355: break;
356: }
357: case FILE_MODE_WRITE:
358: PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, plist_create_id, plist_access_id));
359: break;
360: case FILE_MODE_UNDEFINED:
361: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
362: default:
363: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
364: }
365: PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
366: PetscCallHDF5(H5Pclose, (plist_access_id));
367: PetscCallHDF5(H5Pclose, (plist_create_id));
368: PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
369: PetscFunctionReturn(PETSC_SUCCESS);
370: }
372: static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
373: {
374: PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
376: PetscFunctionBegin;
377: *name = vhdf5->filename;
378: PetscFunctionReturn(PETSC_SUCCESS);
379: }
381: static PetscErrorCode PetscViewerSetUp_HDF5(PETSC_UNUSED PetscViewer viewer)
382: {
383: PetscFunctionBegin;
384: PetscFunctionReturn(PETSC_SUCCESS);
385: }
387: static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
388: {
389: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
391: PetscFunctionBegin;
392: hdf5->defTimestepping = flg;
393: PetscFunctionReturn(PETSC_SUCCESS);
394: }
396: static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
397: {
398: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
400: PetscFunctionBegin;
401: *flg = hdf5->defTimestepping;
402: PetscFunctionReturn(PETSC_SUCCESS);
403: }
405: static PetscErrorCode PetscViewerHDF5SetCompress_HDF5(PetscViewer viewer, PetscBool flg)
406: {
407: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
409: PetscFunctionBegin;
410: hdf5->compress = flg;
411: PetscFunctionReturn(PETSC_SUCCESS);
412: }
414: static PetscErrorCode PetscViewerHDF5GetCompress_HDF5(PetscViewer viewer, PetscBool *flg)
415: {
416: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
418: PetscFunctionBegin;
419: *flg = hdf5->compress;
420: PetscFunctionReturn(PETSC_SUCCESS);
421: }
423: /*@C
424: PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
426: Not Collective
428: Input Parameter:
429: . viewer - the `PetscViewer`
431: Output Parameter:
432: . file_id - The file id
434: Level: intermediate
436: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
437: @*/
438: PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
439: {
440: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
442: PetscFunctionBegin;
444: if (file_id) *file_id = hdf5->file_id;
445: PetscFunctionReturn(PETSC_SUCCESS);
446: }
448: static PetscErrorCode PetscViewerHDF5PushGroup_HDF5(PetscViewer viewer, const char name[])
449: {
450: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
451: PetscViewerHDF5GroupList *groupNode;
452: size_t i, len;
453: char buf[PETSC_MAX_PATH_LEN];
454: const char *gname;
456: PetscFunctionBegin;
457: PetscCall(PetscStrlen(name, &len));
458: gname = NULL;
459: if (len) {
460: if (len == 1 && name[0] == '.') {
461: /* use current name */
462: gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
463: } else if (name[0] == '/') {
464: /* absolute */
465: for (i = 1; i < len; i++) {
466: if (name[i] != '/') {
467: gname = name;
468: break;
469: }
470: }
471: } else {
472: /* relative */
473: const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
474: PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
475: gname = buf;
476: }
477: }
478: PetscCall(PetscNew(&groupNode));
479: PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
480: groupNode->next = hdf5->groups;
481: hdf5->groups = groupNode;
482: PetscFunctionReturn(PETSC_SUCCESS);
483: }
485: static PetscErrorCode PetscViewerHDF5PopGroup_HDF5(PetscViewer viewer)
486: {
487: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
488: PetscViewerHDF5GroupList *groupNode;
490: PetscFunctionBegin;
491: PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
492: groupNode = hdf5->groups;
493: hdf5->groups = hdf5->groups->next;
494: PetscCall(PetscFree(groupNode->name));
495: PetscCall(PetscFree(groupNode));
496: PetscFunctionReturn(PETSC_SUCCESS);
497: }
499: /*@C
500: PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
501: and return this group's ID and file ID.
502: If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
504: Not Collective
506: Input Parameters:
507: + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
508: - path - (Optional) The path relative to the pushed group
510: Output Parameters:
511: + fileId - The HDF5 file ID
512: - groupId - The HDF5 group ID
514: Level: intermediate
516: Note:
517: If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
518: `NULL` or empty path means the current pushed group.
520: If the viewer is writable, the group is created if it doesn't exist yet.
522: .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
523: @*/
524: PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
525: {
526: hid_t file_id;
527: H5O_type_t type;
528: const char *fileName = NULL;
529: const char *groupName = NULL;
530: PetscBool writable, has;
532: PetscFunctionBegin;
534: if (path) PetscAssertPointer(path, 2);
535: PetscAssertPointer(fileId, 3);
536: PetscAssertPointer(groupId, 4);
537: PetscCall(PetscViewerWritable(viewer, &writable));
538: PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
539: PetscCall(PetscViewerFileGetName(viewer, &fileName));
540: PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
541: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
542: if (!has) {
543: PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
544: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
545: }
546: PetscCheck(type == H5O_TYPE_GROUP, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName);
547: PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
548: PetscCall(PetscFree(groupName));
549: *fileId = file_id;
550: PetscFunctionReturn(PETSC_SUCCESS);
551: }
553: static PetscErrorCode PetscViewerHDF5WriteGroup_HDF5(PetscViewer viewer, const char path[])
554: {
555: hid_t fileId, groupId;
557: PetscFunctionBegin;
558: PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
559: PetscCallHDF5(H5Gclose, (groupId));
560: PetscFunctionReturn(PETSC_SUCCESS);
561: }
563: static PetscErrorCode PetscViewerHDF5PushTimestepping_HDF5(PetscViewer viewer)
564: {
565: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
567: PetscFunctionBegin;
568: PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
569: hdf5->timestepping = PETSC_TRUE;
570: if (hdf5->timestep < 0) hdf5->timestep = 0;
571: PetscFunctionReturn(PETSC_SUCCESS);
572: }
574: static PetscErrorCode PetscViewerHDF5PopTimestepping_HDF5(PetscViewer viewer)
575: {
576: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
578: PetscFunctionBegin;
579: PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
580: hdf5->timestepping = PETSC_FALSE;
581: PetscFunctionReturn(PETSC_SUCCESS);
582: }
584: static PetscErrorCode PetscViewerHDF5IsTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
585: {
586: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
588: PetscFunctionBegin;
589: *flg = hdf5->timestepping;
590: PetscFunctionReturn(PETSC_SUCCESS);
591: }
593: static PetscErrorCode PetscViewerHDF5IncrementTimestep_HDF5(PetscViewer viewer)
594: {
595: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
597: PetscFunctionBegin;
598: PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
599: ++hdf5->timestep;
600: PetscFunctionReturn(PETSC_SUCCESS);
601: }
603: static PetscErrorCode PetscViewerHDF5SetTimestep_HDF5(PetscViewer viewer, PetscInt timestep)
604: {
605: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
607: PetscFunctionBegin;
608: PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
609: PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
610: hdf5->timestep = timestep;
611: PetscFunctionReturn(PETSC_SUCCESS);
612: }
614: static PetscErrorCode PetscViewerHDF5GetTimestep_HDF5(PetscViewer viewer, PetscInt *timestep)
615: {
616: PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
618: PetscFunctionBegin;
619: PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
620: *timestep = hdf5->timestep;
621: PetscFunctionReturn(PETSC_SUCCESS);
622: }
624: /*@C
625: PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
627: Not Collective
629: Input Parameter:
630: . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
632: Output Parameter:
633: . htype - the HDF5 datatype
635: Level: advanced
637: .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
638: @*/
639: PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
640: {
641: PetscFunctionBegin;
642: if (ptype == PETSC_INT) *htype = PetscDefined(USE_64BIT_INDICES) ? H5T_NATIVE_LLONG : H5T_NATIVE_INT;
643: else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
644: else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
645: else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
646: else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
647: else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_HBOOL;
648: else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
649: else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
650: else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
651: else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
652: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
653: PetscFunctionReturn(PETSC_SUCCESS);
654: }
656: /*@C
657: PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
659: Not Collective
661: Input Parameter:
662: . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
664: Output Parameter:
665: . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
667: Level: advanced
669: .seealso: [](sec_viewers), `PetscDataType`
670: @*/
671: PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
672: {
673: PetscFunctionBegin;
674: if (htype == H5T_NATIVE_INT) *ptype = PetscDefined(USE_64BIT_INDICES) ? PETSC_LONG : PETSC_INT;
675: #if defined(PETSC_USE_64BIT_INDICES)
676: else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
677: #endif
678: else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
679: else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
680: else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
681: else if (htype == H5T_NATIVE_HBOOL) *ptype = PETSC_BOOL;
682: else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
683: else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
684: else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
685: else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
686: else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
687: PetscFunctionReturn(PETSC_SUCCESS);
688: }
690: static PetscErrorCode PetscViewerHDF5WriteAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
691: {
692: const char *parentAbsPath;
693: hid_t h5, dataspace, obj, attribute, dtype;
694: PetscBool has;
696: PetscFunctionBegin;
697: PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
698: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
699: PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
700: PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
701: if (datatype == PETSC_STRING) {
702: size_t len;
703: PetscCall(PetscStrlen((const char *)value, &len));
704: PetscCallHDF5(H5Tset_size, (dtype, len + 1));
705: }
706: PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
707: PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
708: PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
709: if (has) {
710: PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
711: } else {
712: PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
713: }
714: PetscCallHDF5(H5Awrite, (attribute, dtype, value));
715: if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
716: PetscCallHDF5(H5Aclose, (attribute));
717: PetscCallHDF5(H5Oclose, (obj));
718: PetscCallHDF5(H5Sclose, (dataspace));
719: PetscCall(PetscFree(parentAbsPath));
720: PetscFunctionReturn(PETSC_SUCCESS);
721: }
723: static PetscErrorCode PetscViewerHDF5ReadAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
724: {
725: const char *parentAbsPath;
726: hid_t h5, obj, attribute, dtype;
727: PetscBool has;
729: PetscFunctionBegin;
730: PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
731: PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
732: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
733: if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
734: if (!has) {
735: if (defaultValue) {
736: if (defaultValue != value) {
737: if (datatype == PETSC_STRING) {
738: PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
739: } else {
740: size_t len;
741: PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
742: PetscCall(PetscMemcpy(value, defaultValue, len));
743: }
744: }
745: PetscCall(PetscFree(parentAbsPath));
746: PetscFunctionReturn(PETSC_SUCCESS);
747: } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
748: }
749: PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
750: PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
751: PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
752: if (datatype == PETSC_STRING) {
753: size_t len;
754: hid_t atype;
755: PetscCallHDF5Return(atype, H5Aget_type, (attribute));
756: PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
757: PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
758: PetscCallHDF5(H5Tset_size, (dtype, len + 1));
759: PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
760: } else {
761: PetscCallHDF5(H5Aread, (attribute, dtype, value));
762: }
763: PetscCallHDF5(H5Aclose, (attribute));
764: /* H5Oclose can be used to close groups, datasets, or committed datatypes */
765: PetscCallHDF5(H5Oclose, (obj));
766: PetscCall(PetscFree(parentAbsPath));
767: PetscFunctionReturn(PETSC_SUCCESS);
768: }
770: static PetscErrorCode PetscViewerHDF5HasGroup_HDF5(PetscViewer viewer, const char path[], PetscBool *has)
771: {
772: H5O_type_t type;
773: const char *abspath;
775: PetscFunctionBegin;
776: PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
777: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
778: *has = (PetscBool)(type == H5O_TYPE_GROUP);
779: PetscCall(PetscFree(abspath));
780: PetscFunctionReturn(PETSC_SUCCESS);
781: }
783: static PetscErrorCode PetscViewerHDF5HasDataset_HDF5(PetscViewer viewer, const char path[], PetscBool *has)
784: {
785: H5O_type_t type;
786: const char *abspath;
788: PetscFunctionBegin;
789: PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
790: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
791: *has = (PetscBool)(type == H5O_TYPE_DATASET);
792: PetscCall(PetscFree(abspath));
793: PetscFunctionReturn(PETSC_SUCCESS);
794: }
796: static PetscErrorCode PetscViewerHDF5HasAttribute_HDF5(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
797: {
798: const char *parentAbsPath;
800: PetscFunctionBegin;
801: PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
802: PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
803: if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
804: PetscCall(PetscFree(parentAbsPath));
805: PetscFunctionReturn(PETSC_SUCCESS);
806: }
808: /*MC
809: PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
811: Level: beginner
813: .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
814: `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
815: `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
816: `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
817: M*/
819: PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
820: {
821: PetscViewer_HDF5 *hdf5;
823: PetscFunctionBegin;
824: #if !defined(H5_HAVE_PARALLEL)
825: {
826: PetscMPIInt size;
827: PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
828: PetscCheck(size <= 1, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)");
829: }
830: #endif
832: PetscCall(PetscNew(&hdf5));
834: v->data = (void *)hdf5;
835: v->ops->destroy = PetscViewerDestroy_HDF5;
836: v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
837: v->ops->setup = PetscViewerSetUp_HDF5;
838: v->ops->view = PetscViewerView_HDF5;
839: v->ops->flush = PetscViewerFlush_HDF5;
840: hdf5->btype = FILE_MODE_UNDEFINED;
841: hdf5->filename = NULL;
842: hdf5->timestep = -1;
843: hdf5->groups = NULL;
845: PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
847: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
848: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
849: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
850: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
851: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
852: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetBaseDimension2_C", PetscViewerHDF5GetBaseDimension2_HDF5));
853: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
854: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
855: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
856: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
857: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
858: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCompress_C", PetscViewerHDF5GetCompress_HDF5));
859: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCompress_C", PetscViewerHDF5SetCompress_HDF5));
860: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5WriteGroup_C", PetscViewerHDF5WriteGroup_HDF5));
861: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PushGroup_C", PetscViewerHDF5PushGroup_HDF5));
862: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PopGroup_C", PetscViewerHDF5PopGroup_HDF5));
863: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetGroup_C", PetscViewerHDF5GetGroup_HDF5));
864: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PushTimestepping_C", PetscViewerHDF5PushTimestepping_HDF5));
865: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5PopTimestepping_C", PetscViewerHDF5PopTimestepping_HDF5));
866: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5IsTimestepping_C", PetscViewerHDF5IsTimestepping_HDF5));
867: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5IncrementTimestep_C", PetscViewerHDF5IncrementTimestep_HDF5));
868: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetTimestep_C", PetscViewerHDF5SetTimestep_HDF5));
869: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetTimestep_C", PetscViewerHDF5GetTimestep_HDF5));
870: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5WriteAttribute_C", PetscViewerHDF5WriteAttribute_HDF5));
871: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5ReadAttribute_C", PetscViewerHDF5ReadAttribute_HDF5));
872: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasGroup_C", PetscViewerHDF5HasGroup_HDF5));
873: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasDataset_C", PetscViewerHDF5HasDataset_HDF5));
874: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5HasAttribute_C", PetscViewerHDF5HasAttribute_HDF5));
875: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetSPOutput_C", PetscViewerHDF5GetSPOutput_HDF5));
876: PetscFunctionReturn(PETSC_SUCCESS);
877: }