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