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