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