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