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