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