xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision af0996ce37bc06907c37d8d91773840993d61e62)
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 } PetscViewer_HDF5;
17 
18 #undef __FUNCT__
19 #define __FUNCT__ "PetscViewerSetFromOptions_HDF5"
20 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptions *PetscOptionsObject,PetscViewer v)
21 {
22   PetscErrorCode   ierr;
23   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
24 
25   PetscFunctionBegin;
26   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
27   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
28   ierr = PetscOptionsTail();CHKERRQ(ierr);
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "PetscViewerFileClose_HDF5"
34 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
35 {
36   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
37   PetscErrorCode   ierr;
38 
39   PetscFunctionBegin;
40   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
41   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
42   PetscFunctionReturn(0);
43 }
44 
45 #undef __FUNCT__
46 #define __FUNCT__ "PetscViewerDestroy_HDF5"
47 PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
48 {
49   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
50   PetscErrorCode   ierr;
51 
52   PetscFunctionBegin;
53   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
54   while (hdf5->groups) {
55     GroupList *tmp = hdf5->groups->next;
56 
57     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
58     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
59     hdf5->groups = tmp;
60   }
61   ierr = PetscFree(hdf5);CHKERRQ(ierr);
62   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
63   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "PetscViewerFileSetMode_HDF5"
69 PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
70 {
71   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
72 
73   PetscFunctionBegin;
74   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
75   hdf5->btype = type;
76   PetscFunctionReturn(0);
77 }
78 
79 #undef __FUNCT__
80 #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2_HDF5"
81 PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
82 {
83   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
84 
85   PetscFunctionBegin;
86   hdf5->basedimension2 = flg;
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "PetscViewerHDF5SetBaseDimension2"
92 /*@C
93      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
94        dimension of 2.
95 
96     Logically Collective on PetscViewer
97 
98   Input Parameters:
99 +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
100 -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
101 
102   Options Database:
103 .  -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
104 
105 
106   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
107          of one when the dimension is lower. Others think the option is crazy.
108 
109   Level: intermediate
110 
111 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
112 
113 @*/
114 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
115 {
116   PetscErrorCode ierr;
117 
118   PetscFunctionBegin;
119   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
120   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
121   PetscFunctionReturn(0);
122 }
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "PetscViewerHDF5GetBaseDimension2"
126 /*@C
127      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
128        dimension of 2.
129 
130     Logically Collective on PetscViewer
131 
132   Input Parameter:
133 .  viewer - the PetscViewer, must be of type HDF5
134 
135   Output Parameter:
136 .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
137 
138   Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
139          of one when the dimension is lower. Others think the option is crazy.
140 
141   Level: intermediate
142 
143 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
144 
145 @*/
146 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
147 {
148   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
149 
150   PetscFunctionBegin;
151   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
152   *flg = hdf5->basedimension2;
153   PetscFunctionReturn(0);
154 }
155 
156 #undef __FUNCT__
157 #define __FUNCT__ "PetscViewerFileSetName_HDF5"
158 PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
159 {
160   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
161 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
162   MPI_Info          info = MPI_INFO_NULL;
163 #endif
164   hid_t             plist_id;
165   PetscErrorCode    ierr;
166 
167   PetscFunctionBegin;
168   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
169   /* Set up file access property list with parallel I/O access */
170   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
171 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
172   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
173 #endif
174   /* Create or open the file collectively */
175   switch (hdf5->btype) {
176   case FILE_MODE_READ:
177     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
178     break;
179   case FILE_MODE_APPEND:
180     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
181     break;
182   case FILE_MODE_WRITE:
183     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
184     break;
185   default:
186     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
187   }
188   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
189   PetscStackCallHDF5(H5Pclose,(plist_id));
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "PetscViewerCreate_HDF5"
195 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
196 {
197   PetscViewer_HDF5 *hdf5;
198   PetscErrorCode   ierr;
199 
200   PetscFunctionBegin;
201   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
202 
203   v->data                = (void*) hdf5;
204   v->ops->destroy        = PetscViewerDestroy_HDF5;
205   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
206   v->ops->flush          = 0;
207   hdf5->btype            = (PetscFileMode) -1;
208   hdf5->filename         = 0;
209   hdf5->timestep         = -1;
210   hdf5->groups           = NULL;
211 
212   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
213   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
214   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 #undef __FUNCT__
219 #define __FUNCT__ "PetscViewerHDF5Open"
220 /*@C
221    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
222 
223    Collective on MPI_Comm
224 
225    Input Parameters:
226 +  comm - MPI communicator
227 .  name - name of file
228 -  type - type of file
229 $    FILE_MODE_WRITE - create new file for binary output
230 $    FILE_MODE_READ - open existing file for binary input
231 $    FILE_MODE_APPEND - open existing file for binary output
232 
233    Output Parameter:
234 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
235 
236   Options Database:
237 .  -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
238 
239    Level: beginner
240 
241    Note:
242    This PetscViewer should be destroyed with PetscViewerDestroy().
243 
244    Concepts: HDF5 files
245    Concepts: PetscViewerHDF5^creating
246 
247 .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), PetscViewerHDF5GetBaseDimension2()
248           VecView(), MatView(), VecLoad(), MatLoad(),
249           PetscFileMode, PetscViewer
250 @*/
251 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
252 {
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
257   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
258   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
259   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
260   PetscFunctionReturn(0);
261 }
262 
263 #undef __FUNCT__
264 #define __FUNCT__ "PetscViewerHDF5GetFileId"
265 /*@C
266   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
267 
268   Not collective
269 
270   Input Parameter:
271 . viewer - the PetscViewer
272 
273   Output Parameter:
274 . file_id - The file id
275 
276   Level: intermediate
277 
278 .seealso: PetscViewerHDF5Open()
279 @*/
280 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
281 {
282   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
283 
284   PetscFunctionBegin;
285   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
286   if (file_id) *file_id = hdf5->file_id;
287   PetscFunctionReturn(0);
288 }
289 
290 #undef __FUNCT__
291 #define __FUNCT__ "PetscViewerHDF5PushGroup"
292 /*@C
293   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
294 
295   Not collective
296 
297   Input Parameters:
298 + viewer - the PetscViewer
299 - name - The group name
300 
301   Level: intermediate
302 
303 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
304 @*/
305 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
306 {
307   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
308   GroupList        *groupNode;
309   PetscErrorCode   ierr;
310 
311   PetscFunctionBegin;
312   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
313   PetscValidCharPointer(name,2);
314   ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr);
315   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
316 
317   groupNode->next = hdf5->groups;
318   hdf5->groups    = groupNode;
319   PetscFunctionReturn(0);
320 }
321 
322 #undef __FUNCT__
323 #define __FUNCT__ "PetscViewerHDF5PopGroup"
324 /*@
325   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
326 
327   Not collective
328 
329   Input Parameter:
330 . viewer - the PetscViewer
331 
332   Level: intermediate
333 
334 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
335 @*/
336 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
337 {
338   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
339   GroupList        *groupNode;
340   PetscErrorCode   ierr;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
344   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
345   groupNode    = hdf5->groups;
346   hdf5->groups = hdf5->groups->next;
347   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
348   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
349   PetscFunctionReturn(0);
350 }
351 
352 #undef __FUNCT__
353 #define __FUNCT__ "PetscViewerHDF5GetGroup"
354 /*@C
355   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
356 
357   Not collective
358 
359   Input Parameter:
360 . viewer - the PetscViewer
361 
362   Output Parameter:
363 . name - The group name
364 
365   Level: intermediate
366 
367 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
368 @*/
369 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
370 {
371   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
372 
373   PetscFunctionBegin;
374   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
375   PetscValidPointer(name,2);
376   if (hdf5->groups) *name = hdf5->groups->name;
377   else *name = NULL;
378   PetscFunctionReturn(0);
379 }
380 
381 #undef __FUNCT__
382 #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
383 /*@
384   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
385 
386   Not collective
387 
388   Input Parameter:
389 . viewer - the PetscViewer
390 
391   Level: intermediate
392 
393 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
394 @*/
395 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
396 {
397   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
398 
399   PetscFunctionBegin;
400   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
401   ++hdf5->timestep;
402   PetscFunctionReturn(0);
403 }
404 
405 #undef __FUNCT__
406 #define __FUNCT__ "PetscViewerHDF5SetTimestep"
407 /*@
408   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
409   of -1 disables blocking with timesteps.
410 
411   Not collective
412 
413   Input Parameters:
414 + viewer - the PetscViewer
415 - timestep - The timestep number
416 
417   Level: intermediate
418 
419 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
420 @*/
421 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
422 {
423   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
424 
425   PetscFunctionBegin;
426   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
427   hdf5->timestep = timestep;
428   PetscFunctionReturn(0);
429 }
430 
431 #undef __FUNCT__
432 #define __FUNCT__ "PetscViewerHDF5GetTimestep"
433 /*@
434   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
435 
436   Not collective
437 
438   Input Parameter:
439 . viewer - the PetscViewer
440 
441   Output Parameter:
442 . timestep - The timestep number
443 
444   Level: intermediate
445 
446 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
447 @*/
448 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
449 {
450   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
451 
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
454   PetscValidPointer(timestep,2);
455   *timestep = hdf5->timestep;
456   PetscFunctionReturn(0);
457 }
458 
459 #undef __FUNCT__
460 #define __FUNCT__ "PetscDataTypeToHDF5DataType"
461 /*@C
462   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
463 
464   Not collective
465 
466   Input Parameter:
467 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
468 
469   Output Parameter:
470 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
471 
472   Level: advanced
473 
474 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
475 @*/
476 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
477 {
478   PetscFunctionBegin;
479   if (ptype == PETSC_INT)
480 #if defined(PETSC_USE_64BIT_INDICES)
481                                        *htype = H5T_NATIVE_LLONG;
482 #else
483                                        *htype = H5T_NATIVE_INT;
484 #endif
485   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
486   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
487   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
488   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
489   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
490   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
491   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
492   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
493   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
494   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
495   PetscFunctionReturn(0);
496 }
497 
498 #undef __FUNCT__
499 #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
500 /*@C
501   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
502 
503   Not collective
504 
505   Input Parameter:
506 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
507 
508   Output Parameter:
509 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
510 
511   Level: advanced
512 
513 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
514 @*/
515 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
516 {
517   PetscFunctionBegin;
518 #if defined(PETSC_USE_64BIT_INDICES)
519   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
520   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
521 #else
522   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
523 #endif
524   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
525   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
526   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
527   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
528   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
529   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
530   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
531   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
532   PetscFunctionReturn(0);
533 }
534 
535 #undef __FUNCT__
536 #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
537 /*@
538  PetscViewerHDF5WriteAttribute - Write a scalar attribute
539 
540   Input Parameters:
541 + viewer - The HDF5 viewer
542 . parent - The parent name
543 . name   - The attribute name
544 . datatype - The attribute type
545 - value    - The attribute value
546 
547   Level: advanced
548 
549 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
550 @*/
551 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
552 {
553   hid_t          h5, dataspace, dataset, attribute, dtype;
554   PetscErrorCode ierr;
555 
556   PetscFunctionBegin;
557   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
558   PetscValidPointer(parent, 2);
559   PetscValidPointer(name, 3);
560   PetscValidPointer(value, 4);
561   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
562   if (datatype == PETSC_STRING) {
563     size_t len;
564     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
565     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
566   }
567   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
568   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
569 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
570   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
571   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
572 #else
573   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
574   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
575 #endif
576   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
577   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
578   PetscStackCallHDF5(H5Aclose,(attribute));
579   PetscStackCallHDF5(H5Dclose,(dataset));
580   PetscStackCallHDF5(H5Sclose,(dataspace));
581   PetscFunctionReturn(0);
582 }
583 
584 #undef __FUNCT__
585 #define __FUNCT__ "PetscViewerHDF5ReadAttribute"
586 /*@
587  PetscViewerHDF5ReadAttribute - Read a scalar attribute
588 
589   Input Parameters:
590 + viewer - The HDF5 viewer
591 . parent - The parent name
592 . name   - The attribute name
593 - datatype - The attribute type
594 
595   Output Parameter:
596 . value    - The attribute value
597 
598   Level: advanced
599 
600 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
601 @*/
602 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
603 {
604   hid_t          h5, dataspace, dataset, attribute, dtype;
605   PetscErrorCode ierr;
606 
607   PetscFunctionBegin;
608   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
609   PetscValidPointer(parent, 2);
610   PetscValidPointer(name, 3);
611   PetscValidPointer(value, 4);
612   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
613   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
614   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
615 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
616   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
617 #else
618   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
619 #endif
620   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
621   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
622   PetscStackCallHDF5(H5Aclose,(attribute));
623   PetscStackCallHDF5(H5Dclose,(dataset));
624   PetscStackCallHDF5(H5Sclose,(dataspace));
625   PetscFunctionReturn(0);
626 }
627 
628 #undef __FUNCT__
629 #define __FUNCT__ "PetscViewerHDF5HasObject"
630 static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
631 {
632   hid_t          h5;
633   PetscErrorCode ierr;
634 
635   PetscFunctionBegin;
636   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
637   PetscValidPointer(name, 2);
638   PetscValidPointer(has, 3);
639   *has = PETSC_FALSE;
640   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
641   if (H5Lexists(h5, name, H5P_DEFAULT)) {
642     H5O_info_t info;
643     hid_t      obj;
644 
645     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
646     PetscStackCallHDF5(H5Oget_info,(obj, &info));
647     if (otype == info.type) *has = PETSC_TRUE;
648     PetscStackCallHDF5(H5Oclose,(obj));
649   }
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "PetscViewerHDF5HasAttribute"
655 /*@
656  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
657 
658   Input Parameters:
659 + viewer - The HDF5 viewer
660 . parent - The parent name
661 - name   - The attribute name
662 
663   Output Parameter:
664 . has    - Flag for attribute existence
665 
666   Level: advanced
667 
668 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
669 @*/
670 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
671 {
672   hid_t          h5, dataset;
673   htri_t         hhas;
674   PetscBool      exists;
675   PetscErrorCode ierr;
676 
677   PetscFunctionBegin;
678   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
679   PetscValidPointer(parent, 2);
680   PetscValidPointer(name, 3);
681   PetscValidPointer(has, 4);
682   *has = PETSC_FALSE;
683   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
684   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
685   if (exists) {
686 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
687     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
688 #else
689     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
690 #endif
691     if (dataset < 0) PetscFunctionReturn(0);
692     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
693     if (hhas < 0) {
694       PetscStackCallHDF5(H5Dclose,(dataset));
695       PetscFunctionReturn(0);
696     }
697     PetscStackCallHDF5(H5Dclose,(dataset));
698     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
699   }
700   PetscFunctionReturn(0);
701 }
702 
703 /*
704   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
705   is attached to a communicator, in this case the attribute is a PetscViewer.
706 */
707 static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
708 
709 #undef __FUNCT__
710 #define __FUNCT__ "PETSC_VIEWER_HDF5_"
711 /*@C
712   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
713 
714   Collective on MPI_Comm
715 
716   Input Parameter:
717 . comm - the MPI communicator to share the HDF5 PetscViewer
718 
719   Level: intermediate
720 
721   Options Database Keys:
722 . -viewer_hdf5_filename <name>
723 
724   Environmental variables:
725 . PETSC_VIEWER_HDF5_FILENAME
726 
727   Notes:
728   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
729   an error code.  The HDF5 PetscViewer is usually used in the form
730 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
731 
732 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
733 @*/
734 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
735 {
736   PetscErrorCode ierr;
737   PetscBool      flg;
738   PetscViewer    viewer;
739   char           fname[PETSC_MAX_PATH_LEN];
740   MPI_Comm       ncomm;
741 
742   PetscFunctionBegin;
743   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
744   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
745     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
746     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
747   }
748   ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
749   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
750   if (!flg) { /* PetscViewer not yet created */
751     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
752     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
753     if (!flg) {
754       ierr = PetscStrcpy(fname,"output.h5");
755       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
756     }
757     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
758     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
759     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
760     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
761     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
762     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
763   }
764   ierr = PetscCommDestroy(&ncomm);
765   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
766   PetscFunctionReturn(viewer);
767 }
768