xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 2bf68e3e0f2a61f71e7c65bee250bfa1c8ce0cdb)
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 /*MC
289    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
290 
291 
292 .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
293            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
294            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
295            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
296 
297   Level: beginner
298 M*/
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "PetscViewerCreate_HDF5"
302 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
303 {
304   PetscViewer_HDF5 *hdf5;
305   PetscErrorCode   ierr;
306 
307   PetscFunctionBegin;
308   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
309 
310   v->data                = (void*) hdf5;
311   v->ops->destroy        = PetscViewerDestroy_HDF5;
312   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
313   v->ops->flush          = 0;
314   hdf5->btype            = (PetscFileMode) -1;
315   hdf5->filename         = 0;
316   hdf5->timestep         = -1;
317   hdf5->groups           = NULL;
318 
319   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
320   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
321   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
322   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
323   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
324   PetscFunctionReturn(0);
325 }
326 
327 #undef __FUNCT__
328 #define __FUNCT__ "PetscViewerHDF5Open"
329 /*@C
330    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
331 
332    Collective on MPI_Comm
333 
334    Input Parameters:
335 +  comm - MPI communicator
336 .  name - name of file
337 -  type - type of file
338 $    FILE_MODE_WRITE - create new file for binary output
339 $    FILE_MODE_READ - open existing file for binary input
340 $    FILE_MODE_APPEND - open existing file for binary output
341 
342    Output Parameter:
343 .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
344 
345   Options Database:
346 .  -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
347 .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
348 
349    Level: beginner
350 
351    Note:
352    This PetscViewer should be destroyed with PetscViewerDestroy().
353 
354    Concepts: HDF5 files
355    Concepts: PetscViewerHDF5^creating
356 
357 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
358           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
359           MatLoad(), PetscFileMode, PetscViewer
360 @*/
361 PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
362 {
363   PetscErrorCode ierr;
364 
365   PetscFunctionBegin;
366   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
367   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
368   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
369   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373 #undef __FUNCT__
374 #define __FUNCT__ "PetscViewerHDF5GetFileId"
375 /*@C
376   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
377 
378   Not collective
379 
380   Input Parameter:
381 . viewer - the PetscViewer
382 
383   Output Parameter:
384 . file_id - The file id
385 
386   Level: intermediate
387 
388 .seealso: PetscViewerHDF5Open()
389 @*/
390 PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
391 {
392   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
393 
394   PetscFunctionBegin;
395   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
396   if (file_id) *file_id = hdf5->file_id;
397   PetscFunctionReturn(0);
398 }
399 
400 #undef __FUNCT__
401 #define __FUNCT__ "PetscViewerHDF5PushGroup"
402 /*@C
403   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
404 
405   Not collective
406 
407   Input Parameters:
408 + viewer - the PetscViewer
409 - name - The group name
410 
411   Level: intermediate
412 
413 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
414 @*/
415 PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
416 {
417   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
418   GroupList        *groupNode;
419   PetscErrorCode   ierr;
420 
421   PetscFunctionBegin;
422   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
423   PetscValidCharPointer(name,2);
424   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
425   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
426 
427   groupNode->next = hdf5->groups;
428   hdf5->groups    = groupNode;
429   PetscFunctionReturn(0);
430 }
431 
432 #undef __FUNCT__
433 #define __FUNCT__ "PetscViewerHDF5PopGroup"
434 /*@
435   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
436 
437   Not collective
438 
439   Input Parameter:
440 . viewer - the PetscViewer
441 
442   Level: intermediate
443 
444 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup()
445 @*/
446 PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
447 {
448   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
449   GroupList        *groupNode;
450   PetscErrorCode   ierr;
451 
452   PetscFunctionBegin;
453   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
454   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
455   groupNode    = hdf5->groups;
456   hdf5->groups = hdf5->groups->next;
457   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
458   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
459   PetscFunctionReturn(0);
460 }
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "PetscViewerHDF5GetGroup"
464 /*@C
465   PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL.
466 
467   Not collective
468 
469   Input Parameter:
470 . viewer - the PetscViewer
471 
472   Output Parameter:
473 . name - The group name
474 
475   Level: intermediate
476 
477 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup()
478 @*/
479 PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
480 {
481   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
482 
483   PetscFunctionBegin;
484   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
485   PetscValidPointer(name,2);
486   if (hdf5->groups) *name = hdf5->groups->name;
487   else *name = NULL;
488   PetscFunctionReturn(0);
489 }
490 
491 #undef __FUNCT__
492 #define __FUNCT__ "PetscViewerHDF5IncrementTimestep"
493 /*@
494   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
495 
496   Not collective
497 
498   Input Parameter:
499 . viewer - the PetscViewer
500 
501   Level: intermediate
502 
503 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
504 @*/
505 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
506 {
507   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
508 
509   PetscFunctionBegin;
510   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
511   ++hdf5->timestep;
512   PetscFunctionReturn(0);
513 }
514 
515 #undef __FUNCT__
516 #define __FUNCT__ "PetscViewerHDF5SetTimestep"
517 /*@
518   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
519   of -1 disables blocking with timesteps.
520 
521   Not collective
522 
523   Input Parameters:
524 + viewer - the PetscViewer
525 - timestep - The timestep number
526 
527   Level: intermediate
528 
529 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
530 @*/
531 PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
532 {
533   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
534 
535   PetscFunctionBegin;
536   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
537   hdf5->timestep = timestep;
538   PetscFunctionReturn(0);
539 }
540 
541 #undef __FUNCT__
542 #define __FUNCT__ "PetscViewerHDF5GetTimestep"
543 /*@
544   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
545 
546   Not collective
547 
548   Input Parameter:
549 . viewer - the PetscViewer
550 
551   Output Parameter:
552 . timestep - The timestep number
553 
554   Level: intermediate
555 
556 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
557 @*/
558 PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
559 {
560   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
561 
562   PetscFunctionBegin;
563   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
564   PetscValidPointer(timestep,2);
565   *timestep = hdf5->timestep;
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "PetscDataTypeToHDF5DataType"
571 /*@C
572   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
573 
574   Not collective
575 
576   Input Parameter:
577 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
578 
579   Output Parameter:
580 . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
581 
582   Level: advanced
583 
584 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
585 @*/
586 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
587 {
588   PetscFunctionBegin;
589   if (ptype == PETSC_INT)
590 #if defined(PETSC_USE_64BIT_INDICES)
591                                        *htype = H5T_NATIVE_LLONG;
592 #else
593                                        *htype = H5T_NATIVE_INT;
594 #endif
595   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
596   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
597   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
598   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_DOUBLE;
599   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_DOUBLE;
600   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
601   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
602   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
603   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
604   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
605   PetscFunctionReturn(0);
606 }
607 
608 #undef __FUNCT__
609 #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType"
610 /*@C
611   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
612 
613   Not collective
614 
615   Input Parameter:
616 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
617 
618   Output Parameter:
619 . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
620 
621   Level: advanced
622 
623 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
624 @*/
625 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
626 {
627   PetscFunctionBegin;
628 #if defined(PETSC_USE_64BIT_INDICES)
629   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
630   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
631 #else
632   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
633 #endif
634   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
635   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
636   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
637   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
638   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
639   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
640   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
641   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
642   PetscFunctionReturn(0);
643 }
644 
645 #undef __FUNCT__
646 #define __FUNCT__ "PetscViewerHDF5WriteAttribute"
647 /*@C
648  PetscViewerHDF5WriteAttribute - Write a scalar attribute
649 
650   Input Parameters:
651 + viewer - The HDF5 viewer
652 . parent - The parent name
653 . name   - The attribute name
654 . datatype - The attribute type
655 - value    - The attribute value
656 
657   Level: advanced
658 
659 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute()
660 @*/
661 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
662 {
663   hid_t          h5, dataspace, dataset, attribute, dtype;
664   PetscErrorCode ierr;
665 
666   PetscFunctionBegin;
667   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
668   PetscValidPointer(parent, 2);
669   PetscValidPointer(name, 3);
670   PetscValidPointer(value, 4);
671   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
672   if (datatype == PETSC_STRING) {
673     size_t len;
674     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
675     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
676   }
677   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
678   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
679 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
680   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
681   PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
682 #else
683   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
684   PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT));
685 #endif
686   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
687   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
688   PetscStackCallHDF5(H5Aclose,(attribute));
689   PetscStackCallHDF5(H5Dclose,(dataset));
690   PetscStackCallHDF5(H5Sclose,(dataspace));
691   PetscFunctionReturn(0);
692 }
693 
694 #undef __FUNCT__
695 #define __FUNCT__ "PetscViewerHDF5ReadAttribute"
696 /*@C
697  PetscViewerHDF5ReadAttribute - Read a scalar attribute
698 
699   Input Parameters:
700 + viewer - The HDF5 viewer
701 . parent - The parent name
702 . name   - The attribute name
703 - datatype - The attribute type
704 
705   Output Parameter:
706 . value    - The attribute value
707 
708   Level: advanced
709 
710 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute()
711 @*/
712 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
713 {
714   hid_t          h5, dataspace, dataset, attribute, dtype;
715   PetscErrorCode ierr;
716 
717   PetscFunctionBegin;
718   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
719   PetscValidPointer(parent, 2);
720   PetscValidPointer(name, 3);
721   PetscValidPointer(value, 4);
722   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
723   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
724   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
725 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
726   PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT));
727 #else
728   PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent));
729 #endif
730   PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name));
731   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
732   PetscStackCallHDF5(H5Aclose,(attribute));
733   PetscStackCallHDF5(H5Dclose,(dataset));
734   PetscStackCallHDF5(H5Sclose,(dataspace));
735   PetscFunctionReturn(0);
736 }
737 
738 #undef __FUNCT__
739 #define __FUNCT__ "PetscViewerHDF5HasObject"
740 static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has)
741 {
742   hid_t          h5;
743   PetscErrorCode ierr;
744 
745   PetscFunctionBegin;
746   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
747   PetscValidPointer(name, 2);
748   PetscValidPointer(has, 3);
749   *has = PETSC_FALSE;
750   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
751   if (H5Lexists(h5, name, H5P_DEFAULT)) {
752     H5O_info_t info;
753     hid_t      obj;
754 
755     PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT));
756     PetscStackCallHDF5(H5Oget_info,(obj, &info));
757     if (otype == info.type) *has = PETSC_TRUE;
758     PetscStackCallHDF5(H5Oclose,(obj));
759   }
760   PetscFunctionReturn(0);
761 }
762 
763 #undef __FUNCT__
764 #define __FUNCT__ "PetscViewerHDF5HasAttribute"
765 /*@C
766  PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists
767 
768   Input Parameters:
769 + viewer - The HDF5 viewer
770 . parent - The parent name
771 - name   - The attribute name
772 
773   Output Parameter:
774 . has    - Flag for attribute existence
775 
776   Level: advanced
777 
778 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute()
779 @*/
780 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
781 {
782   hid_t          h5, dataset;
783   htri_t         hhas;
784   PetscBool      exists;
785   PetscErrorCode ierr;
786 
787   PetscFunctionBegin;
788   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
789   PetscValidPointer(parent, 2);
790   PetscValidPointer(name, 3);
791   PetscValidPointer(has, 4);
792   *has = PETSC_FALSE;
793   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
794   ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr);
795   if (exists) {
796 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800)
797     PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT));
798 #else
799     PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent));
800 #endif
801     if (dataset < 0) PetscFunctionReturn(0);
802     PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name));
803     if (hhas < 0) {
804       PetscStackCallHDF5(H5Dclose,(dataset));
805       PetscFunctionReturn(0);
806     }
807     PetscStackCallHDF5(H5Dclose,(dataset));
808     *has = hhas ? PETSC_TRUE : PETSC_FALSE;
809   }
810   PetscFunctionReturn(0);
811 }
812 
813 /*
814   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
815   is attached to a communicator, in this case the attribute is a PetscViewer.
816 */
817 static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
818 
819 #undef __FUNCT__
820 #define __FUNCT__ "PETSC_VIEWER_HDF5_"
821 /*@C
822   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
823 
824   Collective on MPI_Comm
825 
826   Input Parameter:
827 . comm - the MPI communicator to share the HDF5 PetscViewer
828 
829   Level: intermediate
830 
831   Options Database Keys:
832 . -viewer_hdf5_filename <name>
833 
834   Environmental variables:
835 . PETSC_VIEWER_HDF5_FILENAME
836 
837   Notes:
838   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
839   an error code.  The HDF5 PetscViewer is usually used in the form
840 $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
841 
842 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
843 @*/
844 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
845 {
846   PetscErrorCode ierr;
847   PetscBool      flg;
848   PetscViewer    viewer;
849   char           fname[PETSC_MAX_PATH_LEN];
850   MPI_Comm       ncomm;
851 
852   PetscFunctionBegin;
853   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
854   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
855     ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
856     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
857   }
858   ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
859   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
860   if (!flg) { /* PetscViewer not yet created */
861     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
862     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
863     if (!flg) {
864       ierr = PetscStrcpy(fname,"output.h5");
865       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
866     }
867     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
868     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
869     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
870     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
871     ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
872     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
873   }
874   ierr = PetscCommDestroy(&ncomm);
875   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
876   PetscFunctionReturn(viewer);
877 }
878