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