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