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