xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision d083f849a86f1f43e18d534ee43954e2786cb29a)
1bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h>
2d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
3ded06c3fSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800)
4ded06c3fSVaclav Hapla #error "PETSc needs HDF5 version >= 1.8.0"
5ded06c3fSVaclav Hapla #endif
65c6c1daeSBarry Smith 
7bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
906db490cSVaclav Hapla 
105c6c1daeSBarry Smith typedef struct GroupList {
115c6c1daeSBarry Smith   const char       *name;
125c6c1daeSBarry Smith   struct GroupList *next;
135c6c1daeSBarry Smith } GroupList;
145c6c1daeSBarry Smith 
155c6c1daeSBarry Smith typedef struct {
165c6c1daeSBarry Smith   char          *filename;
175c6c1daeSBarry Smith   PetscFileMode btype;
185c6c1daeSBarry Smith   hid_t         file_id;
195c6c1daeSBarry Smith   PetscInt      timestep;
205c6c1daeSBarry Smith   GroupList     *groups;
2182ea9b62SBarry Smith   PetscBool     basedimension2;  /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */
229a0502c6SHåkon Strandenes   PetscBool     spoutput;  /* write data in single precision even if PETSc is compiled with double precision PetscReal */
235c6c1daeSBarry Smith } PetscViewer_HDF5;
245c6c1daeSBarry Smith 
256c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
266c132bc1SVaclav Hapla {
276c132bc1SVaclav Hapla   const char *group;
286c132bc1SVaclav Hapla   char buf[PETSC_MAX_PATH_LEN]="";
296c132bc1SVaclav Hapla   PetscErrorCode ierr;
306c132bc1SVaclav Hapla 
316c132bc1SVaclav Hapla   PetscFunctionBegin;
326c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
336c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
346c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
356c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
366c132bc1SVaclav Hapla   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
376c132bc1SVaclav Hapla   PetscFunctionReturn(0);
386c132bc1SVaclav Hapla }
396c132bc1SVaclav Hapla 
40577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
41577e0e04SVaclav Hapla {
42577e0e04SVaclav Hapla   PetscBool has;
43577e0e04SVaclav Hapla   const char *group;
44577e0e04SVaclav Hapla   PetscErrorCode ierr;
45577e0e04SVaclav Hapla 
46577e0e04SVaclav Hapla   PetscFunctionBegin;
47577e0e04SVaclav Hapla   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
48577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
49577e0e04SVaclav Hapla   if (!has) {
50577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
51577e0e04SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
52577e0e04SVaclav Hapla   }
53577e0e04SVaclav Hapla   PetscFunctionReturn(0);
54577e0e04SVaclav Hapla }
55577e0e04SVaclav Hapla 
564416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
5782ea9b62SBarry Smith {
5882ea9b62SBarry Smith   PetscErrorCode   ierr;
5982ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
6082ea9b62SBarry Smith 
6182ea9b62SBarry Smith   PetscFunctionBegin;
6282ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
6382ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
649a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
6582ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
6682ea9b62SBarry Smith   PetscFunctionReturn(0);
6782ea9b62SBarry Smith }
6882ea9b62SBarry Smith 
691b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
701b793a25SVaclav Hapla {
711b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
721b793a25SVaclav Hapla   PetscErrorCode    ierr;
731b793a25SVaclav Hapla 
741b793a25SVaclav Hapla   PetscFunctionBegin;
751b793a25SVaclav Hapla   if (hdf5->filename) {
761b793a25SVaclav Hapla     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
771b793a25SVaclav Hapla   }
781b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
791b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
801b793a25SVaclav Hapla   PetscFunctionReturn(0);
811b793a25SVaclav Hapla }
821b793a25SVaclav Hapla 
835c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
845c6c1daeSBarry Smith {
855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
865c6c1daeSBarry Smith   PetscErrorCode   ierr;
875c6c1daeSBarry Smith 
885c6c1daeSBarry Smith   PetscFunctionBegin;
895c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
90729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
915c6c1daeSBarry Smith   PetscFunctionReturn(0);
925c6c1daeSBarry Smith }
935c6c1daeSBarry Smith 
949b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
955c6c1daeSBarry Smith {
965c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
975c6c1daeSBarry Smith   PetscErrorCode   ierr;
985c6c1daeSBarry Smith 
995c6c1daeSBarry Smith   PetscFunctionBegin;
1005c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
1015c6c1daeSBarry Smith   while (hdf5->groups) {
1025c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
1035c6c1daeSBarry Smith 
1045c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
1055c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
1065c6c1daeSBarry Smith     hdf5->groups = tmp;
1075c6c1daeSBarry Smith   }
1085c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
1090b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
110d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1110b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
112058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
113058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1145c6c1daeSBarry Smith   PetscFunctionReturn(0);
1155c6c1daeSBarry Smith }
1165c6c1daeSBarry Smith 
1179b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1185c6c1daeSBarry Smith {
1195c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1205c6c1daeSBarry Smith 
1215c6c1daeSBarry Smith   PetscFunctionBegin;
1225c6c1daeSBarry Smith   hdf5->btype = type;
1235c6c1daeSBarry Smith   PetscFunctionReturn(0);
1245c6c1daeSBarry Smith }
1255c6c1daeSBarry Smith 
1260b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1270b62783dSVaclav Hapla {
1280b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1290b62783dSVaclav Hapla 
1300b62783dSVaclav Hapla   PetscFunctionBegin;
1310b62783dSVaclav Hapla   *type = hdf5->btype;
1320b62783dSVaclav Hapla   PetscFunctionReturn(0);
1330b62783dSVaclav Hapla }
1340b62783dSVaclav Hapla 
1359b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
13682ea9b62SBarry Smith {
13782ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
13882ea9b62SBarry Smith 
13982ea9b62SBarry Smith   PetscFunctionBegin;
14082ea9b62SBarry Smith   hdf5->basedimension2 = flg;
14182ea9b62SBarry Smith   PetscFunctionReturn(0);
14282ea9b62SBarry Smith }
14382ea9b62SBarry Smith 
144df863907SAlex Fikl /*@
14582ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
14682ea9b62SBarry Smith        dimension of 2.
14782ea9b62SBarry Smith 
14882ea9b62SBarry Smith     Logically Collective on PetscViewer
14982ea9b62SBarry Smith 
15082ea9b62SBarry Smith   Input Parameters:
15182ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
15282ea9b62SBarry Smith -  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith   Options Database:
15582ea9b62SBarry Smith .  -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
15682ea9b62SBarry Smith 
15782ea9b62SBarry Smith 
15895452b02SPatrick Sanan   Notes:
15995452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
16082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
16182ea9b62SBarry Smith 
16282ea9b62SBarry Smith   Level: intermediate
16382ea9b62SBarry Smith 
16482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
16582ea9b62SBarry Smith 
16682ea9b62SBarry Smith @*/
16782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
16882ea9b62SBarry Smith {
16982ea9b62SBarry Smith   PetscErrorCode ierr;
17082ea9b62SBarry Smith 
17182ea9b62SBarry Smith   PetscFunctionBegin;
17282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
17382ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
17482ea9b62SBarry Smith   PetscFunctionReturn(0);
17582ea9b62SBarry Smith }
17682ea9b62SBarry Smith 
177df863907SAlex Fikl /*@
17882ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
17982ea9b62SBarry Smith        dimension of 2.
18082ea9b62SBarry Smith 
18182ea9b62SBarry Smith     Logically Collective on PetscViewer
18282ea9b62SBarry Smith 
18382ea9b62SBarry Smith   Input Parameter:
18482ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
18582ea9b62SBarry Smith 
18682ea9b62SBarry Smith   Output Parameter:
18782ea9b62SBarry Smith .  flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1
18882ea9b62SBarry Smith 
18995452b02SPatrick Sanan   Notes:
19095452b02SPatrick Sanan     Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof
19182ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
19282ea9b62SBarry Smith 
19382ea9b62SBarry Smith   Level: intermediate
19482ea9b62SBarry Smith 
19582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
19682ea9b62SBarry Smith 
19782ea9b62SBarry Smith @*/
19882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
19982ea9b62SBarry Smith {
20082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
20182ea9b62SBarry Smith 
20282ea9b62SBarry Smith   PetscFunctionBegin;
20382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
20482ea9b62SBarry Smith   *flg = hdf5->basedimension2;
20582ea9b62SBarry Smith   PetscFunctionReturn(0);
20682ea9b62SBarry Smith }
20782ea9b62SBarry Smith 
2089b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2099a0502c6SHåkon Strandenes {
2109a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2119a0502c6SHåkon Strandenes 
2129a0502c6SHåkon Strandenes   PetscFunctionBegin;
2139a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2149a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2159a0502c6SHåkon Strandenes }
2169a0502c6SHåkon Strandenes 
217df863907SAlex Fikl /*@
2189a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2199a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2209a0502c6SHåkon Strandenes 
2219a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes   Input Parameters:
2249a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2259a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2269a0502c6SHåkon Strandenes 
2279a0502c6SHåkon Strandenes   Options Database:
2289a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2299a0502c6SHåkon Strandenes 
2309a0502c6SHåkon Strandenes 
23195452b02SPatrick Sanan   Notes:
23295452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2339a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2349a0502c6SHåkon Strandenes 
2359a0502c6SHåkon Strandenes   Level: intermediate
2369a0502c6SHåkon Strandenes 
2379a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2389a0502c6SHåkon Strandenes           PetscReal
2399a0502c6SHåkon Strandenes 
2409a0502c6SHåkon Strandenes @*/
2419a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2429a0502c6SHåkon Strandenes {
2439a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2449a0502c6SHåkon Strandenes 
2459a0502c6SHåkon Strandenes   PetscFunctionBegin;
2469a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2479a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2489a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2499a0502c6SHåkon Strandenes }
2509a0502c6SHåkon Strandenes 
251df863907SAlex Fikl /*@
2529a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2539a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2549a0502c6SHåkon Strandenes 
2559a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2569a0502c6SHåkon Strandenes 
2579a0502c6SHåkon Strandenes   Input Parameter:
2589a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2599a0502c6SHåkon Strandenes 
2609a0502c6SHåkon Strandenes   Output Parameter:
2619a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2629a0502c6SHåkon Strandenes 
26395452b02SPatrick Sanan   Notes:
26495452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2659a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2669a0502c6SHåkon Strandenes 
2679a0502c6SHåkon Strandenes   Level: intermediate
2689a0502c6SHåkon Strandenes 
2699a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2709a0502c6SHåkon Strandenes           PetscReal
2719a0502c6SHåkon Strandenes 
2729a0502c6SHåkon Strandenes @*/
2739a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2749a0502c6SHåkon Strandenes {
2759a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2769a0502c6SHåkon Strandenes 
2779a0502c6SHåkon Strandenes   PetscFunctionBegin;
2789a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2799a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2809a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2819a0502c6SHåkon Strandenes }
2829a0502c6SHåkon Strandenes 
2839b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2845c6c1daeSBarry Smith {
2855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2865c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2875c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2885c6c1daeSBarry Smith #endif
2895c6c1daeSBarry Smith   hid_t             plist_id;
2905c6c1daeSBarry Smith   PetscErrorCode    ierr;
2915c6c1daeSBarry Smith 
2925c6c1daeSBarry Smith   PetscFunctionBegin;
293f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
294f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2955c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2965c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
297729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2985c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
299729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
3005c6c1daeSBarry Smith #endif
3015c6c1daeSBarry Smith   /* Create or open the file collectively */
3025c6c1daeSBarry Smith   switch (hdf5->btype) {
3035c6c1daeSBarry Smith   case FILE_MODE_READ:
304729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3055c6c1daeSBarry Smith     break;
3065c6c1daeSBarry Smith   case FILE_MODE_APPEND:
307729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
3085c6c1daeSBarry Smith     break;
3095c6c1daeSBarry Smith   case FILE_MODE_WRITE:
310729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
3115c6c1daeSBarry Smith     break;
3125c6c1daeSBarry Smith   default:
3135c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3145c6c1daeSBarry Smith   }
3155c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
316729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
3175c6c1daeSBarry Smith   PetscFunctionReturn(0);
3185c6c1daeSBarry Smith }
3195c6c1daeSBarry Smith 
320d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
321d1232d7fSToby Isaac {
322d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
323d1232d7fSToby Isaac 
324d1232d7fSToby Isaac   PetscFunctionBegin;
325d1232d7fSToby Isaac   *name = vhdf5->filename;
326d1232d7fSToby Isaac   PetscFunctionReturn(0);
327d1232d7fSToby Isaac }
328d1232d7fSToby Isaac 
329b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
330b723ab35SVaclav Hapla {
3315dc64a97SVaclav Hapla   /*
332b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
333b723ab35SVaclav Hapla   PetscErrorCode   ierr;
3345dc64a97SVaclav Hapla   */
335b723ab35SVaclav Hapla 
336b723ab35SVaclav Hapla   PetscFunctionBegin;
337b723ab35SVaclav Hapla   PetscFunctionReturn(0);
338b723ab35SVaclav Hapla }
339b723ab35SVaclav Hapla 
3408556b5ebSBarry Smith /*MC
3418556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
3428556b5ebSBarry Smith 
3438556b5ebSBarry Smith 
3448556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
3458556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
3468556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
3478556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
3488556b5ebSBarry Smith 
3491b266c99SBarry Smith   Level: beginner
3508556b5ebSBarry Smith M*/
351d1232d7fSToby Isaac 
3528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
3535c6c1daeSBarry Smith {
3545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
3555c6c1daeSBarry Smith   PetscErrorCode   ierr;
3565c6c1daeSBarry Smith 
3575c6c1daeSBarry Smith   PetscFunctionBegin;
358b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
3595c6c1daeSBarry Smith 
3605c6c1daeSBarry Smith   v->data                = (void*) hdf5;
3615c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
36282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
363b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
3641b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
3655c6c1daeSBarry Smith   v->ops->flush          = 0;
3665c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
3675c6c1daeSBarry Smith   hdf5->filename         = 0;
3685c6c1daeSBarry Smith   hdf5->timestep         = -1;
3690298fd71SBarry Smith   hdf5->groups           = NULL;
3705c6c1daeSBarry Smith 
3710b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
372d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
3730b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
3740b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
37582ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3769a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3775c6c1daeSBarry Smith   PetscFunctionReturn(0);
3785c6c1daeSBarry Smith }
3795c6c1daeSBarry Smith 
3805c6c1daeSBarry Smith /*@C
3815c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3825c6c1daeSBarry Smith 
383*d083f849SBarry Smith    Collective
3845c6c1daeSBarry Smith 
3855c6c1daeSBarry Smith    Input Parameters:
3865c6c1daeSBarry Smith +  comm - MPI communicator
3875c6c1daeSBarry Smith .  name - name of file
3885c6c1daeSBarry Smith -  type - type of file
3895c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3905c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3915c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3925c6c1daeSBarry Smith 
3935c6c1daeSBarry Smith    Output Parameter:
3945c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3955c6c1daeSBarry Smith 
39682ea9b62SBarry Smith   Options Database:
39782ea9b62SBarry Smith .  -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
3989a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
39982ea9b62SBarry Smith 
4005c6c1daeSBarry Smith    Level: beginner
4015c6c1daeSBarry Smith 
4025c6c1daeSBarry Smith    Note:
4035c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
4045c6c1daeSBarry Smith 
4055c6c1daeSBarry Smith 
4066a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
4079a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
408a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
4095c6c1daeSBarry Smith @*/
4105c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
4115c6c1daeSBarry Smith {
4125c6c1daeSBarry Smith   PetscErrorCode ierr;
4135c6c1daeSBarry Smith 
4145c6c1daeSBarry Smith   PetscFunctionBegin;
4155c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
4165c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
4175c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
4185c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
419b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
4205c6c1daeSBarry Smith   PetscFunctionReturn(0);
4215c6c1daeSBarry Smith }
4225c6c1daeSBarry Smith 
4235c6c1daeSBarry Smith /*@C
4245c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
4255c6c1daeSBarry Smith 
4265c6c1daeSBarry Smith   Not collective
4275c6c1daeSBarry Smith 
4285c6c1daeSBarry Smith   Input Parameter:
4295c6c1daeSBarry Smith . viewer - the PetscViewer
4305c6c1daeSBarry Smith 
4315c6c1daeSBarry Smith   Output Parameter:
4325c6c1daeSBarry Smith . file_id - The file id
4335c6c1daeSBarry Smith 
4345c6c1daeSBarry Smith   Level: intermediate
4355c6c1daeSBarry Smith 
4365c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
4375c6c1daeSBarry Smith @*/
4385c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
4395c6c1daeSBarry Smith {
4405c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4415c6c1daeSBarry Smith 
4425c6c1daeSBarry Smith   PetscFunctionBegin;
4435c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4445c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
4455c6c1daeSBarry Smith   PetscFunctionReturn(0);
4465c6c1daeSBarry Smith }
4475c6c1daeSBarry Smith 
4485c6c1daeSBarry Smith /*@C
4495c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
4505c6c1daeSBarry Smith 
4515c6c1daeSBarry Smith   Not collective
4525c6c1daeSBarry Smith 
4535c6c1daeSBarry Smith   Input Parameters:
4545c6c1daeSBarry Smith + viewer - the PetscViewer
4555c6c1daeSBarry Smith - name - The group name
4565c6c1daeSBarry Smith 
4575c6c1daeSBarry Smith   Level: intermediate
4585c6c1daeSBarry Smith 
459874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4605c6c1daeSBarry Smith @*/
4615c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
4625c6c1daeSBarry Smith {
4635c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4645c6c1daeSBarry Smith   GroupList        *groupNode;
4655c6c1daeSBarry Smith   PetscErrorCode   ierr;
4665c6c1daeSBarry Smith 
4675c6c1daeSBarry Smith   PetscFunctionBegin;
4685c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4695c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
47095dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
4715c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
472a297a907SKarl Rupp 
4735c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4745c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4755c6c1daeSBarry Smith   PetscFunctionReturn(0);
4765c6c1daeSBarry Smith }
4775c6c1daeSBarry Smith 
4783ef9c667SSatish Balay /*@
4795c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4805c6c1daeSBarry Smith 
4815c6c1daeSBarry Smith   Not collective
4825c6c1daeSBarry Smith 
4835c6c1daeSBarry Smith   Input Parameter:
4845c6c1daeSBarry Smith . viewer - the PetscViewer
4855c6c1daeSBarry Smith 
4865c6c1daeSBarry Smith   Level: intermediate
4875c6c1daeSBarry Smith 
488874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4895c6c1daeSBarry Smith @*/
4905c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4915c6c1daeSBarry Smith {
4925c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4935c6c1daeSBarry Smith   GroupList        *groupNode;
4945c6c1daeSBarry Smith   PetscErrorCode   ierr;
4955c6c1daeSBarry Smith 
4965c6c1daeSBarry Smith   PetscFunctionBegin;
4975c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
49882f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4995c6c1daeSBarry Smith   groupNode    = hdf5->groups;
5005c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
5015c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
5025c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
5035c6c1daeSBarry Smith   PetscFunctionReturn(0);
5045c6c1daeSBarry Smith }
5055c6c1daeSBarry Smith 
5065c6c1daeSBarry Smith /*@C
507874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
508874270d9SVaclav Hapla   If none has been assigned, returns NULL.
5095c6c1daeSBarry Smith 
5105c6c1daeSBarry Smith   Not collective
5115c6c1daeSBarry Smith 
5125c6c1daeSBarry Smith   Input Parameter:
5135c6c1daeSBarry Smith . viewer - the PetscViewer
5145c6c1daeSBarry Smith 
5155c6c1daeSBarry Smith   Output Parameter:
5165c6c1daeSBarry Smith . name - The group name
5175c6c1daeSBarry Smith 
5185c6c1daeSBarry Smith   Level: intermediate
5195c6c1daeSBarry Smith 
520874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
5215c6c1daeSBarry Smith @*/
5225c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
5235c6c1daeSBarry Smith {
5245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
5255c6c1daeSBarry Smith 
5265c6c1daeSBarry Smith   PetscFunctionBegin;
5275c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
528c959eef4SJed Brown   PetscValidPointer(name,2);
529a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
5300298fd71SBarry Smith   else *name = NULL;
5315c6c1daeSBarry Smith   PetscFunctionReturn(0);
5325c6c1daeSBarry Smith }
5335c6c1daeSBarry Smith 
5348c557b5aSVaclav Hapla /*@
535874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
536874270d9SVaclav Hapla   and return this group's ID and file ID.
537874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
538874270d9SVaclav Hapla 
539874270d9SVaclav Hapla   Not collective
540874270d9SVaclav Hapla 
541874270d9SVaclav Hapla   Input Parameter:
542874270d9SVaclav Hapla . viewer - the PetscViewer
543874270d9SVaclav Hapla 
544874270d9SVaclav Hapla   Output Parameter:
545874270d9SVaclav Hapla + fileId - The HDF5 file ID
546874270d9SVaclav Hapla - groupId - The HDF5 group ID
547874270d9SVaclav Hapla 
548e74616beSVaclav Hapla   Notes:
549e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
550e74616beSVaclav Hapla 
551874270d9SVaclav Hapla   Level: intermediate
552874270d9SVaclav Hapla 
553874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
554874270d9SVaclav Hapla @*/
55554dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
55654dbf706SVaclav Hapla {
55790e03537SVaclav Hapla   hid_t          file_id;
55890e03537SVaclav Hapla   H5O_type_t     type;
55954dbf706SVaclav Hapla   const char     *groupName = NULL;
560e74616beSVaclav Hapla   PetscBool      create;
56154dbf706SVaclav Hapla   PetscErrorCode ierr;
56254dbf706SVaclav Hapla 
56354dbf706SVaclav Hapla   PetscFunctionBegin;
564e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
56554dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
56654dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
567e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
56890e03537SVaclav Hapla   if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName);
56990e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
57054dbf706SVaclav Hapla   *fileId  = file_id;
57154dbf706SVaclav Hapla   PetscFunctionReturn(0);
57254dbf706SVaclav Hapla }
57354dbf706SVaclav Hapla 
5743ef9c667SSatish Balay /*@
5755c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
5765c6c1daeSBarry Smith 
5775c6c1daeSBarry Smith   Not collective
5785c6c1daeSBarry Smith 
5795c6c1daeSBarry Smith   Input Parameter:
5805c6c1daeSBarry Smith . viewer - the PetscViewer
5815c6c1daeSBarry Smith 
5825c6c1daeSBarry Smith   Level: intermediate
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
5855c6c1daeSBarry Smith @*/
5865c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
5875c6c1daeSBarry Smith {
5885c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5895c6c1daeSBarry Smith 
5905c6c1daeSBarry Smith   PetscFunctionBegin;
5915c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5925c6c1daeSBarry Smith   ++hdf5->timestep;
5935c6c1daeSBarry Smith   PetscFunctionReturn(0);
5945c6c1daeSBarry Smith }
5955c6c1daeSBarry Smith 
5963ef9c667SSatish Balay /*@
5975c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5985c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith   Not collective
6015c6c1daeSBarry Smith 
6025c6c1daeSBarry Smith   Input Parameters:
6035c6c1daeSBarry Smith + viewer - the PetscViewer
6045c6c1daeSBarry Smith - timestep - The timestep number
6055c6c1daeSBarry Smith 
6065c6c1daeSBarry Smith   Level: intermediate
6075c6c1daeSBarry Smith 
6085c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
6095c6c1daeSBarry Smith @*/
6105c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
6115c6c1daeSBarry Smith {
6125c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6135c6c1daeSBarry Smith 
6145c6c1daeSBarry Smith   PetscFunctionBegin;
6155c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6165c6c1daeSBarry Smith   hdf5->timestep = timestep;
6175c6c1daeSBarry Smith   PetscFunctionReturn(0);
6185c6c1daeSBarry Smith }
6195c6c1daeSBarry Smith 
6203ef9c667SSatish Balay /*@
6215c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
6225c6c1daeSBarry Smith 
6235c6c1daeSBarry Smith   Not collective
6245c6c1daeSBarry Smith 
6255c6c1daeSBarry Smith   Input Parameter:
6265c6c1daeSBarry Smith . viewer - the PetscViewer
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith   Output Parameter:
6295c6c1daeSBarry Smith . timestep - The timestep number
6305c6c1daeSBarry Smith 
6315c6c1daeSBarry Smith   Level: intermediate
6325c6c1daeSBarry Smith 
6335c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
6345c6c1daeSBarry Smith @*/
6355c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
6365c6c1daeSBarry Smith {
6375c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6385c6c1daeSBarry Smith 
6395c6c1daeSBarry Smith   PetscFunctionBegin;
6405c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6415c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
6425c6c1daeSBarry Smith   *timestep = hdf5->timestep;
6435c6c1daeSBarry Smith   PetscFunctionReturn(0);
6445c6c1daeSBarry Smith }
6455c6c1daeSBarry Smith 
64636bce6e8SMatthew G. Knepley /*@C
64736bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
64836bce6e8SMatthew G. Knepley 
64936bce6e8SMatthew G. Knepley   Not collective
65036bce6e8SMatthew G. Knepley 
65136bce6e8SMatthew G. Knepley   Input Parameter:
65236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
65336bce6e8SMatthew G. Knepley 
65436bce6e8SMatthew G. Knepley   Output Parameter:
65536bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
65636bce6e8SMatthew G. Knepley 
65736bce6e8SMatthew G. Knepley   Level: advanced
65836bce6e8SMatthew G. Knepley 
65936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
66036bce6e8SMatthew G. Knepley @*/
66136bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
66236bce6e8SMatthew G. Knepley {
66336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
66436bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
66536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
66636bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
66736bce6e8SMatthew G. Knepley #else
66836bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
66936bce6e8SMatthew G. Knepley #endif
67036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
67136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
67236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
673c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
674de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
67536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
67636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
67736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
6787e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
67936bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
68036bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
68136bce6e8SMatthew G. Knepley }
68236bce6e8SMatthew G. Knepley 
68336bce6e8SMatthew G. Knepley /*@C
68436bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
68536bce6e8SMatthew G. Knepley 
68636bce6e8SMatthew G. Knepley   Not collective
68736bce6e8SMatthew G. Knepley 
68836bce6e8SMatthew G. Knepley   Input Parameter:
68936bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
69036bce6e8SMatthew G. Knepley 
69136bce6e8SMatthew G. Knepley   Output Parameter:
69236bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
69336bce6e8SMatthew G. Knepley 
69436bce6e8SMatthew G. Knepley   Level: advanced
69536bce6e8SMatthew G. Knepley 
69636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
69736bce6e8SMatthew G. Knepley @*/
69836bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
69936bce6e8SMatthew G. Knepley {
70036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
70136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
70236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
70336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
70436bce6e8SMatthew G. Knepley #else
70536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
70636bce6e8SMatthew G. Knepley #endif
70736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
70836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
70936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
71036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
71136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
71236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
7137e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
71436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
71536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
71636bce6e8SMatthew G. Knepley }
71736bce6e8SMatthew G. Knepley 
718df863907SAlex Fikl /*@C
719b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
72036bce6e8SMatthew G. Knepley 
72136bce6e8SMatthew G. Knepley   Input Parameters:
72236bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
72357d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
72436bce6e8SMatthew G. Knepley . name   - The attribute name
72536bce6e8SMatthew G. Knepley . datatype - The attribute type
72636bce6e8SMatthew G. Knepley - value    - The attribute value
72736bce6e8SMatthew G. Knepley 
72836bce6e8SMatthew G. Knepley   Level: advanced
72936bce6e8SMatthew G. Knepley 
730577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
73136bce6e8SMatthew G. Knepley @*/
73257d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
73336bce6e8SMatthew G. Knepley {
73457d22f7dSVaclav Hapla   char           *parent;
73560568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
736080f144cSVaclav Hapla   PetscBool      has;
73736bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
73836bce6e8SMatthew G. Knepley 
73936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
7405cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
74157d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
742c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
743b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
74457d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
745bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
746b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
74736bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
7487e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
7497e97c476SMatthew G. Knepley     size_t len;
7507e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
751729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
7527e97c476SMatthew G. Knepley   }
75336bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
754729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
75560568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
756080f144cSVaclav Hapla   if (has) {
757080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
758080f144cSVaclav Hapla   } else {
75960568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
760080f144cSVaclav Hapla   }
761729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
762729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
763729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
76460568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
765729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
76657d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
76736bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
76836bce6e8SMatthew G. Knepley }
76936bce6e8SMatthew G. Knepley 
770df863907SAlex Fikl /*@C
771577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
772577e0e04SVaclav Hapla 
773577e0e04SVaclav Hapla   Input Parameters:
774577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
775577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
776577e0e04SVaclav Hapla . name     - The attribute name
777577e0e04SVaclav Hapla . datatype - The attribute type
778577e0e04SVaclav Hapla - value    - The attribute value
779577e0e04SVaclav Hapla 
780577e0e04SVaclav Hapla   Notes:
781577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
782577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
783577e0e04SVaclav Hapla 
784577e0e04SVaclav Hapla   Level: advanced
785577e0e04SVaclav Hapla 
786577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
787577e0e04SVaclav Hapla @*/
788577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
789577e0e04SVaclav Hapla {
790577e0e04SVaclav Hapla   PetscErrorCode ierr;
791577e0e04SVaclav Hapla 
792577e0e04SVaclav Hapla   PetscFunctionBegin;
793577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
794577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
795577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
796b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
797577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
798577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
799577e0e04SVaclav Hapla   PetscFunctionReturn(0);
800577e0e04SVaclav Hapla }
801577e0e04SVaclav Hapla 
802577e0e04SVaclav Hapla /*@C
803b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
804ad0c61feSMatthew G. Knepley 
805ad0c61feSMatthew G. Knepley   Input Parameters:
806ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
80757d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
808ad0c61feSMatthew G. Knepley . name   - The attribute name
809ad0c61feSMatthew G. Knepley - datatype - The attribute type
810ad0c61feSMatthew G. Knepley 
811ad0c61feSMatthew G. Knepley   Output Parameter:
812ad0c61feSMatthew G. Knepley . value    - The attribute value
813ad0c61feSMatthew G. Knepley 
814ad0c61feSMatthew G. Knepley   Level: advanced
815ad0c61feSMatthew G. Knepley 
816577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
817ad0c61feSMatthew G. Knepley @*/
81857d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
819ad0c61feSMatthew G. Knepley {
82057d22f7dSVaclav Hapla   char           *parent;
821f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
822969748fdSVaclav Hapla   PetscBool      has;
823ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
824ad0c61feSMatthew G. Knepley 
825ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
8265cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
82757d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
828c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
829b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
83057d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
831969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
832969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
833969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
834ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
835ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
83660568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
83760568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
838f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
839f0b58479SMatthew G. Knepley     size_t len;
840e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
84115b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
842f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
843f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
844f0b58479SMatthew G. Knepley   }
84570efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
846729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
847e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
848e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
84957d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
850ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
851ad0c61feSMatthew G. Knepley }
852ad0c61feSMatthew G. Knepley 
853577e0e04SVaclav Hapla /*@C
854577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
855577e0e04SVaclav Hapla 
856577e0e04SVaclav Hapla   Input Parameters:
857577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
858577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
859577e0e04SVaclav Hapla . name     - The attribute name
860577e0e04SVaclav Hapla - datatype - The attribute type
861577e0e04SVaclav Hapla 
862577e0e04SVaclav Hapla   Output Parameter:
863577e0e04SVaclav Hapla . value    - The attribute value
864577e0e04SVaclav Hapla 
865577e0e04SVaclav Hapla   Notes:
866577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
867577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
868577e0e04SVaclav Hapla 
869577e0e04SVaclav Hapla   Level: advanced
870577e0e04SVaclav Hapla 
871577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
872577e0e04SVaclav Hapla @*/
873577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
874577e0e04SVaclav Hapla {
875577e0e04SVaclav Hapla   PetscErrorCode ierr;
876577e0e04SVaclav Hapla 
877577e0e04SVaclav Hapla   PetscFunctionBegin;
878577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
879577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
880577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
881b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
882577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
883577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
884577e0e04SVaclav Hapla   PetscFunctionReturn(0);
885577e0e04SVaclav Hapla }
886577e0e04SVaclav Hapla 
887a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
888a75c3fd4SVaclav Hapla {
889a75c3fd4SVaclav Hapla   htri_t exists;
890a75c3fd4SVaclav Hapla   hid_t group;
891a75c3fd4SVaclav Hapla 
892a75c3fd4SVaclav Hapla   PetscFunctionBegin;
893a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
894a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
895a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
896a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
897a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
898a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
899a75c3fd4SVaclav Hapla   }
900a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
901a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
902a75c3fd4SVaclav Hapla }
903a75c3fd4SVaclav Hapla 
904bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
9055cdeb108SMatthew G. Knepley {
90690e03537SVaclav Hapla   const char     rootGroupName[] = "/";
9075cdeb108SMatthew G. Knepley   hid_t          h5;
908e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
90915b861d2SVaclav Hapla   PetscInt       i;
91015b861d2SVaclav Hapla   int            n;
91185688ae2SVaclav Hapla   char           **hierarchy;
91285688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
9135cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
9145cdeb108SMatthew G. Knepley 
9155cdeb108SMatthew G. Knepley   PetscFunctionBegin;
9165cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
91790e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
91890e03537SVaclav Hapla   else name = rootGroupName;
919ccd66a83SVaclav Hapla   if (has) {
92056cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
9215cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
922ccd66a83SVaclav Hapla   }
923ccd66a83SVaclav Hapla   if (otype) {
924ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
92556cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
926ccd66a83SVaclav Hapla   }
927c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
92885688ae2SVaclav Hapla 
92985688ae2SVaclav Hapla   /*
93085688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
93185688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
93285688ae2SVaclav Hapla      1) whether it's a valid link
93385688ae2SVaclav Hapla      2) whether this link resolves to an object
93485688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
93585688ae2SVaclav Hapla   */
93685688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
93785688ae2SVaclav Hapla   if (!n) {
93885688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
939ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
940ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
94185688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
94285688ae2SVaclav Hapla     PetscFunctionReturn(0);
94385688ae2SVaclav Hapla   }
94485688ae2SVaclav Hapla   for (i=0; i<n; i++) {
94585688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
94685688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
947a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
948a75c3fd4SVaclav Hapla     if (!exists) break;
94985688ae2SVaclav Hapla   }
95085688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
95185688ae2SVaclav Hapla 
95285688ae2SVaclav Hapla   /* If the object exists, get its type */
953ccd66a83SVaclav Hapla   if (exists && otype) {
9545cdeb108SMatthew G. Knepley     H5O_info_t info;
9555cdeb108SMatthew G. Knepley 
95676276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
95704633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
95856cc0592SVaclav Hapla     *otype = info.type;
9595cdeb108SMatthew G. Knepley   }
960ccd66a83SVaclav Hapla   if (has) *has = exists;
9615cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
9625cdeb108SMatthew G. Knepley }
9635cdeb108SMatthew G. Knepley 
964ecce1506SVaclav Hapla /*@
9650a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
9660a9f38efSVaclav Hapla 
9670a9f38efSVaclav Hapla   Input Parameters:
9680a9f38efSVaclav Hapla . viewer - The HDF5 viewer
9690a9f38efSVaclav Hapla 
9700a9f38efSVaclav Hapla   Output Parameter:
9710a9f38efSVaclav Hapla . has    - Flag for group existence
9720a9f38efSVaclav Hapla 
9730a9f38efSVaclav Hapla   Notes:
9740a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
9750a9f38efSVaclav Hapla 
9760a9f38efSVaclav Hapla   Level: advanced
9770a9f38efSVaclav Hapla 
9780a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
9790a9f38efSVaclav Hapla @*/
9800a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
9810a9f38efSVaclav Hapla {
9820a9f38efSVaclav Hapla   H5O_type_t type;
9830a9f38efSVaclav Hapla   const char *name;
9840a9f38efSVaclav Hapla   PetscErrorCode ierr;
9850a9f38efSVaclav Hapla 
9860a9f38efSVaclav Hapla   PetscFunctionBegin;
9870a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
9880a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
9890a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
9900a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
9910a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
9920a9f38efSVaclav Hapla   PetscFunctionReturn(0);
9930a9f38efSVaclav Hapla }
9940a9f38efSVaclav Hapla 
9950a9f38efSVaclav Hapla /*@
996e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
997ecce1506SVaclav Hapla 
998ecce1506SVaclav Hapla   Input Parameters:
999ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1000ecce1506SVaclav Hapla - obj    - The named object
1001ecce1506SVaclav Hapla 
1002ecce1506SVaclav Hapla   Output Parameter:
1003ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1004ecce1506SVaclav Hapla 
1005e3f143f7SVaclav Hapla   Notes:
1006e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1007e3f143f7SVaclav Hapla 
1008ecce1506SVaclav Hapla   Level: advanced
1009ecce1506SVaclav Hapla 
1010e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1011ecce1506SVaclav Hapla @*/
1012ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1013ecce1506SVaclav Hapla {
101456cc0592SVaclav Hapla   H5O_type_t type;
10156c132bc1SVaclav Hapla   char *path;
1016ecce1506SVaclav Hapla   PetscErrorCode ierr;
1017ecce1506SVaclav Hapla 
1018ecce1506SVaclav Hapla   PetscFunctionBegin;
1019c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1020c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1021c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1022ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1023e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
10246c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
10256c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
102656cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
10276c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1028ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1029ecce1506SVaclav Hapla }
1030ecce1506SVaclav Hapla 
1031df863907SAlex Fikl /*@C
1032ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1033e7bdbf8eSMatthew G. Knepley 
1034e7bdbf8eSMatthew G. Knepley   Input Parameters:
1035e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
103610e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1037e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1038e7bdbf8eSMatthew G. Knepley 
1039e7bdbf8eSMatthew G. Knepley   Output Parameter:
1040e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1041e7bdbf8eSMatthew G. Knepley 
1042e7bdbf8eSMatthew G. Knepley   Level: advanced
1043e7bdbf8eSMatthew G. Knepley 
1044577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1045e7bdbf8eSMatthew G. Knepley @*/
104610e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1047e7bdbf8eSMatthew G. Knepley {
104810e69e8fSVaclav Hapla   char           *parent;
1049e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1050e7bdbf8eSMatthew G. Knepley 
1051e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
10525cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
105310e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1054c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1055c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
105610e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1057bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
105810e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
105910e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
106006db490cSVaclav Hapla   PetscFunctionReturn(0);
106106db490cSVaclav Hapla }
106206db490cSVaclav Hapla 
1063577e0e04SVaclav Hapla /*@C
1064577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1065577e0e04SVaclav Hapla 
1066577e0e04SVaclav Hapla   Input Parameters:
1067577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1068577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1069577e0e04SVaclav Hapla - name   - The attribute name
1070577e0e04SVaclav Hapla 
1071577e0e04SVaclav Hapla   Output Parameter:
1072577e0e04SVaclav Hapla . has    - Flag for attribute existence
1073577e0e04SVaclav Hapla 
1074577e0e04SVaclav Hapla   Notes:
1075577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1076577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1077577e0e04SVaclav Hapla 
1078577e0e04SVaclav Hapla   Level: advanced
1079577e0e04SVaclav Hapla 
1080577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1081577e0e04SVaclav Hapla @*/
1082577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1083577e0e04SVaclav Hapla {
1084577e0e04SVaclav Hapla   PetscErrorCode ierr;
1085577e0e04SVaclav Hapla 
1086577e0e04SVaclav Hapla   PetscFunctionBegin;
1087577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1088577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1089577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1090577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1091577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1092577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1093577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1094577e0e04SVaclav Hapla }
1095577e0e04SVaclav Hapla 
109606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
109706db490cSVaclav Hapla {
109806db490cSVaclav Hapla   hid_t          h5;
109906db490cSVaclav Hapla   htri_t         hhas;
110006db490cSVaclav Hapla   PetscErrorCode ierr;
110106db490cSVaclav Hapla 
110206db490cSVaclav Hapla   PetscFunctionBegin;
110306db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
11042f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
11055cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1106e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1107e7bdbf8eSMatthew G. Knepley }
1108e7bdbf8eSMatthew G. Knepley 
1109a75e6a4aSMatthew G. Knepley /*
1110a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1111a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1112a75e6a4aSMatthew G. Knepley */
1113d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1114a75e6a4aSMatthew G. Knepley 
1115a75e6a4aSMatthew G. Knepley /*@C
1116a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1117a75e6a4aSMatthew G. Knepley 
1118*d083f849SBarry Smith   Collective
1119a75e6a4aSMatthew G. Knepley 
1120a75e6a4aSMatthew G. Knepley   Input Parameter:
1121a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1122a75e6a4aSMatthew G. Knepley 
1123a75e6a4aSMatthew G. Knepley   Level: intermediate
1124a75e6a4aSMatthew G. Knepley 
1125a75e6a4aSMatthew G. Knepley   Options Database Keys:
1126a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1127a75e6a4aSMatthew G. Knepley 
1128a75e6a4aSMatthew G. Knepley   Environmental variables:
1129a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1130a75e6a4aSMatthew G. Knepley 
1131a75e6a4aSMatthew G. Knepley   Notes:
1132a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1133a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1134a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1135a75e6a4aSMatthew G. Knepley 
1136a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1137a75e6a4aSMatthew G. Knepley @*/
1138a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1139a75e6a4aSMatthew G. Knepley {
1140a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1141a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1142a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1143a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1144a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1145a75e6a4aSMatthew G. Knepley 
1146a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1147a75e6a4aSMatthew G. Knepley   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1148a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
114912801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1150a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1151a75e6a4aSMatthew G. Knepley   }
115247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1153a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1154a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1155a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1156a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1157a75e6a4aSMatthew G. Knepley     if (!flg) {
1158a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1159a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1160a75e6a4aSMatthew G. Knepley     }
1161a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1162a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1163a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1164a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
116547435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1166a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1167a75e6a4aSMatthew G. Knepley   }
1168a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
1169a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1170a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1171a75e6a4aSMatthew G. Knepley }
1172