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