xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision a041468a873ecb0eaa683d63e3e590db31f0074c)
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   Level: intermediate
666 
667 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
668 @*/
669 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
670 {
671   hid_t          file_id;
672   H5O_type_t     type;
673   const char     *groupName = NULL;
674   PetscErrorCode ierr;
675 
676   PetscFunctionBegin;
677   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
678   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
679   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, PETSC_TRUE, NULL, &type);CHKERRQ(ierr);
680   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
681   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
682   *fileId  = file_id;
683   PetscFunctionReturn(0);
684 }
685 
686 /*@
687   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
688 
689   Not collective
690 
691   Input Parameter:
692 . viewer - the PetscViewer
693 
694   Level: intermediate
695 
696 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
697 @*/
698 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
699 {
700   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
704   ++hdf5->timestep;
705   PetscFunctionReturn(0);
706 }
707 
708 /*@
709   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
710   of -1 disables blocking with timesteps.
711 
712   Not collective
713 
714   Input Parameters:
715 + viewer - the PetscViewer
716 - timestep - The timestep number
717 
718   Level: intermediate
719 
720 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
721 @*/
722 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
723 {
724   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
725 
726   PetscFunctionBegin;
727   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
728   hdf5->timestep = timestep;
729   PetscFunctionReturn(0);
730 }
731 
732 /*@
733   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
734 
735   Not collective
736 
737   Input Parameter:
738 . viewer - the PetscViewer
739 
740   Output Parameter:
741 . timestep - The timestep number
742 
743   Level: intermediate
744 
745 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
746 @*/
747 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
748 {
749   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
750 
751   PetscFunctionBegin;
752   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
753   PetscValidPointer(timestep,2);
754   *timestep = hdf5->timestep;
755   PetscFunctionReturn(0);
756 }
757 
758 /*@C
759   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
760 
761   Not collective
762 
763   Input Parameter:
764 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
765 
766   Output Parameter:
767 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
768 
769   Level: advanced
770 
771 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
772 @*/
773 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
774 {
775   PetscFunctionBegin;
776   if (ptype == PETSC_INT)
777 #if defined(PETSC_USE_64BIT_INDICES)
778                                        *htype = H5T_NATIVE_LLONG;
779 #else
780                                        *htype = H5T_NATIVE_INT;
781 #endif
782   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
783   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
784   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
785   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
786   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
787   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
788   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
789   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
790   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
791   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
792   PetscFunctionReturn(0);
793 }
794 
795 /*@C
796   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
797 
798   Not collective
799 
800   Input Parameter:
801 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
802 
803   Output Parameter:
804 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
805 
806   Level: advanced
807 
808 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
809 @*/
810 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
811 {
812   PetscFunctionBegin;
813 #if defined(PETSC_USE_64BIT_INDICES)
814   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
815   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
816 #else
817   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
818 #endif
819   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
820   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
821   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
822   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
823   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
824   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
825   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
826   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
827   PetscFunctionReturn(0);
828 }
829 
830 /*@C
831  PetscViewerHDF5WriteAttribute - Write an attribute
832 
833   Input Parameters:
834 + viewer - The HDF5 viewer
835 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
836 . name   - The attribute name
837 . datatype - The attribute type
838 - value    - The attribute value
839 
840   Level: advanced
841 
842 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
843 @*/
844 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
845 {
846   char           *parent;
847   hid_t          h5, dataspace, obj, attribute, dtype;
848   PetscBool      has;
849   PetscErrorCode ierr;
850 
851   PetscFunctionBegin;
852   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
853   if (dataset) PetscValidCharPointer(dataset, 2);
854   PetscValidCharPointer(name, 3);
855   PetscValidPointer(value, 5);
856   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
857   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
858   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
859   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
860   if (datatype == PETSC_STRING) {
861     size_t len;
862     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
863     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
864   }
865   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
866   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
867   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
868   if (has) {
869     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
870   } else {
871     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
872   }
873   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
874   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
875   PetscStackCallHDF5(H5Aclose,(attribute));
876   PetscStackCallHDF5(H5Oclose,(obj));
877   PetscStackCallHDF5(H5Sclose,(dataspace));
878   ierr = PetscFree(parent);CHKERRQ(ierr);
879   PetscFunctionReturn(0);
880 }
881 
882 /*@C
883  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
884 
885   Input Parameters:
886 + viewer   - The HDF5 viewer
887 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
888 . name     - The attribute name
889 . datatype - The attribute type
890 - value    - The attribute value
891 
892   Notes:
893   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
894   You might want to check first if it does using PetscViewerHDF5HasObject().
895 
896   Level: advanced
897 
898 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
899 @*/
900 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
901 {
902   PetscErrorCode ierr;
903 
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
906   PetscValidHeader(obj,2);
907   PetscValidCharPointer(name,3);
908   PetscValidPointer(value,5);
909   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
910   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
911   PetscFunctionReturn(0);
912 }
913 
914 /*@C
915  PetscViewerHDF5ReadAttribute - Read an attribute
916 
917   Input Parameters:
918 + viewer - The HDF5 viewer
919 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
920 . name   - The attribute name
921 - datatype - The attribute type
922 
923   Output Parameter:
924 . value    - The attribute value
925 
926   Level: advanced
927 
928 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
929 @*/
930 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
931 {
932   char           *parent;
933   hid_t          h5, obj, attribute, atype, dtype;
934   PetscBool      has;
935   PetscErrorCode ierr;
936 
937   PetscFunctionBegin;
938   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
939   if (dataset) PetscValidCharPointer(dataset, 2);
940   PetscValidCharPointer(name, 3);
941   PetscValidPointer(value, 5);
942   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
943   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
944   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
945   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
946   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
947   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
948   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
949   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
950   if (datatype == PETSC_STRING) {
951     size_t len;
952     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
953     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
954     PetscStackCallHDF5(H5Tclose,(atype));
955     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
956   }
957   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
958   PetscStackCallHDF5(H5Aclose,(attribute));
959   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
960   PetscStackCallHDF5(H5Oclose,(obj));
961   ierr = PetscFree(parent);CHKERRQ(ierr);
962   PetscFunctionReturn(0);
963 }
964 
965 /*@C
966  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
967 
968   Input Parameters:
969 + viewer   - The HDF5 viewer
970 . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
971 . name     - The attribute name
972 - datatype - The attribute type
973 
974   Output Parameter:
975 . value    - The attribute value
976 
977   Notes:
978   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
979   You might want to check first if it does using PetscViewerHDF5HasObject().
980 
981   Level: advanced
982 
983 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
984 @*/
985 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
986 {
987   PetscErrorCode ierr;
988 
989   PetscFunctionBegin;
990   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
991   PetscValidHeader(obj,2);
992   PetscValidCharPointer(name,3);
993   PetscValidPointer(value, 5);
994   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
995   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
996   PetscFunctionReturn(0);
997 }
998 
999 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1000 {
1001   htri_t exists;
1002   hid_t group;
1003 
1004   PetscFunctionBegin;
1005   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1006   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1007   if (!exists && createGroup) {
1008     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1009     PetscStackCallHDF5(H5Gclose,(group));
1010     exists = PETSC_TRUE;
1011   }
1012   *exists_ = (PetscBool) exists;
1013   PetscFunctionReturn(0);
1014 }
1015 
1016 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1017 {
1018   const char     rootGroupName[] = "/";
1019   hid_t          h5;
1020   PetscBool      exists=PETSC_FALSE;
1021   PetscInt       i;
1022   int            n;
1023   char           **hierarchy;
1024   char           buf[PETSC_MAX_PATH_LEN]="";
1025   PetscErrorCode ierr;
1026 
1027   PetscFunctionBegin;
1028   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1029   if (name) PetscValidCharPointer(name, 2);
1030   else name = rootGroupName;
1031   if (has) {
1032     PetscValidIntPointer(has, 3);
1033     *has = PETSC_FALSE;
1034   }
1035   if (otype) {
1036     PetscValidIntPointer(otype, 4);
1037     *otype = H5O_TYPE_UNKNOWN;
1038   }
1039   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1040 
1041   /*
1042      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
1043      Hence, each of them needs to be tested separately:
1044      1) whether it's a valid link
1045      2) whether this link resolves to an object
1046      See H5Oexists_by_name() documentation.
1047   */
1048   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
1049   if (!n) {
1050     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1051     if (has)   *has   = PETSC_TRUE;
1052     if (otype) *otype = H5O_TYPE_GROUP;
1053     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1054     PetscFunctionReturn(0);
1055   }
1056   for (i=0; i<n; i++) {
1057     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
1058     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1059     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1060     if (!exists) break;
1061   }
1062   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
1063 
1064   /* If the object exists, get its type */
1065   if (exists && otype) {
1066     H5O_info_t info;
1067 
1068     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1069     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
1070     *otype = info.type;
1071   }
1072   if (has) *has = exists;
1073   PetscFunctionReturn(0);
1074 }
1075 
1076 /*@
1077  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
1078 
1079   Input Parameters:
1080 . viewer - The HDF5 viewer
1081 
1082   Output Parameter:
1083 . has    - Flag for group existence
1084 
1085   Notes:
1086   If the path exists but is not a group, this returns PETSC_FALSE as well.
1087 
1088   Level: advanced
1089 
1090 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
1091 @*/
1092 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
1093 {
1094   H5O_type_t type;
1095   const char *name;
1096   PetscErrorCode ierr;
1097 
1098   PetscFunctionBegin;
1099   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1100   PetscValidIntPointer(has,2);
1101   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
1102   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1103   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
1104   PetscFunctionReturn(0);
1105 }
1106 
1107 /*@
1108  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1109 
1110   Input Parameters:
1111 + viewer - The HDF5 viewer
1112 - obj    - The named object
1113 
1114   Output Parameter:
1115 . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1116 
1117   Notes:
1118   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1119 
1120   Level: advanced
1121 
1122 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1123 @*/
1124 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1125 {
1126   H5O_type_t type;
1127   char *path;
1128   PetscErrorCode ierr;
1129 
1130   PetscFunctionBegin;
1131   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1132   PetscValidHeader(obj,2);
1133   PetscValidIntPointer(has,3);
1134   *has = PETSC_FALSE;
1135   if (!obj->name) PetscFunctionReturn(0);
1136   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
1137   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
1138   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
1139   ierr = PetscFree(path);CHKERRQ(ierr);
1140   PetscFunctionReturn(0);
1141 }
1142 
1143 /*@C
1144  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1145 
1146   Input Parameters:
1147 + viewer - The HDF5 viewer
1148 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1149 - name   - The attribute name
1150 
1151   Output Parameter:
1152 . has    - Flag for attribute existence
1153 
1154   Level: advanced
1155 
1156 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1157 @*/
1158 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1159 {
1160   char           *parent;
1161   PetscErrorCode ierr;
1162 
1163   PetscFunctionBegin;
1164   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1165   if (dataset) PetscValidCharPointer(dataset,2);
1166   PetscValidCharPointer(name,3);
1167   PetscValidIntPointer(has,4);
1168   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1169   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
1170   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
1171   ierr = PetscFree(parent);CHKERRQ(ierr);
1172   PetscFunctionReturn(0);
1173 }
1174 
1175 /*@C
1176  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1177 
1178   Input Parameters:
1179 + viewer - The HDF5 viewer
1180 . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1181 - name   - The attribute name
1182 
1183   Output Parameter:
1184 . has    - Flag for attribute existence
1185 
1186   Notes:
1187   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1188   You might want to check first if it does using PetscViewerHDF5HasObject().
1189 
1190   Level: advanced
1191 
1192 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1193 @*/
1194 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1195 {
1196   PetscErrorCode ierr;
1197 
1198   PetscFunctionBegin;
1199   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1200   PetscValidHeader(obj,2);
1201   PetscValidCharPointer(name,3);
1202   PetscValidIntPointer(has,4);
1203   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1204   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1209 {
1210   hid_t          h5;
1211   htri_t         hhas;
1212   PetscErrorCode ierr;
1213 
1214   PetscFunctionBegin;
1215   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
1216   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
1217   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1218   PetscFunctionReturn(0);
1219 }
1220 
1221 static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx)
1222 {
1223   HDF5ReadCtx    h=NULL;
1224   PetscErrorCode ierr;
1225 
1226   PetscFunctionBegin;
1227   ierr = PetscNew(&h);CHKERRQ(ierr);
1228   ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr);
1229   PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT));
1230   PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset));
1231   ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr);
1232   ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr);
1233   if (h->complexVal) {ierr = PetscViewerHDF5ReadAttribute(viewer,name,"complex",PETSC_BOOL,&h->complexVal);CHKERRQ(ierr);}
1234   /* MATLAB stores column vectors horizontally */
1235   ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr);
1236   /* Create property list for collective dataset read */
1237   PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER));
1238 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
1239   PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE));
1240 #endif
1241   /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */
1242   *ctx = h;
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx)
1247 {
1248   HDF5ReadCtx    h;
1249   PetscErrorCode ierr;
1250 
1251   PetscFunctionBegin;
1252   h = *ctx;
1253   PetscStackCallHDF5(H5Pclose,(h->plist));
1254   PetscStackCallHDF5(H5Gclose,(h->group));
1255   PetscStackCallHDF5(H5Sclose,(h->dataspace));
1256   PetscStackCallHDF5(H5Dclose,(h->dataset));
1257   ierr = PetscFree(*ctx);CHKERRQ(ierr);
1258   PetscFunctionReturn(0);
1259 }
1260 
1261 static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_)
1262 {
1263   int            rdim, dim;
1264   hsize_t        dims[4];
1265   PetscInt       bsInd, lenInd, bs, len, N;
1266   PetscLayout    map;
1267   PetscErrorCode ierr;
1268 
1269   PetscFunctionBegin;
1270   if (!(*map_)) {
1271     ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr);
1272   }
1273   map = *map_;
1274   /* calculate expected number of dimensions */
1275   dim = 0;
1276   if (ctx->timestep >= 0) ++dim;
1277   ++dim; /* length in blocks */
1278   if (ctx->complexVal) ++dim;
1279   /* get actual number of dimensions in dataset */
1280   PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL));
1281   /* calculate expected dimension indices */
1282   lenInd = 0;
1283   if (ctx->timestep >= 0) ++lenInd;
1284   bsInd = lenInd + 1;
1285   ctx->dim2 = PETSC_FALSE;
1286   if (rdim == dim) {
1287     bs = 1; /* support vectors stored as 1D array */
1288   } else if (rdim == dim+1) {
1289     bs = (PetscInt) dims[bsInd];
1290     if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */
1291   } else {
1292     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim);
1293   }
1294   len = dims[lenInd];
1295   if (ctx->horizontal) {
1296     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.");
1297     len = bs;
1298     bs = 1;
1299   }
1300   N = (PetscInt) len*bs;
1301 
1302   /* Set Vec sizes,blocksize,and type if not already set */
1303   if (map->bs < 0) {
1304     ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr);
1305   } 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);
1306   if (map->N < 0) {
1307     ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr);
1308   } 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);
1309   PetscFunctionReturn(0);
1310 }
1311 
1312 static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace)
1313 {
1314   hsize_t        count[4], offset[4];
1315   int            dim;
1316   PetscInt       bs, n, low;
1317   PetscErrorCode ierr;
1318 
1319   PetscFunctionBegin;
1320   /* Compute local size and ownership range */
1321   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1322   ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr);
1323   ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr);
1324   ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr);
1325 
1326   /* Each process defines a dataset and reads it from the hyperslab in the file */
1327   dim  = 0;
1328   if (ctx->timestep >= 0) {
1329     count[dim]  = 1;
1330     offset[dim] = ctx->timestep;
1331     ++dim;
1332   }
1333   if (ctx->horizontal) {
1334     count[dim]  = 1;
1335     offset[dim] = 0;
1336     ++dim;
1337   }
1338   {
1339     ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr);
1340     ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr);
1341     ++dim;
1342   }
1343   if (bs > 1 || ctx->dim2) {
1344     if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1");
1345     count[dim]  = bs;
1346     offset[dim] = 0;
1347     ++dim;
1348   }
1349   if (ctx->complexVal) {
1350     count[dim]  = 2;
1351     offset[dim] = 0;
1352     ++dim;
1353   }
1354   PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL));
1355   PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL));
1356   PetscFunctionReturn(0);
1357 }
1358 
1359 static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr)
1360 {
1361   PetscFunctionBegin;
1362   PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr));
1363   PetscFunctionReturn(0);
1364 }
1365 
1366 PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr)
1367 {
1368   HDF5ReadCtx     h=NULL;
1369   hid_t           memspace=0;
1370   size_t          unitsize;
1371   void            *arr;
1372   PetscErrorCode  ierr;
1373 
1374   PetscFunctionBegin;
1375   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1376 #if defined(PETSC_USE_COMPLEX)
1377   if (!h->complexVal) {
1378     H5T_class_t clazz = H5Tget_class(datatype);
1379     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.");
1380   }
1381 #else
1382   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.");
1383 #endif
1384 
1385   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1386   ierr = PetscLayoutSetUp(map);CHKERRQ(ierr);
1387   ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr);
1388 
1389   unitsize = H5Tget_size(datatype);
1390   if (h->complexVal) unitsize *= 2;
1391   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);
1392   ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr);
1393 
1394   ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr);
1395   PetscStackCallHDF5(H5Sclose,(memspace));
1396   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
1397   *newarr = arr;
1398   PetscFunctionReturn(0);
1399 }
1400 
1401 /*@C
1402  PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file.
1403 
1404   Input Parameters:
1405 + viewer - The HDF5 viewer
1406 - name   - The vector name
1407 
1408   Output Parameter:
1409 + bs     - block size
1410 - N      - global size
1411 
1412   Note:
1413   A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order:
1414   1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex).
1415 
1416   A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2().
1417 
1418   Level: advanced
1419 
1420 .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2()
1421 @*/
1422 PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N)
1423 {
1424   HDF5ReadCtx    h=NULL;
1425   PetscLayout    map=NULL;
1426   PetscErrorCode ierr;
1427 
1428   PetscFunctionBegin;
1429   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1430   ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr);
1431   ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr);
1432   ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr);
1433   if (bs) *bs = map->bs;
1434   if (N) *N = map->N;
1435   ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr);
1436   PetscFunctionReturn(0);
1437 }
1438 
1439 /*
1440   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1441   is attached to a communicator, in this case the attribute is a PetscViewer.
1442 */
1443 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1444 
1445 /*@C
1446   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1447 
1448   Collective on MPI_Comm
1449 
1450   Input Parameter:
1451 . comm - the MPI communicator to share the HDF5 PetscViewer
1452 
1453   Level: intermediate
1454 
1455   Options Database Keys:
1456 . -viewer_hdf5_filename <name>
1457 
1458   Environmental variables:
1459 . PETSC_VIEWER_HDF5_FILENAME
1460 
1461   Notes:
1462   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1463   an error code.  The HDF5 PetscViewer is usually used in the form
1464 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1465 
1466 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1467 @*/
1468 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1469 {
1470   PetscErrorCode ierr;
1471   PetscBool      flg;
1472   PetscViewer    viewer;
1473   char           fname[PETSC_MAX_PATH_LEN];
1474   MPI_Comm       ncomm;
1475 
1476   PetscFunctionBegin;
1477   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1478   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1479     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1480     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1481   }
1482   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1483   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1484   if (!flg) { /* PetscViewer not yet created */
1485     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1486     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1487     if (!flg) {
1488       ierr = PetscStrcpy(fname,"output.h5");
1489       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1490     }
1491     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1492     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1493     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1494     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1495     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1496     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1497   }
1498   ierr = PetscCommDestroy(&ncomm);
1499   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1500   PetscFunctionReturn(viewer);
1501 }
1502