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