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