xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision c43546b0ae3b8d6e5c6ba2880b5f8967dfd80092)
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 /*@
667   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
668 
669   Not collective
670 
671   Input Parameter:
672 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
673 
674   Level: intermediate
675 
676 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
677 @*/
678 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) {
679   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
680   PetscViewerHDF5GroupList *groupNode;
681 
682   PetscFunctionBegin;
683   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
684   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
685   groupNode    = hdf5->groups;
686   hdf5->groups = hdf5->groups->next;
687   PetscCall(PetscFree(groupNode->name));
688   PetscCall(PetscFree(groupNode));
689   PetscFunctionReturn(0);
690 }
691 
692 /*@C
693   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
694   If none has been assigned, returns NULL.
695 
696   Not collective
697 
698   Input Parameter:
699 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
700 
701   Output Parameter:
702 . name - The group name
703 
704   Level: intermediate
705 
706 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`
707 @*/
708 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) {
709   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
710 
711   PetscFunctionBegin;
712   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
713   PetscValidPointer(name, 2);
714   if (hdf5->groups) *name = hdf5->groups->name;
715   else *name = NULL;
716   PetscFunctionReturn(0);
717 }
718 
719 /*@
720   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
721   and return this group's ID and file ID.
722   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
723 
724   Not collective
725 
726   Input Parameter:
727 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
728 
729   Output Parameters:
730 + fileId - The HDF5 file ID
731 - groupId - The HDF5 group ID
732 
733   Note:
734   If the viewer is writable, the group is created if it doesn't exist yet.
735 
736   Level: intermediate
737 
738 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
739 @*/
740 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) {
741   hid_t       file_id;
742   H5O_type_t  type;
743   const char *groupName = NULL, *fileName = NULL;
744   PetscBool   writable, has;
745 
746   PetscFunctionBegin;
747   PetscCall(PetscViewerWritable(viewer, &writable));
748   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
749   PetscCall(PetscViewerFileGetName(viewer, &fileName));
750   PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName));
751   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
752   if (!has) {
753     PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
754     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
755   }
756   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);
757   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName ? groupName : "/", H5P_DEFAULT));
758   *fileId = file_id;
759   PetscFunctionReturn(0);
760 }
761 
762 /*@
763   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
764 
765   Not collective
766 
767   Input Parameter:
768 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
769 
770   Level: intermediate
771 
772   Notes:
773   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
774   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
775   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. `VecView()`].
776   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
777   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
778 
779   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
780   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
781 
782   Developer note:
783   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
784 
785 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
786 @*/
787 PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) {
788   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
789 
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
792   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
793   hdf5->timestepping = PETSC_TRUE;
794   if (hdf5->timestep < 0) hdf5->timestep = 0;
795   PetscFunctionReturn(0);
796 }
797 
798 /*@
799   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
800 
801   Not collective
802 
803   Input Parameter:
804 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
805 
806   Level: intermediate
807 
808   Note:
809   See `PetscViewerHDF5PushTimestepping()` for details.
810 
811 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
812 @*/
813 PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) {
814   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
815 
816   PetscFunctionBegin;
817   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
818   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
819   hdf5->timestepping = PETSC_FALSE;
820   PetscFunctionReturn(0);
821 }
822 
823 /*@
824   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
825 
826   Not collective
827 
828   Input Parameter:
829 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
830 
831   Output Parameter:
832 . flg - is timestepping active?
833 
834   Level: intermediate
835 
836   Note:
837   See `PetscViewerHDF5PushTimestepping()` for details.
838 
839 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
840 @*/
841 PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) {
842   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
843 
844   PetscFunctionBegin;
845   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
846   *flg = hdf5->timestepping;
847   PetscFunctionReturn(0);
848 }
849 
850 /*@
851   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
852 
853   Not collective
854 
855   Input Parameter:
856 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
857 
858   Level: intermediate
859 
860   Note:
861   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
862 
863 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
864 @*/
865 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) {
866   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
867 
868   PetscFunctionBegin;
869   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
870   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
871   ++hdf5->timestep;
872   PetscFunctionReturn(0);
873 }
874 
875 /*@
876   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
877 
878   Logically collective
879 
880   Input Parameters:
881 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
882 - timestep - The timestep
883 
884   Level: intermediate
885 
886   Note:
887   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
888 
889 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
890 @*/
891 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) {
892   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
893 
894   PetscFunctionBegin;
895   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
896   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
897   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
898   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
899   hdf5->timestep = timestep;
900   PetscFunctionReturn(0);
901 }
902 
903 /*@
904   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
905 
906   Not collective
907 
908   Input Parameter:
909 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
910 
911   Output Parameter:
912 . timestep - The timestep
913 
914   Level: intermediate
915 
916   Note:
917   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
918 
919 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
920 @*/
921 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) {
922   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
923 
924   PetscFunctionBegin;
925   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
926   PetscValidIntPointer(timestep, 2);
927   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
928   *timestep = hdf5->timestep;
929   PetscFunctionReturn(0);
930 }
931 
932 /*@C
933   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
934 
935   Not collective
936 
937   Input Parameter:
938 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
939 
940   Output Parameter:
941 . mtype - the HDF5  datatype
942 
943   Level: advanced
944 
945 .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
946 @*/
947 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) {
948   PetscFunctionBegin;
949   if (ptype == PETSC_INT)
950 #if defined(PETSC_USE_64BIT_INDICES)
951     *htype = H5T_NATIVE_LLONG;
952 #else
953     *htype = H5T_NATIVE_INT;
954 #endif
955   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
956   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
957   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
958   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
959   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
960   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
961   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
962   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
963   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
964   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
965   PetscFunctionReturn(0);
966 }
967 
968 /*@C
969   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
970 
971   Not collective
972 
973   Input Parameter:
974 . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
975 
976   Output Parameter:
977 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
978 
979   Level: advanced
980 
981 .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
982 @*/
983 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) {
984   PetscFunctionBegin;
985 #if defined(PETSC_USE_64BIT_INDICES)
986   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
987   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
988 #else
989   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
990 #endif
991   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
992   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
993   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
994   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
995   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
996   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
997   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
998   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
999   PetscFunctionReturn(0);
1000 }
1001 
1002 /*@C
1003  PetscViewerHDF5WriteAttribute - Write an attribute
1004 
1005   Collective
1006 
1007   Input Parameters:
1008 + viewer - The `PETSCVIEWERHDF5` viewer
1009 . parent - The parent dataset/group name
1010 . name   - The attribute name
1011 . datatype - The attribute type
1012 - value    - The attribute value
1013 
1014   Level: advanced
1015 
1016   Note:
1017   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.
1018 
1019 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1020           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1021 @*/
1022 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) {
1023   char     *parentAbsPath;
1024   hid_t     h5, dataspace, obj, attribute, dtype;
1025   PetscBool has;
1026 
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1029   if (parent) PetscValidCharPointer(parent, 2);
1030   PetscValidCharPointer(name, 3);
1031   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1032   PetscValidPointer(value, 5);
1033   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
1034   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
1035   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1036   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1037   if (datatype == PETSC_STRING) {
1038     size_t len;
1039     PetscCall(PetscStrlen((const char *)value, &len));
1040     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1041   }
1042   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1043   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1044   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1045   if (has) {
1046     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1047   } else {
1048     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1049   }
1050   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1051   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1052   PetscCallHDF5(H5Aclose, (attribute));
1053   PetscCallHDF5(H5Oclose, (obj));
1054   PetscCallHDF5(H5Sclose, (dataspace));
1055   PetscCall(PetscFree(parentAbsPath));
1056   PetscFunctionReturn(0);
1057 }
1058 
1059 /*@C
1060  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1061 
1062   Collective
1063 
1064   Input Parameters:
1065 + viewer   - The `PETSCVIEWERHDF5` viewer
1066 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1067 . name     - The attribute name
1068 . datatype - The attribute type
1069 - value    - The attribute value
1070 
1071   Note:
1072   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).
1073   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1074 
1075   Level: advanced
1076 
1077 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1078           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1079 @*/
1080 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) {
1081   PetscFunctionBegin;
1082   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1083   PetscValidHeader(obj, 2);
1084   PetscValidCharPointer(name, 3);
1085   PetscValidPointer(value, 5);
1086   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1087   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 /*@C
1092  PetscViewerHDF5ReadAttribute - Read an attribute
1093 
1094   Collective
1095 
1096   Input Parameters:
1097 + viewer - The `PETSCVIEWERHDF5` viewer
1098 . parent - The parent dataset/group name
1099 . name   - The attribute name
1100 . datatype - The attribute type
1101 - defaultValue - The pointer to the default value
1102 
1103   Output Parameter:
1104 . value    - The pointer to the read HDF5 attribute value
1105 
1106   Notes:
1107   If defaultValue is NULL and the attribute is not found, an error occurs.
1108   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1109   The pointers defaultValue and value can be the same; for instance
1110 $  flg = `PETSC_FALSE`;
1111 $  PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1112   is valid, but make sure the default value is initialized.
1113 
1114   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1115 
1116   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.
1117 
1118   Level: advanced
1119 
1120 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1121 @*/
1122 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) {
1123   char     *parentAbsPath;
1124   hid_t     h5, obj, attribute, dtype;
1125   PetscBool has;
1126 
1127   PetscFunctionBegin;
1128   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1129   if (parent) PetscValidCharPointer(parent, 2);
1130   PetscValidCharPointer(name, 3);
1131   if (defaultValue) PetscValidPointer(defaultValue, 5);
1132   PetscValidPointer(value, 6);
1133   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
1134   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
1135   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
1136   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1137   if (!has) {
1138     if (defaultValue) {
1139       if (defaultValue != value) {
1140         if (datatype == PETSC_STRING) {
1141           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1142         } else {
1143           size_t len;
1144           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
1145           PetscCall(PetscMemcpy(value, defaultValue, len));
1146         }
1147       }
1148       PetscCall(PetscFree(parentAbsPath));
1149       PetscFunctionReturn(0);
1150     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1151   }
1152   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1153   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1154   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1155   if (datatype == PETSC_STRING) {
1156     size_t len;
1157     hid_t  atype;
1158     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1159     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
1160     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1161     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1162     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1163   } else {
1164     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1165   }
1166   PetscCallHDF5(H5Aclose, (attribute));
1167   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1168   PetscCallHDF5(H5Oclose, (obj));
1169   PetscCall(PetscFree(parentAbsPath));
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 /*@C
1174  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1175 
1176   Collective
1177 
1178   Input Parameters:
1179 + viewer   - The `PETSCVIEWERHDF5` viewer
1180 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1181 . name     - The attribute name
1182 - datatype - The attribute type
1183 
1184   Output Parameter:
1185 . value    - The attribute value
1186 
1187   Note:
1188   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1189   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1190 
1191   Level: advanced
1192 
1193 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1194 @*/
1195 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) {
1196   PetscFunctionBegin;
1197   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1198   PetscValidHeader(obj, 2);
1199   PetscValidCharPointer(name, 3);
1200   PetscValidPointer(value, 6);
1201   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1202   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) {
1207   htri_t exists;
1208   hid_t  group;
1209 
1210   PetscFunctionBegin;
1211   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1212   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1213   if (!exists && createGroup) {
1214     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1215     PetscCallHDF5(H5Gclose, (group));
1216     exists = PETSC_TRUE;
1217   }
1218   *exists_ = (PetscBool)exists;
1219   PetscFunctionReturn(0);
1220 }
1221 
1222 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) {
1223   const char rootGroupName[] = "/";
1224   hid_t      h5;
1225   PetscBool  exists = PETSC_FALSE;
1226   PetscInt   i;
1227   int        n;
1228   char     **hierarchy;
1229   char       buf[PETSC_MAX_PATH_LEN] = "";
1230 
1231   PetscFunctionBegin;
1232   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1233   if (name) PetscValidCharPointer(name, 2);
1234   else name = rootGroupName;
1235   if (has) {
1236     PetscValidBoolPointer(has, 4);
1237     *has = PETSC_FALSE;
1238   }
1239   if (otype) {
1240     PetscValidIntPointer(otype, 5);
1241     *otype = H5O_TYPE_UNKNOWN;
1242   }
1243   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1244 
1245   /*
1246      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1247      Hence, each of them needs to be tested separately:
1248      1) whether it's a valid link
1249      2) whether this link resolves to an object
1250      See H5Oexists_by_name() documentation.
1251   */
1252   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
1253   if (!n) {
1254     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1255     if (has) *has = PETSC_TRUE;
1256     if (otype) *otype = H5O_TYPE_GROUP;
1257     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
1258     PetscFunctionReturn(0);
1259   }
1260   for (i = 0; i < n; i++) {
1261     PetscCall(PetscStrcat(buf, "/"));
1262     PetscCall(PetscStrcat(buf, hierarchy[i]));
1263     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1264     if (!exists) break;
1265   }
1266   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
1267 
1268   /* If the object exists, get its type */
1269   if (exists && otype) {
1270     H5O_info_t info;
1271 
1272     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1273     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
1274     *otype = info.type;
1275   }
1276   if (has) *has = exists;
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 /*@C
1281  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1282 
1283   Collective
1284 
1285   Input Parameters:
1286 + viewer - The `PETSCVIEWERHDF5` viewer
1287 - path - The group path
1288 
1289   Output Parameter:
1290 . has - Flag for group existence
1291 
1292   Level: advanced
1293 
1294   Notes:
1295   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1296   So NULL or empty path means the current pushed group.
1297 
1298   If path exists but is not a group, `PETSC_FALSE` is returned.
1299 
1300 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
1301 @*/
1302 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) {
1303   H5O_type_t type;
1304   char      *abspath;
1305 
1306   PetscFunctionBegin;
1307   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1308   if (path) PetscValidCharPointer(path, 2);
1309   PetscValidBoolPointer(has, 3);
1310   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
1311   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1312   *has = (PetscBool)(type == H5O_TYPE_GROUP);
1313   PetscCall(PetscFree(abspath));
1314   PetscFunctionReturn(0);
1315 }
1316 
1317 /*@C
1318  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
1319 
1320   Collective
1321 
1322   Input Parameters:
1323 + viewer - The `PETSCVIEWERHDF5` viewer
1324 - path - The dataset path
1325 
1326   Output Parameter:
1327 . has - Flag whether dataset exists
1328 
1329   Level: advanced
1330 
1331   Notes:
1332   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1333 
1334   If path is NULL or empty, has is set to `PETSC_FALSE`.
1335 
1336   If path exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1337 
1338 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1339 @*/
1340 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) {
1341   H5O_type_t type;
1342   char      *abspath;
1343 
1344   PetscFunctionBegin;
1345   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1346   if (path) PetscValidCharPointer(path, 2);
1347   PetscValidBoolPointer(has, 3);
1348   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
1349   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
1350   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1351   PetscCall(PetscFree(abspath));
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 /*@
1356  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1357 
1358   Collective
1359 
1360   Input Parameters:
1361 + viewer - The `PETSCVIEWERHDF5` viewer
1362 - obj    - The named object
1363 
1364   Output Parameter:
1365 . has    - Flag for dataset existence
1366 
1367   Notes:
1368   If the object is unnamed, an error occurs.
1369 
1370   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1371 
1372   Level: advanced
1373 
1374 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1375 @*/
1376 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) {
1377   size_t len;
1378 
1379   PetscFunctionBegin;
1380   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1381   PetscValidHeader(obj, 2);
1382   PetscValidBoolPointer(has, 3);
1383   PetscCall(PetscStrlen(obj->name, &len));
1384   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1385   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1386   PetscFunctionReturn(0);
1387 }
1388 
1389 /*@C
1390  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1391 
1392   Collective
1393 
1394   Input Parameters:
1395 + viewer - The `PETSCVIEWERHDF5` viewer
1396 . parent - The parent dataset/group name
1397 - name   - The attribute name
1398 
1399   Output Parameter:
1400 . has    - Flag for attribute existence
1401 
1402   Level: advanced
1403 
1404   Note:
1405   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.
1406 
1407 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1408 @*/
1409 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
1410   char *parentAbsPath;
1411 
1412   PetscFunctionBegin;
1413   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1414   if (parent) PetscValidCharPointer(parent, 2);
1415   PetscValidCharPointer(name, 3);
1416   PetscValidBoolPointer(has, 4);
1417   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
1418   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
1419   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
1420   PetscCall(PetscFree(parentAbsPath));
1421   PetscFunctionReturn(0);
1422 }
1423 
1424 /*@C
1425  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1426 
1427   Collective
1428 
1429   Input Parameters:
1430 + viewer - The `PETSCVIEWERHDF5` viewer
1431 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1432 - name   - The attribute name
1433 
1434   Output Parameter:
1435 . has    - Flag for attribute existence
1436 
1437   Note:
1438   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1439   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1440 
1441   Level: advanced
1442 
1443 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1444 @*/
1445 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) {
1446   PetscFunctionBegin;
1447   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1448   PetscValidHeader(obj, 2);
1449   PetscValidCharPointer(name, 3);
1450   PetscValidBoolPointer(has, 4);
1451   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
1452   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1453   PetscFunctionReturn(0);
1454 }
1455 
1456 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
1457   hid_t  h5;
1458   htri_t hhas;
1459 
1460   PetscFunctionBegin;
1461   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1462   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
1463   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1464   PetscFunctionReturn(0);
1465 }
1466 
1467 /*
1468   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1469   is attached to a communicator, in this case the attribute is a PetscViewer.
1470 */
1471 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1472 
1473 /*@C
1474   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1475 
1476   Collective
1477 
1478   Input Parameter:
1479 . comm - the MPI communicator to share the HDF5 `PetscViewer`
1480 
1481   Level: intermediate
1482 
1483   Options Database Key:
1484 . -viewer_hdf5_filename <name> - name of the HDF5 file
1485 
1486   Environmental variable:
1487 . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1488 
1489   Note:
1490   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1491   an error code.  The HDF5 `PetscViewer` is usually used in the form
1492 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1493 
1494 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1495 @*/
1496 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) {
1497   PetscErrorCode ierr;
1498   PetscBool      flg;
1499   PetscViewer    viewer;
1500   char           fname[PETSC_MAX_PATH_LEN];
1501   MPI_Comm       ncomm;
1502 
1503   PetscFunctionBegin;
1504   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
1505   if (ierr) {
1506     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1507     PetscFunctionReturn(NULL);
1508   }
1509   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1510     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
1511     if (ierr) {
1512       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1513       PetscFunctionReturn(NULL);
1514     }
1515   }
1516   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
1517   if (ierr) {
1518     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1519     PetscFunctionReturn(NULL);
1520   }
1521   if (!flg) { /* PetscViewer not yet created */
1522     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
1523     if (ierr) {
1524       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1525       PetscFunctionReturn(NULL);
1526     }
1527     if (!flg) {
1528       ierr = PetscStrcpy(fname, "output.h5");
1529       if (ierr) {
1530         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1531         PetscFunctionReturn(NULL);
1532       }
1533     }
1534     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
1535     if (ierr) {
1536       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1537       PetscFunctionReturn(NULL);
1538     }
1539     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1540     if (ierr) {
1541       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1542       PetscFunctionReturn(NULL);
1543     }
1544     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
1545     if (ierr) {
1546       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
1547       PetscFunctionReturn(NULL);
1548     }
1549   }
1550   ierr = PetscCommDestroy(&ncomm);
1551   if (ierr) {
1552     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
1553     PetscFunctionReturn(NULL);
1554   }
1555   PetscFunctionReturn(viewer);
1556 }
1557