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