xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 077c4feaa93e49b43362556d31b04df1b2d3b2d2)
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)
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\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; 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   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
332   H5FD_mpio_xfer_t  mode;
333 
334   PetscFunctionBegin;
335   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
336   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
337   PetscFunctionReturn(0);
338 }
339 
340 /*@
341   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
342 
343   Not Collective
344 
345   Input Parameters:
346 . viewer - the HDF5 PetscViewer
347 
348   Output Parameters:
349 . flg - the flag
350 
351   Level: intermediate
352 
353   Notes:
354   This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned.
355   For more details, see PetscViewerHDF5SetCollective().
356 
357 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
358 
359 @*/
360 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
361 {
362   PetscErrorCode ierr;
363 
364   PetscFunctionBegin;
365   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
366   PetscValidBoolPointer(flg,2);
367 
368   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
369   PetscFunctionReturn(0);
370 }
371 
372 static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
373 {
374   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
375   hid_t             plist_id;
376   PetscErrorCode    ierr;
377 
378   PetscFunctionBegin;
379   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
380   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
381   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
382   /* Set up file access property list with parallel I/O access */
383   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
384   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
385   /* Create or open the file collectively */
386   switch (hdf5->btype) {
387   case FILE_MODE_READ:
388     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
389     break;
390   case FILE_MODE_APPEND:
391     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
392     break;
393   case FILE_MODE_WRITE:
394     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
395     break;
396   default:
397     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
398   }
399   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
400   PetscStackCallHDF5(H5Pclose,(plist_id));
401   PetscFunctionReturn(0);
402 }
403 
404 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
405 {
406   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
407 
408   PetscFunctionBegin;
409   *name = vhdf5->filename;
410   PetscFunctionReturn(0);
411 }
412 
413 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
414 {
415   /*
416   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
417   PetscErrorCode   ierr;
418   */
419 
420   PetscFunctionBegin;
421   PetscFunctionReturn(0);
422 }
423 
424 /*MC
425    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
426 
427 
428 .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
429            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
430            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
431            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
432 
433   Level: beginner
434 M*/
435 
436 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
437 {
438   PetscViewer_HDF5 *hdf5;
439   PetscErrorCode   ierr;
440 
441   PetscFunctionBegin;
442   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
443 
444   v->data                = (void*) hdf5;
445   v->ops->destroy        = PetscViewerDestroy_HDF5;
446   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
447   v->ops->setup          = PetscViewerSetUp_HDF5;
448   v->ops->view           = PetscViewerView_HDF5;
449   v->ops->flush          = 0;
450   hdf5->btype            = (PetscFileMode) -1;
451   hdf5->filename         = 0;
452   hdf5->timestep         = -1;
453   hdf5->groups           = NULL;
454 
455   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
456 
457   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
458   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
459   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
460   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
461   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
462   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
463   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
464   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 /*@C
469    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
470 
471    Collective
472 
473    Input Parameters:
474 +  comm - MPI communicator
475 .  name - name of file
476 -  type - type of file
477 $    FILE_MODE_WRITE - create new file for binary output
478 $    FILE_MODE_READ - open existing file for binary input
479 $    FILE_MODE_APPEND - open existing file for binary output
480 
481    Output Parameter:
482 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
483 
484   Options Database:
485 +  -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
486 -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
487 
488    Level: beginner
489 
490    Note:
491    This PetscViewer should be destroyed with PetscViewerDestroy().
492 
493 
494 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
495           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
496           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
497 @*/
498 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
499 {
500   PetscErrorCode ierr;
501 
502   PetscFunctionBegin;
503   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
504   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
505   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
506   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
507   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
508   PetscFunctionReturn(0);
509 }
510 
511 /*@C
512   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
513 
514   Not collective
515 
516   Input Parameter:
517 . viewer - the PetscViewer
518 
519   Output Parameter:
520 . file_id - The file id
521 
522   Level: intermediate
523 
524 .seealso: PetscViewerHDF5Open()
525 @*/
526 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
527 {
528   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
529 
530   PetscFunctionBegin;
531   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
532   if (file_id) *file_id = hdf5->file_id;
533   PetscFunctionReturn(0);
534 }
535 
536 /*@C
537   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
538 
539   Not collective
540 
541   Input Parameters:
542 + viewer - the PetscViewer
543 - name - The group name
544 
545   Level: intermediate
546 
547 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
548 @*/
549 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
550 {
551   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
552   PetscViewerHDF5GroupList *groupNode;
553   PetscErrorCode   ierr;
554 
555   PetscFunctionBegin;
556   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
557   PetscValidCharPointer(name,2);
558   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
559   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
560 
561   groupNode->next = hdf5->groups;
562   hdf5->groups    = groupNode;
563   PetscFunctionReturn(0);
564 }
565 
566 /*@
567   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
568 
569   Not collective
570 
571   Input Parameter:
572 . viewer - the PetscViewer
573 
574   Level: intermediate
575 
576 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
577 @*/
578 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
579 {
580   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
581   PetscViewerHDF5GroupList *groupNode;
582   PetscErrorCode   ierr;
583 
584   PetscFunctionBegin;
585   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
586   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
587   groupNode    = hdf5->groups;
588   hdf5->groups = hdf5->groups->next;
589   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
590   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
591   PetscFunctionReturn(0);
592 }
593 
594 /*@C
595   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
596   If none has been assigned, returns NULL.
597 
598   Not collective
599 
600   Input Parameter:
601 . viewer - the PetscViewer
602 
603   Output Parameter:
604 . name - The group name
605 
606   Level: intermediate
607 
608 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
609 @*/
610 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
611 {
612   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
613 
614   PetscFunctionBegin;
615   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
616   PetscValidPointer(name,2);
617   if (hdf5->groups) *name = hdf5->groups->name;
618   else *name = NULL;
619   PetscFunctionReturn(0);
620 }
621 
622 /*@
623   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
624   and return this group's ID and file ID.
625   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
626 
627   Not collective
628 
629   Input Parameter:
630 . viewer - the PetscViewer
631 
632   Output Parameter:
633 + fileId - The HDF5 file ID
634 - groupId - The HDF5 group ID
635 
636   Notes:
637   If the viewer is writable, the group is created if it doesn't exist yet.
638 
639   Level: intermediate
640 
641 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
642 @*/
643 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
644 {
645   hid_t          file_id;
646   H5O_type_t     type;
647   const char     *groupName = NULL;
648   PetscBool      create;
649   PetscErrorCode ierr;
650 
651   PetscFunctionBegin;
652   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
653   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
654   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
655   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
656   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
657   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
658   *fileId  = file_id;
659   PetscFunctionReturn(0);
660 }
661 
662 /*@
663   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
664 
665   Not collective
666 
667   Input Parameter:
668 . viewer - the PetscViewer
669 
670   Level: intermediate
671 
672 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
673 @*/
674 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
675 {
676   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
677 
678   PetscFunctionBegin;
679   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
680   ++hdf5->timestep;
681   PetscFunctionReturn(0);
682 }
683 
684 /*@
685   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
686   of -1 disables blocking with timesteps.
687 
688   Not collective
689 
690   Input Parameters:
691 + viewer - the PetscViewer
692 - timestep - The timestep number
693 
694   Level: intermediate
695 
696 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
697 @*/
698 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
699 {
700   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
704   hdf5->timestep = timestep;
705   PetscFunctionReturn(0);
706 }
707 
708 /*@
709   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
710 
711   Not collective
712 
713   Input Parameter:
714 . viewer - the PetscViewer
715 
716   Output Parameter:
717 . timestep - The timestep number
718 
719   Level: intermediate
720 
721 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
722 @*/
723 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
724 {
725   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
726 
727   PetscFunctionBegin;
728   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
729   PetscValidPointer(timestep,2);
730   *timestep = hdf5->timestep;
731   PetscFunctionReturn(0);
732 }
733 
734 /*@C
735   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
736 
737   Not collective
738 
739   Input Parameter:
740 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
741 
742   Output Parameter:
743 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
744 
745   Level: advanced
746 
747 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
748 @*/
749 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
750 {
751   PetscFunctionBegin;
752   if (ptype == PETSC_INT)
753 #if defined(PETSC_USE_64BIT_INDICES)
754                                        *htype = H5T_NATIVE_LLONG;
755 #else
756                                        *htype = H5T_NATIVE_INT;
757 #endif
758   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
759   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
760   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
761   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
762   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
763   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
764   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
765   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
766   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
767   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
768   PetscFunctionReturn(0);
769 }
770 
771 /*@C
772   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
773 
774   Not collective
775 
776   Input Parameter:
777 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
778 
779   Output Parameter:
780 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
781 
782   Level: advanced
783 
784 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
785 @*/
786 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
787 {
788   PetscFunctionBegin;
789 #if defined(PETSC_USE_64BIT_INDICES)
790   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
791   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
792 #else
793   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
794 #endif
795   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
796   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
797   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
798   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
799   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
800   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
801   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
802   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
803   PetscFunctionReturn(0);
804 }
805 
806 /*@C
807  PetscViewerHDF5WriteAttribute - Write an attribute
808 
809   Input Parameters:
810 + viewer - The HDF5 viewer
811 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
812 . name   - The attribute name
813 . datatype - The attribute type
814 - value    - The attribute value
815 
816   Level: advanced
817 
818 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
819 @*/
820 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
821 {
822   char           *parent;
823   hid_t          h5, dataspace, obj, attribute, dtype;
824   PetscBool      has;
825   PetscErrorCode ierr;
826 
827   PetscFunctionBegin;
828   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
829   if (dataset) PetscValidCharPointer(dataset, 2);
830   PetscValidCharPointer(name, 3);
831   PetscValidPointer(value, 5);
832   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
833   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
834   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
835   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
836   if (datatype == PETSC_STRING) {
837     size_t len;
838     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
839     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
840   }
841   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
842   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
843   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
844   if (has) {
845     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
846   } else {
847     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
848   }
849   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
850   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
851   PetscStackCallHDF5(H5Aclose,(attribute));
852   PetscStackCallHDF5(H5Oclose,(obj));
853   PetscStackCallHDF5(H5Sclose,(dataspace));
854   ierr = PetscFree(parent);CHKERRQ(ierr);
855   PetscFunctionReturn(0);
856 }
857 
858 /*@C
859  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
860 
861   Input Parameters:
862 + viewer   - The HDF5 viewer
863 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
864 . name     - The attribute name
865 . datatype - The attribute type
866 - value    - The attribute value
867 
868   Notes:
869   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
870   You might want to check first if it does using PetscViewerHDF5HasObject().
871 
872   Level: advanced
873 
874 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
875 @*/
876 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
877 {
878   PetscErrorCode ierr;
879 
880   PetscFunctionBegin;
881   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
882   PetscValidHeader(obj,2);
883   PetscValidCharPointer(name,3);
884   PetscValidPointer(value,5);
885   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
886   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
887   PetscFunctionReturn(0);
888 }
889 
890 /*@C
891  PetscViewerHDF5ReadAttribute - Read an attribute
892 
893   Input Parameters:
894 + viewer - The HDF5 viewer
895 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
896 . name   - The attribute name
897 - datatype - The attribute type
898 
899   Output Parameter:
900 . value    - The attribute value
901 
902   Level: advanced
903 
904 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
905 @*/
906 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
907 {
908   char           *parent;
909   hid_t          h5, obj, attribute, atype, dtype;
910   PetscBool      has;
911   PetscErrorCode ierr;
912 
913   PetscFunctionBegin;
914   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
915   if (dataset) PetscValidCharPointer(dataset, 2);
916   PetscValidCharPointer(name, 3);
917   PetscValidPointer(value, 5);
918   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
919   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
920   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
921   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
922   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
923   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
924   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
925   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
926   if (datatype == PETSC_STRING) {
927     size_t len;
928     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
929     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
930     PetscStackCallHDF5(H5Tclose,(atype));
931     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
932   }
933   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
934   PetscStackCallHDF5(H5Aclose,(attribute));
935   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
936   PetscStackCallHDF5(H5Oclose,(obj));
937   ierr = PetscFree(parent);CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 /*@C
942  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
943 
944   Input Parameters:
945 + viewer   - The HDF5 viewer
946 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
947 . name     - The attribute name
948 - datatype - The attribute type
949 
950   Output Parameter:
951 . value    - The attribute value
952 
953   Notes:
954   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
955   You might want to check first if it does using PetscViewerHDF5HasObject().
956 
957   Level: advanced
958 
959 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
960 @*/
961 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
962 {
963   PetscErrorCode ierr;
964 
965   PetscFunctionBegin;
966   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
967   PetscValidHeader(obj,2);
968   PetscValidCharPointer(name,3);
969   PetscValidPointer(value, 5);
970   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
971   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
972   PetscFunctionReturn(0);
973 }
974 
975 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
976 {
977   htri_t exists;
978   hid_t group;
979 
980   PetscFunctionBegin;
981   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
982   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
983   if (!exists && createGroup) {
984     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
985     PetscStackCallHDF5(H5Gclose,(group));
986     exists = PETSC_TRUE;
987   }
988   *exists_ = (PetscBool) exists;
989   PetscFunctionReturn(0);
990 }
991 
992 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
993 {
994   const char     rootGroupName[] = "/";
995   hid_t          h5;
996   PetscBool      exists=PETSC_FALSE;
997   PetscInt       i;
998   int            n;
999   char           **hierarchy;
1000   char           buf[PETSC_MAX_PATH_LEN]="";
1001   PetscErrorCode ierr;
1002 
1003   PetscFunctionBegin;
1004   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1005   if (name) PetscValidCharPointer(name, 2);
1006   else name = rootGroupName;
1007   if (has) {
1008     PetscValidIntPointer(has, 3);
1009     *has = PETSC_FALSE;
1010   }
1011   if (otype) {
1012     PetscValidIntPointer(otype, 4);
1013     *otype = H5O_TYPE_UNKNOWN;
1014   }
1015   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1016 
1017   /*
1018      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1019      Hence, each of them needs to be tested separately:
1020      1) whether it's a valid link
1021      2) whether this link resolves to an object
1022      See H5Oexists_by_name() documentation.
1023   */
1024   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
1025   if (!n) {
1026     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1027     if (has)   *has   = PETSC_TRUE;
1028     if (otype) *otype = H5O_TYPE_GROUP;
1029     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1030     PetscFunctionReturn(0);
1031   }
1032   for (i=0; i<n; i++) {
1033     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
1034     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1035     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1036     if (!exists) break;
1037   }
1038   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1039 
1040   /* If the object exists, get its type */
1041   if (exists && otype) {
1042     H5O_info_t info;
1043 
1044     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1045     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1046     *otype = info.type;
1047   }
1048   if (has) *has = exists;
1049   PetscFunctionReturn(0);
1050 }
1051 
1052 /*@
1053  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1054 
1055   Input Parameters:
1056 . viewer - The HDF5 viewer
1057 
1058   Output Parameter:
1059 . has    - Flag for group existence
1060 
1061   Notes:
1062   If the path exists but is not a group, this returns PETSC_FALSE as well.
1063 
1064   Level: advanced
1065 
1066 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1067 @*/
1068 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1069 {
1070   H5O_type_t type;
1071   const char *name;
1072   PetscErrorCode ierr;
1073 
1074   PetscFunctionBegin;
1075   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1076   PetscValidIntPointer(has,2);
1077   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
1078   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1079   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 /*@
1084  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1085 
1086   Input Parameters:
1087 + viewer - The HDF5 viewer
1088 - obj    - The named object
1089 
1090   Output Parameter:
1091 . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1092 
1093   Notes:
1094   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1095 
1096   Level: advanced
1097 
1098 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1099 @*/
1100 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1101 {
1102   H5O_type_t type;
1103   char *path;
1104   PetscErrorCode ierr;
1105 
1106   PetscFunctionBegin;
1107   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1108   PetscValidHeader(obj,2);
1109   PetscValidIntPointer(has,3);
1110   *has = PETSC_FALSE;
1111   if (!obj->name) PetscFunctionReturn(0);
1112   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
1113   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1114   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1115   ierr = PetscFree(path);CHKERRQ(ierr);
1116   PetscFunctionReturn(0);
1117 }
1118 
1119 /*@C
1120  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1121 
1122   Input Parameters:
1123 + viewer - The HDF5 viewer
1124 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1125 - name   - The attribute name
1126 
1127   Output Parameter:
1128 . has    - Flag for attribute existence
1129 
1130   Level: advanced
1131 
1132 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1133 @*/
1134 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1135 {
1136   char           *parent;
1137   PetscErrorCode ierr;
1138 
1139   PetscFunctionBegin;
1140   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1141   if (dataset) PetscValidCharPointer(dataset,2);
1142   PetscValidCharPointer(name,3);
1143   PetscValidIntPointer(has,4);
1144   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1145   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
1146   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
1147   ierr = PetscFree(parent);CHKERRQ(ierr);
1148   PetscFunctionReturn(0);
1149 }
1150 
1151 /*@C
1152  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1153 
1154   Input Parameters:
1155 + viewer - The HDF5 viewer
1156 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1157 - name   - The attribute name
1158 
1159   Output Parameter:
1160 . has    - Flag for attribute existence
1161 
1162   Notes:
1163   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1164   You might want to check first if it does using PetscViewerHDF5HasObject().
1165 
1166   Level: advanced
1167 
1168 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1169 @*/
1170 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1171 {
1172   PetscErrorCode ierr;
1173 
1174   PetscFunctionBegin;
1175   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1176   PetscValidHeader(obj,2);
1177   PetscValidCharPointer(name,3);
1178   PetscValidIntPointer(has,4);
1179   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1180   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1181   PetscFunctionReturn(0);
1182 }
1183 
1184 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1185 {
1186   hid_t          h5;
1187   htri_t         hhas;
1188   PetscErrorCode ierr;
1189 
1190   PetscFunctionBegin;
1191   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1192   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1193   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1194   PetscFunctionReturn(0);
1195 }
1196 
1197 /*
1198   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1199   is attached to a communicator, in this case the attribute is a PetscViewer.
1200 */
1201 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1202 
1203 /*@C
1204   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1205 
1206   Collective
1207 
1208   Input Parameter:
1209 . comm - the MPI communicator to share the HDF5 PetscViewer
1210 
1211   Level: intermediate
1212 
1213   Options Database Keys:
1214 . -viewer_hdf5_filename <name>
1215 
1216   Environmental variables:
1217 . PETSC_VIEWER_HDF5_FILENAME
1218 
1219   Notes:
1220   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1221   an error code.  The HDF5 PetscViewer is usually used in the form
1222 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1223 
1224 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1225 @*/
1226 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1227 {
1228   PetscErrorCode ierr;
1229   PetscBool      flg;
1230   PetscViewer    viewer;
1231   char           fname[PETSC_MAX_PATH_LEN];
1232   MPI_Comm       ncomm;
1233 
1234   PetscFunctionBegin;
1235   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1236   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1237     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1238     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1239   }
1240   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1241   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1242   if (!flg) { /* PetscViewer not yet created */
1243     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1244     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1245     if (!flg) {
1246       ierr = PetscStrcpy(fname,"output.h5");
1247       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1248     }
1249     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1250     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1251     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1252     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1253     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1254     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1255   }
1256   ierr = PetscCommDestroy(&ncomm);
1257   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1258   PetscFunctionReturn(viewer);
1259 }
1260