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