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