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