xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 1b37a2a7cc4a4fb30c3e967db1c694c0a1013f51)
1 #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/
2 
3 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *);
4 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *);
5 static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer, const char *[]);
6 
7 /*@C
8   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
9 
10   Not Collective
11 
12   Input Parameters:
13 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
14 - path   - (Optional) The path relative to the pushed group
15 
16   Output Parameter:
17 . abspath - The absolute HDF5 path (group)
18 
19   Level: intermediate
20 
21   Notes:
22   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
23   `NULL` or empty path means the current pushed group.
24 
25   The output abspath is newly allocated so needs to be freed.
26 
27 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
28 @*/
29 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[])
30 {
31   size_t      len;
32   PetscBool   relative = PETSC_FALSE;
33   const char *group;
34   char        buf[PETSC_MAX_PATH_LEN] = "";
35 
36   PetscFunctionBegin;
37   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
38   if (path) PetscAssertPointer(path, 2);
39   PetscAssertPointer(abspath, 3);
40   PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
41   PetscCall(PetscStrlen(path, &len));
42   relative = (PetscBool)(!len || path[0] != '/');
43   if (relative) {
44     PetscCall(PetscStrncpy(buf, group, sizeof(buf)));
45     if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
46     PetscCall(PetscStrlcat(buf, path, sizeof(buf)));
47     PetscCall(PetscStrallocpy(buf, abspath));
48   } else {
49     PetscCall(PetscStrallocpy(path, abspath));
50   }
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53 
54 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
55 {
56   PetscBool has;
57 
58   PetscFunctionBegin;
59   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
60   if (!has) {
61     char *group;
62     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
63     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group);
64   }
65   PetscFunctionReturn(PETSC_SUCCESS);
66 }
67 
68 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject)
69 {
70   PetscBool         flg  = PETSC_FALSE, set;
71   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
72 
73   PetscFunctionBegin;
74   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
75   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
76   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
77   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
78   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
79   flg = PETSC_FALSE;
80   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
81   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
82   PetscOptionsHeadEnd();
83   PetscFunctionReturn(PETSC_SUCCESS);
84 }
85 
86 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
87 {
88   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
89   PetscBool         flg;
90 
91   PetscFunctionBegin;
92   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
93   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
94   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
95   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
96   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
97   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
98   PetscFunctionReturn(PETSC_SUCCESS);
99 }
100 
101 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
102 {
103   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
104 
105   PetscFunctionBegin;
106   PetscCall(PetscFree(hdf5->filename));
107   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
108   PetscFunctionReturn(PETSC_SUCCESS);
109 }
110 
111 static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
112 {
113   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
114 
115   PetscFunctionBegin;
116   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
117   PetscFunctionReturn(PETSC_SUCCESS);
118 }
119 
120 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
121 {
122   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
123 
124   PetscFunctionBegin;
125   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
126   PetscCall(PetscViewerFileClose_HDF5(viewer));
127   while (hdf5->groups) {
128     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
129 
130     PetscCall(PetscFree(hdf5->groups->name));
131     PetscCall(PetscFree(hdf5->groups));
132     hdf5->groups = tmp;
133   }
134   PetscCall(PetscFree(hdf5));
135   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
136   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
137   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
138   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
139   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
140   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
141   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
142   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
143   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
144   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
149 {
150   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
151 
152   PetscFunctionBegin;
153   hdf5->btype = type;
154   PetscFunctionReturn(PETSC_SUCCESS);
155 }
156 
157 static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
158 {
159   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
160 
161   PetscFunctionBegin;
162   *type = hdf5->btype;
163   PetscFunctionReturn(PETSC_SUCCESS);
164 }
165 
166 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
167 {
168   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
169 
170   PetscFunctionBegin;
171   hdf5->basedimension2 = flg;
172   PetscFunctionReturn(PETSC_SUCCESS);
173 }
174 
175 /*@
176   PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
177   dimension of 2.
178 
179   Logically Collective
180 
181   Input Parameters:
182 + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
183 - flg    - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1
184 
185   Options Database Key:
186 . -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
187 
188   Level: intermediate
189 
190   Note:
191   Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
192   of one when the dimension is lower. Others think the option is crazy.
193 
194 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
195 @*/
196 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
197 {
198   PetscFunctionBegin;
199   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
200   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
201   PetscFunctionReturn(PETSC_SUCCESS);
202 }
203 
204 /*@
205   PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
206   dimension of 2.
207 
208   Logically Collective
209 
210   Input Parameter:
211 . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`
212 
213   Output Parameter:
214 . flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1
215 
216   Level: intermediate
217 
218   Note:
219   Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
220   of one when the dimension is lower. Others think the option is crazy.
221 
222 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
223 @*/
224 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
225 {
226   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
227 
228   PetscFunctionBegin;
229   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
230   *flg = hdf5->basedimension2;
231   PetscFunctionReturn(PETSC_SUCCESS);
232 }
233 
234 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
235 {
236   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
237 
238   PetscFunctionBegin;
239   hdf5->spoutput = flg;
240   PetscFunctionReturn(PETSC_SUCCESS);
241 }
242 
243 /*@
244   PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
245   compiled with double precision `PetscReal`.
246 
247   Logically Collective
248 
249   Input Parameters:
250 + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
251 - flg    - if `PETSC_TRUE` the data will be written to disk with single precision
252 
253   Options Database Key:
254 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
255 
256   Level: intermediate
257 
258   Note:
259   Setting this option does not make any difference if PETSc is compiled with single precision
260   in the first place. It does not affect reading datasets (HDF5 handle this internally).
261 
262 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
263           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
264 @*/
265 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
266 {
267   PetscFunctionBegin;
268   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
269   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
270   PetscFunctionReturn(PETSC_SUCCESS);
271 }
272 
273 /*@
274   PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
275   compiled with double precision `PetscReal`.
276 
277   Logically Collective
278 
279   Input Parameter:
280 . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`
281 
282   Output Parameter:
283 . flg - if `PETSC_TRUE` the data will be written to disk with single precision
284 
285   Level: intermediate
286 
287 .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
288           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
289 @*/
290 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
291 {
292   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
296   *flg = hdf5->spoutput;
297   PetscFunctionReturn(PETSC_SUCCESS);
298 }
299 
300 static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
301 {
302   PetscFunctionBegin;
303   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
304      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
305 #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
306   {
307     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
308     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
309   }
310 #else
311   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"));
312 #endif
313   PetscFunctionReturn(PETSC_SUCCESS);
314 }
315 
316 /*@
317   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
318 
319   Logically Collective; flg must contain common value
320 
321   Input Parameters:
322 + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
323 - flg    - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)
324 
325   Options Database Key:
326 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
327 
328   Level: intermediate
329 
330   Note:
331   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
332   However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions.
333 
334   Developer Notes:
335   In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively.
336   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
337   See HDF5 documentation and MPI-IO documentation for details.
338 
339 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
340 @*/
341 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
342 {
343   PetscFunctionBegin;
344   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
345   PetscValidLogicalCollectiveBool(viewer, flg, 2);
346   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
347   PetscFunctionReturn(PETSC_SUCCESS);
348 }
349 
350 static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
351 {
352 #if defined(H5_HAVE_PARALLEL)
353   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
354   H5FD_mpio_xfer_t  mode;
355 #endif
356 
357   PetscFunctionBegin;
358 #if !defined(H5_HAVE_PARALLEL)
359   *flg = PETSC_FALSE;
360 #else
361   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
362   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
363 #endif
364   PetscFunctionReturn(PETSC_SUCCESS);
365 }
366 
367 /*@
368   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
369 
370   Not Collective
371 
372   Input Parameter:
373 . viewer - the `PETSCVIEWERHDF5` `PetscViewer`
374 
375   Output Parameter:
376 . flg - the flag
377 
378   Level: intermediate
379 
380   Note:
381   This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, `PETSC_FALSE` will be always returned.
382   For more details, see `PetscViewerHDF5SetCollective()`.
383 
384 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
385 @*/
386 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
387 {
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
390   PetscAssertPointer(flg, 2);
391 
392   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
393   PetscFunctionReturn(PETSC_SUCCESS);
394 }
395 
396 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
397 {
398   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
399   hid_t             plist_id;
400 
401   PetscFunctionBegin;
402   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
403   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
404   PetscCall(PetscStrallocpy(name, &hdf5->filename));
405   /* Set up file access property list with parallel I/O access */
406   PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS));
407 #if defined(H5_HAVE_PARALLEL)
408   PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
409 #endif
410   /* Create or open the file collectively */
411   switch (hdf5->btype) {
412   case FILE_MODE_READ:
413     if (PetscDefined(USE_DEBUG)) {
414       PetscMPIInt rank;
415       PetscBool   flg;
416 
417       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
418       if (rank == 0) {
419         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
420         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
421       }
422       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
423     }
424     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id));
425     break;
426   case FILE_MODE_APPEND:
427   case FILE_MODE_UPDATE: {
428     PetscBool flg;
429     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
430     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id));
431     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
432     break;
433   }
434   case FILE_MODE_WRITE:
435     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
436     break;
437   case FILE_MODE_UNDEFINED:
438     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
439   default:
440     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
441   }
442   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
443   PetscCallHDF5(H5Pclose, (plist_id));
444   PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
445   PetscFunctionReturn(PETSC_SUCCESS);
446 }
447 
448 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
449 {
450   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
451 
452   PetscFunctionBegin;
453   *name = vhdf5->filename;
454   PetscFunctionReturn(PETSC_SUCCESS);
455 }
456 
457 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
458 {
459   /*
460   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
461   */
462 
463   PetscFunctionBegin;
464   PetscFunctionReturn(PETSC_SUCCESS);
465 }
466 
467 static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
468 {
469   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
470 
471   PetscFunctionBegin;
472   hdf5->defTimestepping = flg;
473   PetscFunctionReturn(PETSC_SUCCESS);
474 }
475 
476 static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
477 {
478   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
479 
480   PetscFunctionBegin;
481   *flg = hdf5->defTimestepping;
482   PetscFunctionReturn(PETSC_SUCCESS);
483 }
484 
485 /*@
486   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
487 
488   Logically Collective
489 
490   Input Parameters:
491 + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
492 - flg    - if `PETSC_TRUE` we will assume that timestepping is on
493 
494   Options Database Key:
495 . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
496 
497   Level: intermediate
498 
499   Note:
500   If the timestepping attribute is not found for an object, then the default timestepping is used
501 
502 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
503 @*/
504 PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
505 {
506   PetscFunctionBegin;
507   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
508   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
509   PetscFunctionReturn(PETSC_SUCCESS);
510 }
511 
512 /*@
513   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
514 
515   Not Collective
516 
517   Input Parameter:
518 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
519 
520   Output Parameter:
521 . flg - if `PETSC_TRUE` we will assume that timestepping is on
522 
523   Level: intermediate
524 
525 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
526 @*/
527 PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
528 {
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
531   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
532   PetscFunctionReturn(PETSC_SUCCESS);
533 }
534 
535 /*MC
536    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
537 
538   Level: beginner
539 
540 .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
541           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
542           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
543           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
544 M*/
545 
546 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
547 {
548   PetscViewer_HDF5 *hdf5;
549 
550   PetscFunctionBegin;
551 #if !defined(H5_HAVE_PARALLEL)
552   {
553     PetscMPIInt size;
554     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
555     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)");
556   }
557 #endif
558 
559   PetscCall(PetscNew(&hdf5));
560 
561   v->data                = (void *)hdf5;
562   v->ops->destroy        = PetscViewerDestroy_HDF5;
563   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
564   v->ops->setup          = PetscViewerSetUp_HDF5;
565   v->ops->view           = PetscViewerView_HDF5;
566   v->ops->flush          = PetscViewerFlush_HDF5;
567   hdf5->btype            = FILE_MODE_UNDEFINED;
568   hdf5->filename         = NULL;
569   hdf5->timestep         = -1;
570   hdf5->groups           = NULL;
571 
572   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
573 
574   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
575   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
576   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
577   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
578   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
579   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
580   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
581   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
582   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
583   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
584   PetscFunctionReturn(PETSC_SUCCESS);
585 }
586 
587 /*@C
588   PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`
589 
590   Collective
591 
592   Input Parameters:
593 + comm - MPI communicator
594 . name - name of file
595 - type - type of file
596 
597   Output Parameter:
598 . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
599 
600   Options Database Keys:
601 + -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1
602 - -viewer_hdf5_sp_output       - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
603 
604   Level: beginner
605 
606   Notes:
607   Reading is always available, regardless of the mode. Available modes are
608 .vb
609   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
610   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
611   FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
612   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
613 .ve
614 
615   In case of `FILE_MODE_APPEND` / `FILE_MODE_UPDATE`, any stored object (dataset, attribute) can be selectively overwritten if the same fully qualified name (/group/path/to/object) is specified.
616 
617   This PetscViewer should be destroyed with PetscViewerDestroy().
618 
619 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
620           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
621           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
622 @*/
623 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
624 {
625   PetscFunctionBegin;
626   PetscCall(PetscViewerCreate(comm, hdf5v));
627   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
628   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
629   PetscCall(PetscViewerFileSetName(*hdf5v, name));
630   PetscCall(PetscViewerSetFromOptions(*hdf5v));
631   PetscFunctionReturn(PETSC_SUCCESS);
632 }
633 
634 /*@C
635   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
636 
637   Not Collective
638 
639   Input Parameter:
640 . viewer - the `PetscViewer`
641 
642   Output Parameter:
643 . file_id - The file id
644 
645   Level: intermediate
646 
647 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
648 @*/
649 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
650 {
651   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
652 
653   PetscFunctionBegin;
654   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
655   if (file_id) *file_id = hdf5->file_id;
656   PetscFunctionReturn(PETSC_SUCCESS);
657 }
658 
659 /*@C
660   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
661 
662   Not Collective
663 
664   Input Parameters:
665 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
666 - name   - The group name
667 
668   Level: intermediate
669 
670   Notes:
671   This is designed to mnemonically resemble the Unix cd command.
672 .vb
673   If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group.
674   `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
675   "." means the current group is pushed again.
676 .ve
677 
678   Example:
679   Suppose the current group is "/a".
680 .vb
681   If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
682   If name is ".", then the new group will be "/a".
683   If name is "b", then the new group will be "/a/b".
684   If name is "/b", then the new group will be "/b".
685 .ve
686 
687   Developer Notes:
688   The root group "/" is internally stored as `NULL`.
689 
690 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
691 @*/
692 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
693 {
694   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
695   PetscViewerHDF5GroupList *groupNode;
696   size_t                    i, len;
697   char                      buf[PETSC_MAX_PATH_LEN];
698   const char               *gname;
699 
700   PetscFunctionBegin;
701   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
702   if (name) PetscAssertPointer(name, 2);
703   PetscCall(PetscStrlen(name, &len));
704   gname = NULL;
705   if (len) {
706     if (len == 1 && name[0] == '.') {
707       /* use current name */
708       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
709     } else if (name[0] == '/') {
710       /* absolute */
711       for (i = 1; i < len; i++) {
712         if (name[i] != '/') {
713           gname = name;
714           break;
715         }
716       }
717     } else {
718       /* relative */
719       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
720       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
721       gname = buf;
722     }
723   }
724   PetscCall(PetscNew(&groupNode));
725   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
726   groupNode->next = hdf5->groups;
727   hdf5->groups    = groupNode;
728   PetscFunctionReturn(PETSC_SUCCESS);
729 }
730 
731 /*@
732   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
733 
734   Not Collective
735 
736   Input Parameter:
737 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
738 
739   Level: intermediate
740 
741 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
742 @*/
743 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
744 {
745   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
746   PetscViewerHDF5GroupList *groupNode;
747 
748   PetscFunctionBegin;
749   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
750   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
751   groupNode    = hdf5->groups;
752   hdf5->groups = hdf5->groups->next;
753   PetscCall(PetscFree(groupNode->name));
754   PetscCall(PetscFree(groupNode));
755   PetscFunctionReturn(PETSC_SUCCESS);
756 }
757 
758 static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
759 {
760   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
761 
762   PetscFunctionBegin;
763   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
764   PetscAssertPointer(name, 2);
765   if (hdf5->groups) *name = hdf5->groups->name;
766   else *name = NULL;
767   PetscFunctionReturn(PETSC_SUCCESS);
768 }
769 
770 /*@C
771   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
772   and return this group's ID and file ID.
773   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
774 
775   Not Collective
776 
777   Input Parameters:
778 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
779 - path   - (Optional) The path relative to the pushed group
780 
781   Output Parameters:
782 + fileId  - The HDF5 file ID
783 - groupId - The HDF5 group ID
784 
785   Level: intermediate
786 
787   Note:
788   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
789   `NULL` or empty path means the current pushed group.
790 
791   If the viewer is writable, the group is created if it doesn't exist yet.
792 
793 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
794 @*/
795 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
796 {
797   hid_t       file_id;
798   H5O_type_t  type;
799   const char *fileName  = NULL;
800   char       *groupName = NULL;
801   PetscBool   writable, has;
802 
803   PetscFunctionBegin;
804   PetscCall(PetscViewerWritable(viewer, &writable));
805   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
806   PetscCall(PetscViewerFileGetName(viewer, &fileName));
807   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
808   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
809   if (!has) {
810     PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
811     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
812   }
813   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);
814   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
815   PetscCall(PetscFree(groupName));
816   *fileId = file_id;
817   PetscFunctionReturn(PETSC_SUCCESS);
818 }
819 
820 /*@C
821   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file
822 
823   Not Collective
824 
825   Input Parameters:
826 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
827 - path   - (Optional) The path relative to the pushed group
828 
829   Level: intermediate
830 
831   Note:
832   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
833   `NULL` or empty path means the current pushed group.
834 
835   This will fail if the viewer is not writable.
836 
837 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
838 @*/
839 PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
840 {
841   hid_t fileId, groupId;
842 
843   PetscFunctionBegin;
844   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
845   PetscCallHDF5(H5Gclose, (groupId));
846   PetscFunctionReturn(PETSC_SUCCESS);
847 }
848 
849 /*@
850   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
851 
852   Not Collective
853 
854   Input Parameter:
855 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
856 
857   Level: intermediate
858 
859   Notes:
860   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
861   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
862   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. `VecView()`].
863   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
864   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
865 
866   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
867   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
868 
869   Developer Notes:
870   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
871 
872 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
873 @*/
874 PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
875 {
876   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
877 
878   PetscFunctionBegin;
879   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
880   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
881   hdf5->timestepping = PETSC_TRUE;
882   if (hdf5->timestep < 0) hdf5->timestep = 0;
883   PetscFunctionReturn(PETSC_SUCCESS);
884 }
885 
886 /*@
887   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
888 
889   Not Collective
890 
891   Input Parameter:
892 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
893 
894   Level: intermediate
895 
896   Note:
897   See `PetscViewerHDF5PushTimestepping()` for details.
898 
899 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
900 @*/
901 PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
902 {
903   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
904 
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
907   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
908   hdf5->timestepping = PETSC_FALSE;
909   PetscFunctionReturn(PETSC_SUCCESS);
910 }
911 
912 /*@
913   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
914 
915   Not Collective
916 
917   Input Parameter:
918 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
919 
920   Output Parameter:
921 . flg - is timestepping active?
922 
923   Level: intermediate
924 
925   Note:
926   See `PetscViewerHDF5PushTimestepping()` for details.
927 
928 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
929 @*/
930 PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
931 {
932   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
933 
934   PetscFunctionBegin;
935   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
936   *flg = hdf5->timestepping;
937   PetscFunctionReturn(PETSC_SUCCESS);
938 }
939 
940 /*@
941   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
942 
943   Not Collective
944 
945   Input Parameter:
946 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
947 
948   Level: intermediate
949 
950   Note:
951   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
952 
953 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
954 @*/
955 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
956 {
957   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
958 
959   PetscFunctionBegin;
960   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
961   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
962   ++hdf5->timestep;
963   PetscFunctionReturn(PETSC_SUCCESS);
964 }
965 
966 /*@
967   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
968 
969   Logically Collective
970 
971   Input Parameters:
972 + viewer   - the `PetscViewer` of type `PETSCVIEWERHDF5`
973 - timestep - The timestep
974 
975   Level: intermediate
976 
977   Note:
978   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
979 
980 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
981 @*/
982 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
983 {
984   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
985 
986   PetscFunctionBegin;
987   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
988   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
989   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
990   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
991   hdf5->timestep = timestep;
992   PetscFunctionReturn(PETSC_SUCCESS);
993 }
994 
995 /*@
996   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
997 
998   Not Collective
999 
1000   Input Parameter:
1001 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
1002 
1003   Output Parameter:
1004 . timestep - The timestep
1005 
1006   Level: intermediate
1007 
1008   Note:
1009   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1010 
1011 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
1012 @*/
1013 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
1014 {
1015   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1016 
1017   PetscFunctionBegin;
1018   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1019   PetscAssertPointer(timestep, 2);
1020   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
1021   *timestep = hdf5->timestep;
1022   PetscFunctionReturn(PETSC_SUCCESS);
1023 }
1024 
1025 /*@C
1026   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
1027 
1028   Not Collective
1029 
1030   Input Parameter:
1031 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
1032 
1033   Output Parameter:
1034 . htype - the HDF5  datatype
1035 
1036   Level: advanced
1037 
1038 .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
1039 @*/
1040 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
1041 {
1042   PetscFunctionBegin;
1043   if (ptype == PETSC_INT)
1044 #if defined(PETSC_USE_64BIT_INDICES)
1045     *htype = H5T_NATIVE_LLONG;
1046 #else
1047     *htype = H5T_NATIVE_INT;
1048 #endif
1049   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
1050   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
1051   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1052   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
1053   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
1054   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
1055   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
1056   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
1057   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
1058   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
1059   PetscFunctionReturn(PETSC_SUCCESS);
1060 }
1061 
1062 /*@C
1063   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
1064 
1065   Not Collective
1066 
1067   Input Parameter:
1068 . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
1069 
1070   Output Parameter:
1071 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
1072 
1073   Level: advanced
1074 
1075 .seealso: [](sec_viewers), `PetscDataType`
1076 @*/
1077 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
1078 {
1079   PetscFunctionBegin;
1080 #if defined(PETSC_USE_64BIT_INDICES)
1081   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
1082   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
1083 #else
1084   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
1085 #endif
1086   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
1087   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
1088   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
1089   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
1090   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
1091   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
1092   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
1093   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
1094   PetscFunctionReturn(PETSC_SUCCESS);
1095 }
1096 
1097 /*@C
1098   PetscViewerHDF5WriteAttribute - Write an attribute
1099 
1100   Collective
1101 
1102   Input Parameters:
1103 + viewer   - The `PETSCVIEWERHDF5` viewer
1104 . parent   - The parent dataset/group name
1105 . name     - The attribute name
1106 . datatype - The attribute type
1107 - value    - The attribute value
1108 
1109   Level: advanced
1110 
1111   Note:
1112   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
1113 
1114 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1115           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1116 @*/
1117 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1118 {
1119   char     *parentAbsPath;
1120   hid_t     h5, dataspace, obj, attribute, dtype;
1121   PetscBool has;
1122 
1123   PetscFunctionBegin;
1124   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1125   if (parent) PetscAssertPointer(parent, 2);
1126   PetscAssertPointer(name, 3);
1127   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1128   PetscAssertPointer(value, 5);
1129   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1130   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
1131   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1132   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1133   if (datatype == PETSC_STRING) {
1134     size_t len;
1135     PetscCall(PetscStrlen((const char *)value, &len));
1136     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1137   }
1138   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1139   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1140   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1141   if (has) {
1142     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1143   } else {
1144     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1145   }
1146   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1147   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1148   PetscCallHDF5(H5Aclose, (attribute));
1149   PetscCallHDF5(H5Oclose, (obj));
1150   PetscCallHDF5(H5Sclose, (dataspace));
1151   PetscCall(PetscFree(parentAbsPath));
1152   PetscFunctionReturn(PETSC_SUCCESS);
1153 }
1154 
1155 /*@C
1156   PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1157 
1158   Collective
1159 
1160   Input Parameters:
1161 + viewer   - The `PETSCVIEWERHDF5` viewer
1162 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1163 . name     - The attribute name
1164 . datatype - The attribute type
1165 - value    - The attribute value
1166 
1167   Level: advanced
1168 
1169   Note:
1170   This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1171   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1172 
1173 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1174           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1175 @*/
1176 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1177 {
1178   PetscFunctionBegin;
1179   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1180   PetscValidHeader(obj, 2);
1181   PetscAssertPointer(name, 3);
1182   PetscAssertPointer(value, 5);
1183   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1184   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1185   PetscFunctionReturn(PETSC_SUCCESS);
1186 }
1187 
1188 /*@C
1189   PetscViewerHDF5ReadAttribute - Read an attribute
1190 
1191   Collective
1192 
1193   Input Parameters:
1194 + viewer       - The `PETSCVIEWERHDF5` viewer
1195 . parent       - The parent dataset/group name
1196 . name         - The attribute name
1197 . datatype     - The attribute type
1198 - defaultValue - The pointer to the default value
1199 
1200   Output Parameter:
1201 . value - The pointer to the read HDF5 attribute value
1202 
1203   Level: advanced
1204 
1205   Notes:
1206   If defaultValue is `NULL` and the attribute is not found, an error occurs.
1207 
1208   If defaultValue is not `NULL` and the attribute is not found, `defaultValue` is copied to value.
1209 
1210   The pointers `defaultValue` and value can be the same; for instance
1211 .vb
1212   flg = PETSC_FALSE;
1213   PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1214 .ve
1215   is valid, but make sure the default value is initialized.
1216 
1217   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1218 
1219   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. `NULL` means the current pushed group.
1220 
1221 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1222 @*/
1223 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1224 {
1225   char     *parentAbsPath;
1226   hid_t     h5, obj, attribute, dtype;
1227   PetscBool has;
1228 
1229   PetscFunctionBegin;
1230   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1231   if (parent) PetscAssertPointer(parent, 2);
1232   PetscAssertPointer(name, 3);
1233   if (defaultValue) PetscAssertPointer(defaultValue, 5);
1234   PetscAssertPointer(value, 6);
1235   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1236   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1237   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
1238   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1239   if (!has) {
1240     if (defaultValue) {
1241       if (defaultValue != value) {
1242         if (datatype == PETSC_STRING) {
1243           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1244         } else {
1245           size_t len;
1246           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
1247           PetscCall(PetscMemcpy(value, defaultValue, len));
1248         }
1249       }
1250       PetscCall(PetscFree(parentAbsPath));
1251       PetscFunctionReturn(PETSC_SUCCESS);
1252     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1253   }
1254   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1255   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1256   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1257   if (datatype == PETSC_STRING) {
1258     size_t len;
1259     hid_t  atype;
1260     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1261     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
1262     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1263     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1264     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1265   } else {
1266     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1267   }
1268   PetscCallHDF5(H5Aclose, (attribute));
1269   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1270   PetscCallHDF5(H5Oclose, (obj));
1271   PetscCall(PetscFree(parentAbsPath));
1272   PetscFunctionReturn(PETSC_SUCCESS);
1273 }
1274 
1275 /*@C
1276   PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1277 
1278   Collective
1279 
1280   Input Parameters:
1281 + viewer       - The `PETSCVIEWERHDF5` viewer
1282 . obj          - The object whose name is used to lookup the parent dataset, relative to the current group.
1283 . name         - The attribute name
1284 . datatype     - The attribute type
1285 - defaultValue - The default attribute value
1286 
1287   Output Parameter:
1288 . value - The attribute value
1289 
1290   Level: advanced
1291 
1292   Note:
1293   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1294   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1295 
1296 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1297 @*/
1298 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1299 {
1300   PetscFunctionBegin;
1301   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1302   PetscValidHeader(obj, 2);
1303   PetscAssertPointer(name, 3);
1304   PetscAssertPointer(value, 6);
1305   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1306   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1307   PetscFunctionReturn(PETSC_SUCCESS);
1308 }
1309 
1310 static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1311 {
1312   htri_t exists;
1313   hid_t  group;
1314 
1315   PetscFunctionBegin;
1316   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1317   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1318   if (!exists && createGroup) {
1319     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1320     PetscCallHDF5(H5Gclose, (group));
1321     exists = PETSC_TRUE;
1322   }
1323   *exists_ = (PetscBool)exists;
1324   PetscFunctionReturn(PETSC_SUCCESS);
1325 }
1326 
1327 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1328 {
1329   const char rootGroupName[] = "/";
1330   hid_t      h5;
1331   PetscBool  exists = PETSC_FALSE;
1332   PetscInt   i;
1333   int        n;
1334   char     **hierarchy;
1335   char       buf[PETSC_MAX_PATH_LEN] = "";
1336 
1337   PetscFunctionBegin;
1338   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1339   if (name) PetscAssertPointer(name, 2);
1340   else name = rootGroupName;
1341   if (has) {
1342     PetscAssertPointer(has, 4);
1343     *has = PETSC_FALSE;
1344   }
1345   if (otype) {
1346     PetscAssertPointer(otype, 5);
1347     *otype = H5O_TYPE_UNKNOWN;
1348   }
1349   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1350 
1351   /*
1352      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1353      Hence, each of them needs to be tested separately:
1354      1) whether it's a valid link
1355      2) whether this link resolves to an object
1356      See H5Oexists_by_name() documentation.
1357   */
1358   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
1359   if (!n) {
1360     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1361     if (has) *has = PETSC_TRUE;
1362     if (otype) *otype = H5O_TYPE_GROUP;
1363     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
1364     PetscFunctionReturn(PETSC_SUCCESS);
1365   }
1366   for (i = 0; i < n; i++) {
1367     PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
1368     PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf)));
1369     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1370     if (!exists) break;
1371   }
1372   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
1373 
1374   /* If the object exists, get its type */
1375   if (exists && otype) {
1376     H5O_info_t info;
1377 
1378     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1379     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
1380     *otype = info.type;
1381   }
1382   if (has) *has = exists;
1383   PetscFunctionReturn(PETSC_SUCCESS);
1384 }
1385 
1386 /*@C
1387   PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1388 
1389   Collective
1390 
1391   Input Parameters:
1392 + viewer - The `PETSCVIEWERHDF5` viewer
1393 - path   - (Optional) The path relative to the pushed group
1394 
1395   Output Parameter:
1396 . has - Flag for group existence
1397 
1398   Level: advanced
1399 
1400   Notes:
1401   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1402   `NULL` or empty path means the current pushed group.
1403 
1404   If path exists but is not a group, `PETSC_FALSE` is returned.
1405 
1406 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
1407 @*/
1408 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1409 {
1410   H5O_type_t type;
1411   char      *abspath;
1412 
1413   PetscFunctionBegin;
1414   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1415   if (path) PetscAssertPointer(path, 2);
1416   PetscAssertPointer(has, 3);
1417   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
1418   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1419   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1420   PetscCall(PetscFree(abspath));
1421   PetscFunctionReturn(PETSC_SUCCESS);
1422 }
1423 
1424 /*@C
1425   PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
1426 
1427   Collective
1428 
1429   Input Parameters:
1430 + viewer - The `PETSCVIEWERHDF5` viewer
1431 - path   - The dataset path
1432 
1433   Output Parameter:
1434 . has - Flag whether dataset exists
1435 
1436   Level: advanced
1437 
1438   Notes:
1439   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1440 
1441   If `path` is `NULL` or empty, has is set to `PETSC_FALSE`.
1442 
1443   If `path` exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1444 
1445 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1446 @*/
1447 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1448 {
1449   H5O_type_t type;
1450   char      *abspath;
1451 
1452   PetscFunctionBegin;
1453   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1454   if (path) PetscAssertPointer(path, 2);
1455   PetscAssertPointer(has, 3);
1456   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
1457   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1458   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1459   PetscCall(PetscFree(abspath));
1460   PetscFunctionReturn(PETSC_SUCCESS);
1461 }
1462 
1463 /*@
1464   PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1465 
1466   Collective
1467 
1468   Input Parameters:
1469 + viewer - The `PETSCVIEWERHDF5` viewer
1470 - obj    - The named object
1471 
1472   Output Parameter:
1473 . has - Flag for dataset existence
1474 
1475   Level: advanced
1476 
1477   Notes:
1478   If the object is unnamed, an error occurs.
1479 
1480   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1481 
1482 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1483 @*/
1484 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1485 {
1486   size_t len;
1487 
1488   PetscFunctionBegin;
1489   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1490   PetscValidHeader(obj, 2);
1491   PetscAssertPointer(has, 3);
1492   PetscCall(PetscStrlen(obj->name, &len));
1493   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1494   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1495   PetscFunctionReturn(PETSC_SUCCESS);
1496 }
1497 
1498 /*@C
1499   PetscViewerHDF5HasAttribute - Check whether an attribute exists
1500 
1501   Collective
1502 
1503   Input Parameters:
1504 + viewer - The `PETSCVIEWERHDF5` viewer
1505 . parent - The parent dataset/group name
1506 - name   - The attribute name
1507 
1508   Output Parameter:
1509 . has - Flag for attribute existence
1510 
1511   Level: advanced
1512 
1513   Note:
1514   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. `NULL` means the current pushed group.
1515 
1516 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1517 @*/
1518 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1519 {
1520   char *parentAbsPath;
1521 
1522   PetscFunctionBegin;
1523   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1524   if (parent) PetscAssertPointer(parent, 2);
1525   PetscAssertPointer(name, 3);
1526   PetscAssertPointer(has, 4);
1527   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
1528   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
1529   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
1530   PetscCall(PetscFree(parentAbsPath));
1531   PetscFunctionReturn(PETSC_SUCCESS);
1532 }
1533 
1534 /*@C
1535   PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1536 
1537   Collective
1538 
1539   Input Parameters:
1540 + viewer - The `PETSCVIEWERHDF5` viewer
1541 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1542 - name   - The attribute name
1543 
1544   Output Parameter:
1545 . has - Flag for attribute existence
1546 
1547   Level: advanced
1548 
1549   Note:
1550   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1551   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1552 
1553 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1554 @*/
1555 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1556 {
1557   PetscFunctionBegin;
1558   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1559   PetscValidHeader(obj, 2);
1560   PetscAssertPointer(name, 3);
1561   PetscAssertPointer(has, 4);
1562   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1563   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1564   PetscFunctionReturn(PETSC_SUCCESS);
1565 }
1566 
1567 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1568 {
1569   hid_t  h5;
1570   htri_t hhas;
1571 
1572   PetscFunctionBegin;
1573   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1574   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
1575   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1576   PetscFunctionReturn(PETSC_SUCCESS);
1577 }
1578 
1579 /*
1580   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1581   is attached to a communicator, in this case the attribute is a PetscViewer.
1582 */
1583 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1584 
1585 /*@C
1586   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1587 
1588   Collective
1589 
1590   Input Parameter:
1591 . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer`
1592 
1593   Options Database Key:
1594 . -viewer_hdf5_filename <name> - name of the HDF5 file
1595 
1596   Environmental variable:
1597 . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1598 
1599   Level: intermediate
1600 
1601   Note:
1602   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1603   an error code.  The HDF5 `PetscViewer` is usually used in the form
1604 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1605 
1606 .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1607 @*/
1608 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1609 {
1610   PetscErrorCode ierr;
1611   PetscMPIInt    mpi_ierr;
1612   PetscBool      flg;
1613   PetscViewer    viewer;
1614   char           fname[PETSC_MAX_PATH_LEN];
1615   MPI_Comm       ncomm;
1616 
1617   PetscFunctionBegin;
1618   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
1619   if (ierr) {
1620     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1621     PetscFunctionReturn(NULL);
1622   }
1623   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1624     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
1625     if (mpi_ierr) {
1626       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1627       PetscFunctionReturn(NULL);
1628     }
1629   }
1630   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
1631   if (mpi_ierr) {
1632     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1633     PetscFunctionReturn(NULL);
1634   }
1635   if (!flg) { /* PetscViewer not yet created */
1636     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
1637     if (ierr) {
1638       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1639       PetscFunctionReturn(NULL);
1640     }
1641     if (!flg) {
1642       ierr = PetscStrncpy(fname, "output.h5", sizeof(fname));
1643       if (ierr) {
1644         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1645         PetscFunctionReturn(NULL);
1646       }
1647     }
1648     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
1649     if (ierr) {
1650       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1651       PetscFunctionReturn(NULL);
1652     }
1653     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1654     if (ierr) {
1655       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1656       PetscFunctionReturn(NULL);
1657     }
1658     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
1659     if (mpi_ierr) {
1660       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1661       PetscFunctionReturn(NULL);
1662     }
1663   }
1664   ierr = PetscCommDestroy(&ncomm);
1665   if (ierr) {
1666     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1667     PetscFunctionReturn(NULL);
1668   }
1669   PetscFunctionReturn(viewer);
1670 }
1671