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