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