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