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