xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision d28bfa55376163a7647dcf4dd86a4618ebc392bc)
1 #include <petsc/private/viewerimpl.h>
2 #include <petsc/private/viewerhdf5impl.h>
3 #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
4 
5 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
6 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
7 
8 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
9 {
10   const char *group;
11   char buf[PETSC_MAX_PATH_LEN]="";
12   PetscErrorCode ierr;
13 
14   PetscFunctionBegin;
15   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
16   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
17   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
18   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
19   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
20   PetscFunctionReturn(0);
21 }
22 
23 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
24 {
25   PetscBool has;
26   const char *group;
27   PetscErrorCode ierr;
28 
29   PetscFunctionBegin;
30   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
31   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
32   if (!has) {
33     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
34     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
35   }
36   PetscFunctionReturn(0);
37 }
38 
39 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
40 {
41   PetscErrorCode   ierr;
42   PetscBool        flg = PETSC_FALSE, set;
43   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
44 
45   PetscFunctionBegin;
46   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
47   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
48   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
49   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
50   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
51   ierr = PetscOptionsTail();CHKERRQ(ierr);
52   PetscFunctionReturn(0);
53 }
54 
55 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
56 {
57   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
58   PetscErrorCode    ierr;
59 
60   PetscFunctionBegin;
61   if (hdf5->filename) {
62     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
63   }
64   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
65   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
66   PetscFunctionReturn(0);
67 }
68 
69 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
70 {
71   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
72   PetscErrorCode   ierr;
73 
74   PetscFunctionBegin;
75   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
76   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
81 {
82   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
83   PetscErrorCode   ierr;
84 
85   PetscFunctionBegin;
86   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
87   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
88   while (hdf5->groups) {
89     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
90 
91     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
92     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
93     hdf5->groups = tmp;
94   }
95   ierr = PetscFree(hdf5);CHKERRQ(ierr);
96   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
97   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
98   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
99   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
100   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
101   PetscFunctionReturn(0);
102 }
103 
104 static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
105 {
106   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
107 
108   PetscFunctionBegin;
109   hdf5->btype = type;
110   PetscFunctionReturn(0);
111 }
112 
113 static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
114 {
115   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
116 
117   PetscFunctionBegin;
118   *type = hdf5->btype;
119   PetscFunctionReturn(0);
120 }
121 
122 static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
123 {
124   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
125 
126   PetscFunctionBegin;
127   hdf5->basedimension2 = flg;
128   PetscFunctionReturn(0);
129 }
130 
131 /*@
132      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
133        dimension of 2.
134 
135     Logically Collective on PetscViewer
136 
137   Input Parameters:
138 +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
139 -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
140 
141   Options Database:
142 .  -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
143 
144 
145   Notes:
146     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
147          of one when the dimension is lower. Others think the option is crazy.
148 
149   Level: intermediate
150 
151 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
152 
153 @*/
154 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
155 {
156   PetscErrorCode ierr;
157 
158   PetscFunctionBegin;
159   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
160   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
161   PetscFunctionReturn(0);
162 }
163 
164 /*@
165      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
166        dimension of 2.
167 
168     Logically Collective on PetscViewer
169 
170   Input Parameter:
171 .  viewer - the PetscViewer, must be of type HDF5
172 
173   Output Parameter:
174 .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
175 
176   Notes:
177     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
178          of one when the dimension is lower. Others think the option is crazy.
179 
180   Level: intermediate
181 
182 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
183 
184 @*/
185 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
186 {
187   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
188 
189   PetscFunctionBegin;
190   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
191   *flg = hdf5->basedimension2;
192   PetscFunctionReturn(0);
193 }
194 
195 static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
196 {
197   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
198 
199   PetscFunctionBegin;
200   hdf5->spoutput = flg;
201   PetscFunctionReturn(0);
202 }
203 
204 /*@
205      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
206        compiled with double precision PetscReal.
207 
208     Logically Collective on PetscViewer
209 
210   Input Parameters:
211 +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
212 -  flg - if PETSC_TRUE the data will be written to disk with single precision
213 
214   Options Database:
215 .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
216 
217 
218   Notes:
219     Setting this option does not make any difference if PETSc is compiled with single precision
220          in the first place. It does not affect reading datasets (HDF5 handle this internally).
221 
222   Level: intermediate
223 
224 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
225           PetscReal
226 
227 @*/
228 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
229 {
230   PetscErrorCode ierr;
231 
232   PetscFunctionBegin;
233   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
234   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
235   PetscFunctionReturn(0);
236 }
237 
238 /*@
239      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
240        compiled with double precision PetscReal.
241 
242     Logically Collective on PetscViewer
243 
244   Input Parameter:
245 .  viewer - the PetscViewer, must be of type HDF5
246 
247   Output Parameter:
248 .  flg - if PETSC_TRUE the data will be written to disk with single precision
249 
250   Notes:
251     Setting this option does not make any difference if PETSc is compiled with single precision
252          in the first place. It does not affect reading datasets (HDF5 handle this internally).
253 
254   Level: intermediate
255 
256 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
257           PetscReal
258 
259 @*/
260 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
261 {
262   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
263 
264   PetscFunctionBegin;
265   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
266   *flg = hdf5->spoutput;
267   PetscFunctionReturn(0);
268 }
269 
270 static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
271 {
272   PetscFunctionBegin;
273   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
274      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
275 #if H5_VERSION_GE(1,10,3)
276   {
277     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
278     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
279   }
280 #else
281   if (flg) {
282     PetscErrorCode ierr;
283     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3\n");CHKERRQ(ierr);
284   }
285 #endif
286   PetscFunctionReturn(0);
287 }
288 
289 /*@
290   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
291 
292   Logically Collective; flg must contain common value
293 
294   Input Parameters:
295 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
296 - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
297 
298   Options Database:
299 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
300 
301   Notes:
302   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
303   However, this works correctly only since HDF5 1.10.3; hence, we ignore this setting for older versions.
304 
305   Developer notes:
306   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
307   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
308   See HDF5 documentation and MPI-IO documentation for details.
309 
310   Level: intermediate
311 
312 .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
313 
314 @*/
315 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
316 {
317   PetscErrorCode ierr;
318 
319   PetscFunctionBegin;
320   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
321   PetscValidLogicalCollectiveBool(viewer,flg,2);
322   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
323   PetscFunctionReturn(0);
324 }
325 
326 static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
327 {
328   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
329   H5FD_mpio_xfer_t  mode;
330 
331   PetscFunctionBegin;
332   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
333   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
334   PetscFunctionReturn(0);
335 }
336 
337 /*@
338   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
339 
340   Not Collective
341 
342   Input Parameters:
343 . viewer - the HDF5 PetscViewer
344 
345   Output Parameters:
346 . flg - the flag
347 
348   Level: intermediate
349 
350   Notes:
351   This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned.
352   For more details, see PetscViewerHDF5SetCollective().
353 
354 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
355 
356 @*/
357 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
358 {
359   PetscErrorCode ierr;
360 
361   PetscFunctionBegin;
362   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
363   PetscValidPointer(flg,2);
364 
365   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
366   PetscFunctionReturn(0);
367 }
368 
369 static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
370 {
371   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
372   hid_t             plist_id;
373   PetscErrorCode    ierr;
374 
375   PetscFunctionBegin;
376   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
377   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
378   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
379   /* Set up file access property list with parallel I/O access */
380   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
381   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
382   /* Create or open the file collectively */
383   switch (hdf5->btype) {
384   case FILE_MODE_READ:
385     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
386     break;
387   case FILE_MODE_APPEND:
388     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
389     break;
390   case FILE_MODE_WRITE:
391     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
392     break;
393   default:
394     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
395   }
396   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
397   PetscStackCallHDF5(H5Pclose,(plist_id));
398   PetscFunctionReturn(0);
399 }
400 
401 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
402 {
403   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
404 
405   PetscFunctionBegin;
406   *name = vhdf5->filename;
407   PetscFunctionReturn(0);
408 }
409 
410 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
411 {
412   /*
413   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
414   PetscErrorCode   ierr;
415   */
416 
417   PetscFunctionBegin;
418   PetscFunctionReturn(0);
419 }
420 
421 /*MC
422    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
423 
424 
425 .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
426            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
427            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
428            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
429 
430   Level: beginner
431 M*/
432 
433 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
434 {
435   PetscViewer_HDF5 *hdf5;
436   PetscErrorCode   ierr;
437 
438   PetscFunctionBegin;
439   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
440 
441   v->data                = (void*) hdf5;
442   v->ops->destroy        = PetscViewerDestroy_HDF5;
443   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
444   v->ops->setup          = PetscViewerSetUp_HDF5;
445   v->ops->view           = PetscViewerView_HDF5;
446   v->ops->flush          = 0;
447   hdf5->btype            = (PetscFileMode) -1;
448   hdf5->filename         = 0;
449   hdf5->timestep         = -1;
450   hdf5->groups           = NULL;
451 
452   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
453 
454   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
455   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
456   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
457   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
458   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
459   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
460   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
461   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
462   PetscFunctionReturn(0);
463 }
464 
465 /*@C
466    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
467 
468    Collective on MPI_Comm
469 
470    Input Parameters:
471 +  comm - MPI communicator
472 .  name - name of file
473 -  type - type of file
474 $    FILE_MODE_WRITE - create new file for binary output
475 $    FILE_MODE_READ - open existing file for binary input
476 $    FILE_MODE_APPEND - open existing file for binary output
477 
478    Output Parameter:
479 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
480 
481   Options Database:
482 .  -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
483 .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
484 
485    Level: beginner
486 
487    Note:
488    This PetscViewer should be destroyed with PetscViewerDestroy().
489 
490    Concepts: HDF5 files
491    Concepts: PetscViewerHDF5^creating
492 
493 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
494           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
495           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
496 @*/
497 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
498 {
499   PetscErrorCode ierr;
500 
501   PetscFunctionBegin;
502   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
503   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
504   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
505   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
506   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
507   PetscFunctionReturn(0);
508 }
509 
510 /*@C
511   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
512 
513   Not collective
514 
515   Input Parameter:
516 . viewer - the PetscViewer
517 
518   Output Parameter:
519 . file_id - The file id
520 
521   Level: intermediate
522 
523 .seealso: PetscViewerHDF5Open()
524 @*/
525 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
526 {
527   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
528 
529   PetscFunctionBegin;
530   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
531   if (file_id) *file_id = hdf5->file_id;
532   PetscFunctionReturn(0);
533 }
534 
535 /*@C
536   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
537 
538   Not collective
539 
540   Input Parameters:
541 + viewer - the PetscViewer
542 - name - The group name
543 
544   Level: intermediate
545 
546 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
547 @*/
548 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
549 {
550   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
551   PetscViewerHDF5GroupList *groupNode;
552   PetscErrorCode   ierr;
553 
554   PetscFunctionBegin;
555   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
556   PetscValidCharPointer(name,2);
557   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
558   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
559 
560   groupNode->next = hdf5->groups;
561   hdf5->groups    = groupNode;
562   PetscFunctionReturn(0);
563 }
564 
565 /*@
566   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
567 
568   Not collective
569 
570   Input Parameter:
571 . viewer - the PetscViewer
572 
573   Level: intermediate
574 
575 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
576 @*/
577 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
578 {
579   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
580   PetscViewerHDF5GroupList *groupNode;
581   PetscErrorCode   ierr;
582 
583   PetscFunctionBegin;
584   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
585   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
586   groupNode    = hdf5->groups;
587   hdf5->groups = hdf5->groups->next;
588   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
589   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
590   PetscFunctionReturn(0);
591 }
592 
593 /*@C
594   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
595   If none has been assigned, returns NULL.
596 
597   Not collective
598 
599   Input Parameter:
600 . viewer - the PetscViewer
601 
602   Output Parameter:
603 . name - The group name
604 
605   Level: intermediate
606 
607 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
608 @*/
609 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
610 {
611   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
612 
613   PetscFunctionBegin;
614   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
615   PetscValidPointer(name,2);
616   if (hdf5->groups) *name = hdf5->groups->name;
617   else *name = NULL;
618   PetscFunctionReturn(0);
619 }
620 
621 /*@
622   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
623   and return this group's ID and file ID.
624   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
625 
626   Not collective
627 
628   Input Parameter:
629 . viewer - the PetscViewer
630 
631   Output Parameter:
632 + fileId - The HDF5 file ID
633 - groupId - The HDF5 group ID
634 
635   Notes:
636   If the viewer is writable, the group is created if it doesn't exist yet.
637 
638   Level: intermediate
639 
640 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
641 @*/
642 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
643 {
644   hid_t          file_id;
645   H5O_type_t     type;
646   const char     *groupName = NULL;
647   PetscBool      create;
648   PetscErrorCode ierr;
649 
650   PetscFunctionBegin;
651   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
652   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
653   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
654   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
655   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
656   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
657   *fileId  = file_id;
658   PetscFunctionReturn(0);
659 }
660 
661 /*@
662   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
663 
664   Not collective
665 
666   Input Parameter:
667 . viewer - the PetscViewer
668 
669   Level: intermediate
670 
671 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
672 @*/
673 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
674 {
675   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
676 
677   PetscFunctionBegin;
678   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
679   ++hdf5->timestep;
680   PetscFunctionReturn(0);
681 }
682 
683 /*@
684   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
685   of -1 disables blocking with timesteps.
686 
687   Not collective
688 
689   Input Parameters:
690 + viewer - the PetscViewer
691 - timestep - The timestep number
692 
693   Level: intermediate
694 
695 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
696 @*/
697 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
698 {
699   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
700 
701   PetscFunctionBegin;
702   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
703   hdf5->timestep = timestep;
704   PetscFunctionReturn(0);
705 }
706 
707 /*@
708   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
709 
710   Not collective
711 
712   Input Parameter:
713 . viewer - the PetscViewer
714 
715   Output Parameter:
716 . timestep - The timestep number
717 
718   Level: intermediate
719 
720 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
721 @*/
722 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
723 {
724   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
725 
726   PetscFunctionBegin;
727   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
728   PetscValidPointer(timestep,2);
729   *timestep = hdf5->timestep;
730   PetscFunctionReturn(0);
731 }
732 
733 /*@C
734   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
735 
736   Not collective
737 
738   Input Parameter:
739 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
740 
741   Output Parameter:
742 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
743 
744   Level: advanced
745 
746 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
747 @*/
748 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
749 {
750   PetscFunctionBegin;
751   if (ptype == PETSC_INT)
752 #if defined(PETSC_USE_64BIT_INDICES)
753                                        *htype = H5T_NATIVE_LLONG;
754 #else
755                                        *htype = H5T_NATIVE_INT;
756 #endif
757   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
758   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
759   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
760   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
761   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
762   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
763   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
764   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
765   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
766   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
767   PetscFunctionReturn(0);
768 }
769 
770 /*@C
771   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
772 
773   Not collective
774 
775   Input Parameter:
776 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
777 
778   Output Parameter:
779 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
780 
781   Level: advanced
782 
783 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
784 @*/
785 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
786 {
787   PetscFunctionBegin;
788 #if defined(PETSC_USE_64BIT_INDICES)
789   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
790   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
791 #else
792   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
793 #endif
794   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
795   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
796   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
797   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
798   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
799   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
800   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
801   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
802   PetscFunctionReturn(0);
803 }
804 
805 /*@C
806  PetscViewerHDF5WriteAttribute - Write an attribute
807 
808   Input Parameters:
809 + viewer - The HDF5 viewer
810 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
811 . name   - The attribute name
812 . datatype - The attribute type
813 - value    - The attribute value
814 
815   Level: advanced
816 
817 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
818 @*/
819 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
820 {
821   char           *parent;
822   hid_t          h5, dataspace, obj, attribute, dtype;
823   PetscBool      has;
824   PetscErrorCode ierr;
825 
826   PetscFunctionBegin;
827   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
828   if (dataset) PetscValidCharPointer(dataset, 2);
829   PetscValidCharPointer(name, 3);
830   PetscValidPointer(value, 5);
831   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
832   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
833   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
834   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
835   if (datatype == PETSC_STRING) {
836     size_t len;
837     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
838     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
839   }
840   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
841   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
842   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
843   if (has) {
844     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
845   } else {
846     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
847   }
848   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
849   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
850   PetscStackCallHDF5(H5Aclose,(attribute));
851   PetscStackCallHDF5(H5Oclose,(obj));
852   PetscStackCallHDF5(H5Sclose,(dataspace));
853   ierr = PetscFree(parent);CHKERRQ(ierr);
854   PetscFunctionReturn(0);
855 }
856 
857 /*@C
858  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
859 
860   Input Parameters:
861 + viewer   - The HDF5 viewer
862 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
863 . name     - The attribute name
864 . datatype - The attribute type
865 - value    - The attribute value
866 
867   Notes:
868   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
869   You might want to check first if it does using PetscViewerHDF5HasObject().
870 
871   Level: advanced
872 
873 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
874 @*/
875 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
876 {
877   PetscErrorCode ierr;
878 
879   PetscFunctionBegin;
880   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
881   PetscValidHeader(obj,2);
882   PetscValidCharPointer(name,3);
883   PetscValidPointer(value,5);
884   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
885   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
886   PetscFunctionReturn(0);
887 }
888 
889 /*@C
890  PetscViewerHDF5ReadAttribute - Read an attribute
891 
892   Input Parameters:
893 + viewer - The HDF5 viewer
894 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
895 . name   - The attribute name
896 - datatype - The attribute type
897 
898   Output Parameter:
899 . value    - The attribute value
900 
901   Level: advanced
902 
903 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
904 @*/
905 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
906 {
907   char           *parent;
908   hid_t          h5, obj, attribute, atype, dtype;
909   PetscBool      has;
910   PetscErrorCode ierr;
911 
912   PetscFunctionBegin;
913   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
914   if (dataset) PetscValidCharPointer(dataset, 2);
915   PetscValidCharPointer(name, 3);
916   PetscValidPointer(value, 5);
917   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
918   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
919   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
920   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
921   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
922   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
923   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
924   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
925   if (datatype == PETSC_STRING) {
926     size_t len;
927     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
928     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
929     PetscStackCallHDF5(H5Tclose,(atype));
930     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
931   }
932   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
933   PetscStackCallHDF5(H5Aclose,(attribute));
934   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
935   PetscStackCallHDF5(H5Oclose,(obj));
936   ierr = PetscFree(parent);CHKERRQ(ierr);
937   PetscFunctionReturn(0);
938 }
939 
940 /*@C
941  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
942 
943   Input Parameters:
944 + viewer   - The HDF5 viewer
945 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
946 . name     - The attribute name
947 - datatype - The attribute type
948 
949   Output Parameter:
950 . value    - The attribute value
951 
952   Notes:
953   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
954   You might want to check first if it does using PetscViewerHDF5HasObject().
955 
956   Level: advanced
957 
958 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
959 @*/
960 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
961 {
962   PetscErrorCode ierr;
963 
964   PetscFunctionBegin;
965   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
966   PetscValidHeader(obj,2);
967   PetscValidCharPointer(name,3);
968   PetscValidPointer(value, 5);
969   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
970   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
971   PetscFunctionReturn(0);
972 }
973 
974 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
975 {
976   htri_t exists;
977   hid_t group;
978 
979   PetscFunctionBegin;
980   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
981   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
982   if (!exists && createGroup) {
983     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
984     PetscStackCallHDF5(H5Gclose,(group));
985     exists = PETSC_TRUE;
986   }
987   *exists_ = (PetscBool) exists;
988   PetscFunctionReturn(0);
989 }
990 
991 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
992 {
993   const char     rootGroupName[] = "/";
994   hid_t          h5;
995   PetscBool      exists=PETSC_FALSE;
996   PetscInt       i;
997   int            n;
998   char           **hierarchy;
999   char           buf[PETSC_MAX_PATH_LEN]="";
1000   PetscErrorCode ierr;
1001 
1002   PetscFunctionBegin;
1003   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1004   if (name) PetscValidCharPointer(name, 2);
1005   else name = rootGroupName;
1006   if (has) {
1007     PetscValidIntPointer(has, 3);
1008     *has = PETSC_FALSE;
1009   }
1010   if (otype) {
1011     PetscValidIntPointer(otype, 4);
1012     *otype = H5O_TYPE_UNKNOWN;
1013   }
1014   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1015 
1016   /*
1017      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1018      Hence, each of them needs to be tested separately:
1019      1) whether it's a valid link
1020      2) whether this link resolves to an object
1021      See H5Oexists_by_name() documentation.
1022   */
1023   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
1024   if (!n) {
1025     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1026     if (has)   *has   = PETSC_TRUE;
1027     if (otype) *otype = H5O_TYPE_GROUP;
1028     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1029     PetscFunctionReturn(0);
1030   }
1031   for (i=0; i<n; i++) {
1032     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
1033     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1034     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1035     if (!exists) break;
1036   }
1037   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1038 
1039   /* If the object exists, get its type */
1040   if (exists && otype) {
1041     H5O_info_t info;
1042 
1043     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1044     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1045     *otype = info.type;
1046   }
1047   if (has) *has = exists;
1048   PetscFunctionReturn(0);
1049 }
1050 
1051 /*@
1052  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1053 
1054   Input Parameters:
1055 . viewer - The HDF5 viewer
1056 
1057   Output Parameter:
1058 . has    - Flag for group existence
1059 
1060   Notes:
1061   If the path exists but is not a group, this returns PETSC_FALSE as well.
1062 
1063   Level: advanced
1064 
1065 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1066 @*/
1067 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1068 {
1069   H5O_type_t type;
1070   const char *name;
1071   PetscErrorCode ierr;
1072 
1073   PetscFunctionBegin;
1074   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1075   PetscValidIntPointer(has,2);
1076   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
1077   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1078   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1079   PetscFunctionReturn(0);
1080 }
1081 
1082 /*@
1083  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1084 
1085   Input Parameters:
1086 + viewer - The HDF5 viewer
1087 - obj    - The named object
1088 
1089   Output Parameter:
1090 . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1091 
1092   Notes:
1093   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1094 
1095   Level: advanced
1096 
1097 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1098 @*/
1099 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1100 {
1101   H5O_type_t type;
1102   char *path;
1103   PetscErrorCode ierr;
1104 
1105   PetscFunctionBegin;
1106   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1107   PetscValidHeader(obj,2);
1108   PetscValidIntPointer(has,3);
1109   *has = PETSC_FALSE;
1110   if (!obj->name) PetscFunctionReturn(0);
1111   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
1112   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1113   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1114   ierr = PetscFree(path);CHKERRQ(ierr);
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 /*@C
1119  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1120 
1121   Input Parameters:
1122 + viewer - The HDF5 viewer
1123 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1124 - name   - The attribute name
1125 
1126   Output Parameter:
1127 . has    - Flag for attribute existence
1128 
1129   Level: advanced
1130 
1131 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1132 @*/
1133 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1134 {
1135   char           *parent;
1136   PetscErrorCode ierr;
1137 
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1140   if (dataset) PetscValidCharPointer(dataset,2);
1141   PetscValidCharPointer(name,3);
1142   PetscValidIntPointer(has,4);
1143   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1144   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
1145   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
1146   ierr = PetscFree(parent);CHKERRQ(ierr);
1147   PetscFunctionReturn(0);
1148 }
1149 
1150 /*@C
1151  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1152 
1153   Input Parameters:
1154 + viewer - The HDF5 viewer
1155 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1156 - name   - The attribute name
1157 
1158   Output Parameter:
1159 . has    - Flag for attribute existence
1160 
1161   Notes:
1162   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1163   You might want to check first if it does using PetscViewerHDF5HasObject().
1164 
1165   Level: advanced
1166 
1167 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1168 @*/
1169 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1170 {
1171   PetscErrorCode ierr;
1172 
1173   PetscFunctionBegin;
1174   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1175   PetscValidHeader(obj,2);
1176   PetscValidCharPointer(name,3);
1177   PetscValidIntPointer(has,4);
1178   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1179   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1180   PetscFunctionReturn(0);
1181 }
1182 
1183 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1184 {
1185   hid_t          h5;
1186   htri_t         hhas;
1187   PetscErrorCode ierr;
1188 
1189   PetscFunctionBegin;
1190   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1191   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1192   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1193   PetscFunctionReturn(0);
1194 }
1195 
1196 /*
1197   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1198   is attached to a communicator, in this case the attribute is a PetscViewer.
1199 */
1200 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1201 
1202 /*@C
1203   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1204 
1205   Collective on MPI_Comm
1206 
1207   Input Parameter:
1208 . comm - the MPI communicator to share the HDF5 PetscViewer
1209 
1210   Level: intermediate
1211 
1212   Options Database Keys:
1213 . -viewer_hdf5_filename <name>
1214 
1215   Environmental variables:
1216 . PETSC_VIEWER_HDF5_FILENAME
1217 
1218   Notes:
1219   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1220   an error code.  The HDF5 PetscViewer is usually used in the form
1221 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1222 
1223 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1224 @*/
1225 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1226 {
1227   PetscErrorCode ierr;
1228   PetscBool      flg;
1229   PetscViewer    viewer;
1230   char           fname[PETSC_MAX_PATH_LEN];
1231   MPI_Comm       ncomm;
1232 
1233   PetscFunctionBegin;
1234   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1235   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1236     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1237     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1238   }
1239   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1240   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1241   if (!flg) { /* PetscViewer not yet created */
1242     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1243     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1244     if (!flg) {
1245       ierr = PetscStrcpy(fname,"output.h5");
1246       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1247     }
1248     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1249     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1250     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1251     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1252     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1253     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1254   }
1255   ierr = PetscCommDestroy(&ncomm);
1256   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1257   PetscFunctionReturn(viewer);
1258 }
1259