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