xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision b663d2d2a73b43d1022aace8e05429bc2359e7c0)
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 
695c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
705c6c1daeSBarry Smith {
715c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
725c6c1daeSBarry Smith   PetscErrorCode   ierr;
735c6c1daeSBarry Smith 
745c6c1daeSBarry Smith   PetscFunctionBegin;
755c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
76729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
775c6c1daeSBarry Smith   PetscFunctionReturn(0);
785c6c1daeSBarry Smith }
795c6c1daeSBarry Smith 
809b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
815c6c1daeSBarry Smith {
825c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
835c6c1daeSBarry Smith   PetscErrorCode   ierr;
845c6c1daeSBarry Smith 
855c6c1daeSBarry Smith   PetscFunctionBegin;
865c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
875c6c1daeSBarry Smith   while (hdf5->groups) {
885c6c1daeSBarry Smith     GroupList *tmp = hdf5->groups->next;
895c6c1daeSBarry Smith 
905c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
915c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
925c6c1daeSBarry Smith     hdf5->groups = tmp;
935c6c1daeSBarry Smith   }
945c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
950b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
96d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
970b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
98058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
99058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1005c6c1daeSBarry Smith   PetscFunctionReturn(0);
1015c6c1daeSBarry Smith }
1025c6c1daeSBarry Smith 
1039b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1045c6c1daeSBarry Smith {
1055c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1065c6c1daeSBarry Smith 
1075c6c1daeSBarry Smith   PetscFunctionBegin;
1085c6c1daeSBarry Smith   hdf5->btype = type;
1095c6c1daeSBarry Smith   PetscFunctionReturn(0);
1105c6c1daeSBarry Smith }
1115c6c1daeSBarry Smith 
1120b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1130b62783dSVaclav Hapla {
1140b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1150b62783dSVaclav Hapla 
1160b62783dSVaclav Hapla   PetscFunctionBegin;
1170b62783dSVaclav Hapla   *type = hdf5->btype;
1180b62783dSVaclav Hapla   PetscFunctionReturn(0);
1190b62783dSVaclav Hapla }
1200b62783dSVaclav Hapla 
1219b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
12282ea9b62SBarry Smith {
12382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
12482ea9b62SBarry Smith 
12582ea9b62SBarry Smith   PetscFunctionBegin;
12682ea9b62SBarry Smith   hdf5->basedimension2 = flg;
12782ea9b62SBarry Smith   PetscFunctionReturn(0);
12882ea9b62SBarry Smith }
12982ea9b62SBarry Smith 
130df863907SAlex Fikl /*@
13182ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13282ea9b62SBarry Smith        dimension of 2.
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith     Logically Collective on PetscViewer
13582ea9b62SBarry Smith 
13682ea9b62SBarry Smith   Input Parameters:
13782ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
13882ea9b62SBarry 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
13982ea9b62SBarry Smith 
14082ea9b62SBarry Smith   Options Database:
14182ea9b62SBarry 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
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith 
14495452b02SPatrick Sanan   Notes:
14595452b02SPatrick 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
14682ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
14782ea9b62SBarry Smith 
14882ea9b62SBarry Smith   Level: intermediate
14982ea9b62SBarry Smith 
15082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith @*/
15382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
15482ea9b62SBarry Smith {
15582ea9b62SBarry Smith   PetscErrorCode ierr;
15682ea9b62SBarry Smith 
15782ea9b62SBarry Smith   PetscFunctionBegin;
15882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
15982ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
16082ea9b62SBarry Smith   PetscFunctionReturn(0);
16182ea9b62SBarry Smith }
16282ea9b62SBarry Smith 
163df863907SAlex Fikl /*@
16482ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
16582ea9b62SBarry Smith        dimension of 2.
16682ea9b62SBarry Smith 
16782ea9b62SBarry Smith     Logically Collective on PetscViewer
16882ea9b62SBarry Smith 
16982ea9b62SBarry Smith   Input Parameter:
17082ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
17182ea9b62SBarry Smith 
17282ea9b62SBarry Smith   Output Parameter:
17382ea9b62SBarry 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
17482ea9b62SBarry Smith 
17595452b02SPatrick Sanan   Notes:
17695452b02SPatrick 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
17782ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
17882ea9b62SBarry Smith 
17982ea9b62SBarry Smith   Level: intermediate
18082ea9b62SBarry Smith 
18182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
18282ea9b62SBarry Smith 
18382ea9b62SBarry Smith @*/
18482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
18582ea9b62SBarry Smith {
18682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
18782ea9b62SBarry Smith 
18882ea9b62SBarry Smith   PetscFunctionBegin;
18982ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
19082ea9b62SBarry Smith   *flg = hdf5->basedimension2;
19182ea9b62SBarry Smith   PetscFunctionReturn(0);
19282ea9b62SBarry Smith }
19382ea9b62SBarry Smith 
1949b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1959a0502c6SHåkon Strandenes {
1969a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1979a0502c6SHåkon Strandenes 
1989a0502c6SHåkon Strandenes   PetscFunctionBegin;
1999a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2009a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2019a0502c6SHåkon Strandenes }
2029a0502c6SHåkon Strandenes 
203df863907SAlex Fikl /*@
2049a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2059a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2069a0502c6SHåkon Strandenes 
2079a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2089a0502c6SHåkon Strandenes 
2099a0502c6SHåkon Strandenes   Input Parameters:
2109a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2119a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes   Options Database:
2149a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes 
21795452b02SPatrick Sanan   Notes:
21895452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2199a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2209a0502c6SHåkon Strandenes 
2219a0502c6SHåkon Strandenes   Level: intermediate
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2249a0502c6SHåkon Strandenes           PetscReal
2259a0502c6SHåkon Strandenes 
2269a0502c6SHåkon Strandenes @*/
2279a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2289a0502c6SHåkon Strandenes {
2299a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2309a0502c6SHåkon Strandenes 
2319a0502c6SHåkon Strandenes   PetscFunctionBegin;
2329a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2339a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2349a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2359a0502c6SHåkon Strandenes }
2369a0502c6SHåkon Strandenes 
237df863907SAlex Fikl /*@
2389a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2399a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2409a0502c6SHåkon Strandenes 
2419a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2429a0502c6SHåkon Strandenes 
2439a0502c6SHåkon Strandenes   Input Parameter:
2449a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2459a0502c6SHåkon Strandenes 
2469a0502c6SHåkon Strandenes   Output Parameter:
2479a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2489a0502c6SHåkon Strandenes 
24995452b02SPatrick Sanan   Notes:
25095452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2519a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2529a0502c6SHåkon Strandenes 
2539a0502c6SHåkon Strandenes   Level: intermediate
2549a0502c6SHåkon Strandenes 
2559a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2569a0502c6SHåkon Strandenes           PetscReal
2579a0502c6SHåkon Strandenes 
2589a0502c6SHåkon Strandenes @*/
2599a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2609a0502c6SHåkon Strandenes {
2619a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2629a0502c6SHåkon Strandenes 
2639a0502c6SHåkon Strandenes   PetscFunctionBegin;
2649a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2659a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2669a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2679a0502c6SHåkon Strandenes }
2689a0502c6SHåkon Strandenes 
2699b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
2705c6c1daeSBarry Smith {
2715c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2725c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
2735c6c1daeSBarry Smith   MPI_Info          info = MPI_INFO_NULL;
2745c6c1daeSBarry Smith #endif
2755c6c1daeSBarry Smith   hid_t             plist_id;
2765c6c1daeSBarry Smith   PetscErrorCode    ierr;
2775c6c1daeSBarry Smith 
2785c6c1daeSBarry Smith   PetscFunctionBegin;
279f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
280f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
2815c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
2825c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
283729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
2845c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO)
285729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info));
2865c6c1daeSBarry Smith #endif
2875c6c1daeSBarry Smith   /* Create or open the file collectively */
2885c6c1daeSBarry Smith   switch (hdf5->btype) {
2895c6c1daeSBarry Smith   case FILE_MODE_READ:
290729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
2915c6c1daeSBarry Smith     break;
2925c6c1daeSBarry Smith   case FILE_MODE_APPEND:
293729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
2945c6c1daeSBarry Smith     break;
2955c6c1daeSBarry Smith   case FILE_MODE_WRITE:
296729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
2975c6c1daeSBarry Smith     break;
2985c6c1daeSBarry Smith   default:
2995c6c1daeSBarry Smith     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3005c6c1daeSBarry Smith   }
3015c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
302729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
3035c6c1daeSBarry Smith   PetscFunctionReturn(0);
3045c6c1daeSBarry Smith }
3055c6c1daeSBarry Smith 
306d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
307d1232d7fSToby Isaac {
308d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
309d1232d7fSToby Isaac 
310d1232d7fSToby Isaac   PetscFunctionBegin;
311d1232d7fSToby Isaac   *name = vhdf5->filename;
312d1232d7fSToby Isaac   PetscFunctionReturn(0);
313d1232d7fSToby Isaac }
314d1232d7fSToby Isaac 
315b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
316b723ab35SVaclav Hapla {
3175dc64a97SVaclav Hapla   /*
318b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
319b723ab35SVaclav Hapla   PetscErrorCode   ierr;
3205dc64a97SVaclav Hapla   */
321b723ab35SVaclav Hapla 
322b723ab35SVaclav Hapla   PetscFunctionBegin;
323b723ab35SVaclav Hapla   PetscFunctionReturn(0);
324b723ab35SVaclav Hapla }
325b723ab35SVaclav Hapla 
3268556b5ebSBarry Smith /*MC
3278556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
3288556b5ebSBarry Smith 
3298556b5ebSBarry Smith 
3308556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
3318556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
3328556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
3338556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
3348556b5ebSBarry Smith 
3351b266c99SBarry Smith   Level: beginner
3368556b5ebSBarry Smith M*/
337d1232d7fSToby Isaac 
3388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
3395c6c1daeSBarry Smith {
3405c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
3415c6c1daeSBarry Smith   PetscErrorCode   ierr;
3425c6c1daeSBarry Smith 
3435c6c1daeSBarry Smith   PetscFunctionBegin;
344b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
3455c6c1daeSBarry Smith 
3465c6c1daeSBarry Smith   v->data                = (void*) hdf5;
3475c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
34882ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
349b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
3505c6c1daeSBarry Smith   v->ops->flush          = 0;
3515c6c1daeSBarry Smith   hdf5->btype            = (PetscFileMode) -1;
3525c6c1daeSBarry Smith   hdf5->filename         = 0;
3535c6c1daeSBarry Smith   hdf5->timestep         = -1;
3540298fd71SBarry Smith   hdf5->groups           = NULL;
3555c6c1daeSBarry Smith 
3560b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
357d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
3580b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
3590b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
36082ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
3619a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
3625c6c1daeSBarry Smith   PetscFunctionReturn(0);
3635c6c1daeSBarry Smith }
3645c6c1daeSBarry Smith 
3655c6c1daeSBarry Smith /*@C
3665c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
3675c6c1daeSBarry Smith 
3685c6c1daeSBarry Smith    Collective on MPI_Comm
3695c6c1daeSBarry Smith 
3705c6c1daeSBarry Smith    Input Parameters:
3715c6c1daeSBarry Smith +  comm - MPI communicator
3725c6c1daeSBarry Smith .  name - name of file
3735c6c1daeSBarry Smith -  type - type of file
3745c6c1daeSBarry Smith $    FILE_MODE_WRITE - create new file for binary output
3755c6c1daeSBarry Smith $    FILE_MODE_READ - open existing file for binary input
3765c6c1daeSBarry Smith $    FILE_MODE_APPEND - open existing file for binary output
3775c6c1daeSBarry Smith 
3785c6c1daeSBarry Smith    Output Parameter:
3795c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
3805c6c1daeSBarry Smith 
38182ea9b62SBarry Smith   Options Database:
38282ea9b62SBarry 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
3839a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
38482ea9b62SBarry Smith 
3855c6c1daeSBarry Smith    Level: beginner
3865c6c1daeSBarry Smith 
3875c6c1daeSBarry Smith    Note:
3885c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
3895c6c1daeSBarry Smith 
3905c6c1daeSBarry Smith    Concepts: HDF5 files
3915c6c1daeSBarry Smith    Concepts: PetscViewerHDF5^creating
3925c6c1daeSBarry Smith 
3936a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
3949a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
395a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
3965c6c1daeSBarry Smith @*/
3975c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
3985c6c1daeSBarry Smith {
3995c6c1daeSBarry Smith   PetscErrorCode ierr;
4005c6c1daeSBarry Smith 
4015c6c1daeSBarry Smith   PetscFunctionBegin;
4025c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
4035c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
4045c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
4055c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
406*b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
4075c6c1daeSBarry Smith   PetscFunctionReturn(0);
4085c6c1daeSBarry Smith }
4095c6c1daeSBarry Smith 
4105c6c1daeSBarry Smith /*@C
4115c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
4125c6c1daeSBarry Smith 
4135c6c1daeSBarry Smith   Not collective
4145c6c1daeSBarry Smith 
4155c6c1daeSBarry Smith   Input Parameter:
4165c6c1daeSBarry Smith . viewer - the PetscViewer
4175c6c1daeSBarry Smith 
4185c6c1daeSBarry Smith   Output Parameter:
4195c6c1daeSBarry Smith . file_id - The file id
4205c6c1daeSBarry Smith 
4215c6c1daeSBarry Smith   Level: intermediate
4225c6c1daeSBarry Smith 
4235c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
4245c6c1daeSBarry Smith @*/
4255c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
4265c6c1daeSBarry Smith {
4275c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4285c6c1daeSBarry Smith 
4295c6c1daeSBarry Smith   PetscFunctionBegin;
4305c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4315c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
4325c6c1daeSBarry Smith   PetscFunctionReturn(0);
4335c6c1daeSBarry Smith }
4345c6c1daeSBarry Smith 
4355c6c1daeSBarry Smith /*@C
4365c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
4375c6c1daeSBarry Smith 
4385c6c1daeSBarry Smith   Not collective
4395c6c1daeSBarry Smith 
4405c6c1daeSBarry Smith   Input Parameters:
4415c6c1daeSBarry Smith + viewer - the PetscViewer
4425c6c1daeSBarry Smith - name - The group name
4435c6c1daeSBarry Smith 
4445c6c1daeSBarry Smith   Level: intermediate
4455c6c1daeSBarry Smith 
446874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4475c6c1daeSBarry Smith @*/
4485c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name)
4495c6c1daeSBarry Smith {
4505c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4515c6c1daeSBarry Smith   GroupList        *groupNode;
4525c6c1daeSBarry Smith   PetscErrorCode   ierr;
4535c6c1daeSBarry Smith 
4545c6c1daeSBarry Smith   PetscFunctionBegin;
4555c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
4565c6c1daeSBarry Smith   PetscValidCharPointer(name,2);
45795dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
4585c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
459a297a907SKarl Rupp 
4605c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
4615c6c1daeSBarry Smith   hdf5->groups    = groupNode;
4625c6c1daeSBarry Smith   PetscFunctionReturn(0);
4635c6c1daeSBarry Smith }
4645c6c1daeSBarry Smith 
4653ef9c667SSatish Balay /*@
4665c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
4675c6c1daeSBarry Smith 
4685c6c1daeSBarry Smith   Not collective
4695c6c1daeSBarry Smith 
4705c6c1daeSBarry Smith   Input Parameter:
4715c6c1daeSBarry Smith . viewer - the PetscViewer
4725c6c1daeSBarry Smith 
4735c6c1daeSBarry Smith   Level: intermediate
4745c6c1daeSBarry Smith 
475874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
4765c6c1daeSBarry Smith @*/
4775c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
4785c6c1daeSBarry Smith {
4795c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4805c6c1daeSBarry Smith   GroupList        *groupNode;
4815c6c1daeSBarry Smith   PetscErrorCode   ierr;
4825c6c1daeSBarry Smith 
4835c6c1daeSBarry Smith   PetscFunctionBegin;
4845c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
48582f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
4865c6c1daeSBarry Smith   groupNode    = hdf5->groups;
4875c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
4885c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
4895c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
4905c6c1daeSBarry Smith   PetscFunctionReturn(0);
4915c6c1daeSBarry Smith }
4925c6c1daeSBarry Smith 
4935c6c1daeSBarry Smith /*@C
494874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
495874270d9SVaclav Hapla   If none has been assigned, returns NULL.
4965c6c1daeSBarry Smith 
4975c6c1daeSBarry Smith   Not collective
4985c6c1daeSBarry Smith 
4995c6c1daeSBarry Smith   Input Parameter:
5005c6c1daeSBarry Smith . viewer - the PetscViewer
5015c6c1daeSBarry Smith 
5025c6c1daeSBarry Smith   Output Parameter:
5035c6c1daeSBarry Smith . name - The group name
5045c6c1daeSBarry Smith 
5055c6c1daeSBarry Smith   Level: intermediate
5065c6c1daeSBarry Smith 
507874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
5085c6c1daeSBarry Smith @*/
5095c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name)
5105c6c1daeSBarry Smith {
5115c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
5125c6c1daeSBarry Smith 
5135c6c1daeSBarry Smith   PetscFunctionBegin;
5145c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
515c959eef4SJed Brown   PetscValidPointer(name,2);
516a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
5170298fd71SBarry Smith   else *name = NULL;
5185c6c1daeSBarry Smith   PetscFunctionReturn(0);
5195c6c1daeSBarry Smith }
5205c6c1daeSBarry Smith 
5218c557b5aSVaclav Hapla /*@
522874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
523874270d9SVaclav Hapla   and return this group's ID and file ID.
524874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
525874270d9SVaclav Hapla 
526874270d9SVaclav Hapla   Not collective
527874270d9SVaclav Hapla 
528874270d9SVaclav Hapla   Input Parameter:
529874270d9SVaclav Hapla . viewer - the PetscViewer
530874270d9SVaclav Hapla 
531874270d9SVaclav Hapla   Output Parameter:
532874270d9SVaclav Hapla + fileId - The HDF5 file ID
533874270d9SVaclav Hapla - groupId - The HDF5 group ID
534874270d9SVaclav Hapla 
535e74616beSVaclav Hapla   Notes:
536e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
537e74616beSVaclav Hapla 
538874270d9SVaclav Hapla   Level: intermediate
539874270d9SVaclav Hapla 
540874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
541874270d9SVaclav Hapla @*/
54254dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
54354dbf706SVaclav Hapla {
54490e03537SVaclav Hapla   hid_t          file_id;
54590e03537SVaclav Hapla   H5O_type_t     type;
54654dbf706SVaclav Hapla   const char     *groupName = NULL;
547e74616beSVaclav Hapla   PetscBool      create;
54854dbf706SVaclav Hapla   PetscErrorCode ierr;
54954dbf706SVaclav Hapla 
55054dbf706SVaclav Hapla   PetscFunctionBegin;
551e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
55254dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
55354dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
554e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
55590e03537SVaclav 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);
55690e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
55754dbf706SVaclav Hapla   *fileId  = file_id;
55854dbf706SVaclav Hapla   PetscFunctionReturn(0);
55954dbf706SVaclav Hapla }
56054dbf706SVaclav Hapla 
5613ef9c667SSatish Balay /*@
5625c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
5635c6c1daeSBarry Smith 
5645c6c1daeSBarry Smith   Not collective
5655c6c1daeSBarry Smith 
5665c6c1daeSBarry Smith   Input Parameter:
5675c6c1daeSBarry Smith . viewer - the PetscViewer
5685c6c1daeSBarry Smith 
5695c6c1daeSBarry Smith   Level: intermediate
5705c6c1daeSBarry Smith 
5715c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
5725c6c1daeSBarry Smith @*/
5735c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
5745c6c1daeSBarry Smith {
5755c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5765c6c1daeSBarry Smith 
5775c6c1daeSBarry Smith   PetscFunctionBegin;
5785c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5795c6c1daeSBarry Smith   ++hdf5->timestep;
5805c6c1daeSBarry Smith   PetscFunctionReturn(0);
5815c6c1daeSBarry Smith }
5825c6c1daeSBarry Smith 
5833ef9c667SSatish Balay /*@
5845c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
5855c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
5865c6c1daeSBarry Smith 
5875c6c1daeSBarry Smith   Not collective
5885c6c1daeSBarry Smith 
5895c6c1daeSBarry Smith   Input Parameters:
5905c6c1daeSBarry Smith + viewer - the PetscViewer
5915c6c1daeSBarry Smith - timestep - The timestep number
5925c6c1daeSBarry Smith 
5935c6c1daeSBarry Smith   Level: intermediate
5945c6c1daeSBarry Smith 
5955c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
5965c6c1daeSBarry Smith @*/
5975c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
5985c6c1daeSBarry Smith {
5995c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6005c6c1daeSBarry Smith 
6015c6c1daeSBarry Smith   PetscFunctionBegin;
6025c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6035c6c1daeSBarry Smith   hdf5->timestep = timestep;
6045c6c1daeSBarry Smith   PetscFunctionReturn(0);
6055c6c1daeSBarry Smith }
6065c6c1daeSBarry Smith 
6073ef9c667SSatish Balay /*@
6085c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
6095c6c1daeSBarry Smith 
6105c6c1daeSBarry Smith   Not collective
6115c6c1daeSBarry Smith 
6125c6c1daeSBarry Smith   Input Parameter:
6135c6c1daeSBarry Smith . viewer - the PetscViewer
6145c6c1daeSBarry Smith 
6155c6c1daeSBarry Smith   Output Parameter:
6165c6c1daeSBarry Smith . timestep - The timestep number
6175c6c1daeSBarry Smith 
6185c6c1daeSBarry Smith   Level: intermediate
6195c6c1daeSBarry Smith 
6205c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
6215c6c1daeSBarry Smith @*/
6225c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
6235c6c1daeSBarry Smith {
6245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6255c6c1daeSBarry Smith 
6265c6c1daeSBarry Smith   PetscFunctionBegin;
6275c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6285c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
6295c6c1daeSBarry Smith   *timestep = hdf5->timestep;
6305c6c1daeSBarry Smith   PetscFunctionReturn(0);
6315c6c1daeSBarry Smith }
6325c6c1daeSBarry Smith 
63336bce6e8SMatthew G. Knepley /*@C
63436bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
63536bce6e8SMatthew G. Knepley 
63636bce6e8SMatthew G. Knepley   Not collective
63736bce6e8SMatthew G. Knepley 
63836bce6e8SMatthew G. Knepley   Input Parameter:
63936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
64036bce6e8SMatthew G. Knepley 
64136bce6e8SMatthew G. Knepley   Output Parameter:
64236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
64336bce6e8SMatthew G. Knepley 
64436bce6e8SMatthew G. Knepley   Level: advanced
64536bce6e8SMatthew G. Knepley 
64636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
64736bce6e8SMatthew G. Knepley @*/
64836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
64936bce6e8SMatthew G. Knepley {
65036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
65136bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
65236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
65336bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
65436bce6e8SMatthew G. Knepley #else
65536bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
65636bce6e8SMatthew G. Knepley #endif
65736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
65836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
65936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
660c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
661de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
66236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
66336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
66436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
6657e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
66636bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
66736bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
66836bce6e8SMatthew G. Knepley }
66936bce6e8SMatthew G. Knepley 
67036bce6e8SMatthew G. Knepley /*@C
67136bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
67236bce6e8SMatthew G. Knepley 
67336bce6e8SMatthew G. Knepley   Not collective
67436bce6e8SMatthew G. Knepley 
67536bce6e8SMatthew G. Knepley   Input Parameter:
67636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
67736bce6e8SMatthew G. Knepley 
67836bce6e8SMatthew G. Knepley   Output Parameter:
67936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
68036bce6e8SMatthew G. Knepley 
68136bce6e8SMatthew G. Knepley   Level: advanced
68236bce6e8SMatthew G. Knepley 
68336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
68436bce6e8SMatthew G. Knepley @*/
68536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
68636bce6e8SMatthew G. Knepley {
68736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
68836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
68936bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
69036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
69136bce6e8SMatthew G. Knepley #else
69236bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
69336bce6e8SMatthew G. Knepley #endif
69436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
69536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
69636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
69736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
69836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
69936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
7007e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
70136bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
70236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
70336bce6e8SMatthew G. Knepley }
70436bce6e8SMatthew G. Knepley 
705df863907SAlex Fikl /*@C
706b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
70736bce6e8SMatthew G. Knepley 
70836bce6e8SMatthew G. Knepley   Input Parameters:
70936bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
71057d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
71136bce6e8SMatthew G. Knepley . name   - The attribute name
71236bce6e8SMatthew G. Knepley . datatype - The attribute type
71336bce6e8SMatthew G. Knepley - value    - The attribute value
71436bce6e8SMatthew G. Knepley 
71536bce6e8SMatthew G. Knepley   Level: advanced
71636bce6e8SMatthew G. Knepley 
717577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
71836bce6e8SMatthew G. Knepley @*/
71957d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
72036bce6e8SMatthew G. Knepley {
72157d22f7dSVaclav Hapla   char           *parent;
72260568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
723080f144cSVaclav Hapla   PetscBool      has;
72436bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
72536bce6e8SMatthew G. Knepley 
72636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
7275cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
72857d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
729c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
730b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
73157d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
732bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
733b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
73436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
7357e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
7367e97c476SMatthew G. Knepley     size_t len;
7377e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
738729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
7397e97c476SMatthew G. Knepley   }
74036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
741729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
74260568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
743080f144cSVaclav Hapla   if (has) {
744080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
745080f144cSVaclav Hapla   } else {
74660568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
747080f144cSVaclav Hapla   }
748729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
749729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
750729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
75160568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
752729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
75357d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
75436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
75536bce6e8SMatthew G. Knepley }
75636bce6e8SMatthew G. Knepley 
757df863907SAlex Fikl /*@C
758577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
759577e0e04SVaclav Hapla 
760577e0e04SVaclav Hapla   Input Parameters:
761577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
762577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
763577e0e04SVaclav Hapla . name     - The attribute name
764577e0e04SVaclav Hapla . datatype - The attribute type
765577e0e04SVaclav Hapla - value    - The attribute value
766577e0e04SVaclav Hapla 
767577e0e04SVaclav Hapla   Notes:
768577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
769577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
770577e0e04SVaclav Hapla 
771577e0e04SVaclav Hapla   Level: advanced
772577e0e04SVaclav Hapla 
773577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
774577e0e04SVaclav Hapla @*/
775577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
776577e0e04SVaclav Hapla {
777577e0e04SVaclav Hapla   PetscErrorCode ierr;
778577e0e04SVaclav Hapla 
779577e0e04SVaclav Hapla   PetscFunctionBegin;
780577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
781577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
782577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
783b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
784577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
785577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
786577e0e04SVaclav Hapla   PetscFunctionReturn(0);
787577e0e04SVaclav Hapla }
788577e0e04SVaclav Hapla 
789577e0e04SVaclav Hapla /*@C
790b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
791ad0c61feSMatthew G. Knepley 
792ad0c61feSMatthew G. Knepley   Input Parameters:
793ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
79457d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
795ad0c61feSMatthew G. Knepley . name   - The attribute name
796ad0c61feSMatthew G. Knepley - datatype - The attribute type
797ad0c61feSMatthew G. Knepley 
798ad0c61feSMatthew G. Knepley   Output Parameter:
799ad0c61feSMatthew G. Knepley . value    - The attribute value
800ad0c61feSMatthew G. Knepley 
801ad0c61feSMatthew G. Knepley   Level: advanced
802ad0c61feSMatthew G. Knepley 
803577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
804ad0c61feSMatthew G. Knepley @*/
80557d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
806ad0c61feSMatthew G. Knepley {
80757d22f7dSVaclav Hapla   char           *parent;
808f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
809969748fdSVaclav Hapla   PetscBool      has;
810ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
811ad0c61feSMatthew G. Knepley 
812ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
8135cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
81457d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
815c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
816b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
81757d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
818969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
819969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
820969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
821ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
822ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
82360568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
82460568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
825f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
826f0b58479SMatthew G. Knepley     size_t len;
827e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
82815b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
829f0b58479SMatthew G. Knepley     PetscStackCallHDF5(H5Tclose,(atype));
830f0b58479SMatthew G. Knepley     ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr);
831f0b58479SMatthew G. Knepley   }
83270efba33SBarry Smith   PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
833729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
834e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
835e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
83657d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
837ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
838ad0c61feSMatthew G. Knepley }
839ad0c61feSMatthew G. Knepley 
840577e0e04SVaclav Hapla /*@C
841577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
842577e0e04SVaclav Hapla 
843577e0e04SVaclav Hapla   Input Parameters:
844577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
845577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
846577e0e04SVaclav Hapla . name     - The attribute name
847577e0e04SVaclav Hapla - datatype - The attribute type
848577e0e04SVaclav Hapla 
849577e0e04SVaclav Hapla   Output Parameter:
850577e0e04SVaclav Hapla . value    - The attribute value
851577e0e04SVaclav Hapla 
852577e0e04SVaclav Hapla   Notes:
853577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
854577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
855577e0e04SVaclav Hapla 
856577e0e04SVaclav Hapla   Level: advanced
857577e0e04SVaclav Hapla 
858577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
859577e0e04SVaclav Hapla @*/
860577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
861577e0e04SVaclav Hapla {
862577e0e04SVaclav Hapla   PetscErrorCode ierr;
863577e0e04SVaclav Hapla 
864577e0e04SVaclav Hapla   PetscFunctionBegin;
865577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
866577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
867577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
868b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
869577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
870577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
871577e0e04SVaclav Hapla   PetscFunctionReturn(0);
872577e0e04SVaclav Hapla }
873577e0e04SVaclav Hapla 
874a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
875a75c3fd4SVaclav Hapla {
876a75c3fd4SVaclav Hapla   htri_t exists;
877a75c3fd4SVaclav Hapla   hid_t group;
878a75c3fd4SVaclav Hapla 
879a75c3fd4SVaclav Hapla   PetscFunctionBegin;
880a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
881a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
882a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
883a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
884a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
885a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
886a75c3fd4SVaclav Hapla   }
887a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
888a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
889a75c3fd4SVaclav Hapla }
890a75c3fd4SVaclav Hapla 
891bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
8925cdeb108SMatthew G. Knepley {
89390e03537SVaclav Hapla   const char     rootGroupName[] = "/";
8945cdeb108SMatthew G. Knepley   hid_t          h5;
895e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
89615b861d2SVaclav Hapla   PetscInt       i;
89715b861d2SVaclav Hapla   int            n;
89885688ae2SVaclav Hapla   char           **hierarchy;
89985688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
9005cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
9015cdeb108SMatthew G. Knepley 
9025cdeb108SMatthew G. Knepley   PetscFunctionBegin;
9035cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
90490e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
90590e03537SVaclav Hapla   else name = rootGroupName;
906ccd66a83SVaclav Hapla   if (has) {
90756cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
9085cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
909ccd66a83SVaclav Hapla   }
910ccd66a83SVaclav Hapla   if (otype) {
911ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
91256cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
913ccd66a83SVaclav Hapla   }
914c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
91585688ae2SVaclav Hapla 
91685688ae2SVaclav Hapla   /*
91785688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
91885688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
91985688ae2SVaclav Hapla      1) whether it's a valid link
92085688ae2SVaclav Hapla      2) whether this link resolves to an object
92185688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
92285688ae2SVaclav Hapla   */
92385688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
92485688ae2SVaclav Hapla   if (!n) {
92585688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
926ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
927ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
92885688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
92985688ae2SVaclav Hapla     PetscFunctionReturn(0);
93085688ae2SVaclav Hapla   }
93185688ae2SVaclav Hapla   for (i=0; i<n; i++) {
93285688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
93385688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
934a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
935a75c3fd4SVaclav Hapla     if (!exists) break;
93685688ae2SVaclav Hapla   }
93785688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
93885688ae2SVaclav Hapla 
93985688ae2SVaclav Hapla   /* If the object exists, get its type */
940ccd66a83SVaclav Hapla   if (exists && otype) {
9415cdeb108SMatthew G. Knepley     H5O_info_t info;
9425cdeb108SMatthew G. Knepley 
94376276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
94404633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
94556cc0592SVaclav Hapla     *otype = info.type;
9465cdeb108SMatthew G. Knepley   }
947ccd66a83SVaclav Hapla   if (has) *has = exists;
9485cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
9495cdeb108SMatthew G. Knepley }
9505cdeb108SMatthew G. Knepley 
951ecce1506SVaclav Hapla /*@
9520a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
9530a9f38efSVaclav Hapla 
9540a9f38efSVaclav Hapla   Input Parameters:
9550a9f38efSVaclav Hapla . viewer - The HDF5 viewer
9560a9f38efSVaclav Hapla 
9570a9f38efSVaclav Hapla   Output Parameter:
9580a9f38efSVaclav Hapla . has    - Flag for group existence
9590a9f38efSVaclav Hapla 
9600a9f38efSVaclav Hapla   Notes:
9610a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
9620a9f38efSVaclav Hapla 
9630a9f38efSVaclav Hapla   Level: advanced
9640a9f38efSVaclav Hapla 
9650a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
9660a9f38efSVaclav Hapla @*/
9670a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
9680a9f38efSVaclav Hapla {
9690a9f38efSVaclav Hapla   H5O_type_t type;
9700a9f38efSVaclav Hapla   const char *name;
9710a9f38efSVaclav Hapla   PetscErrorCode ierr;
9720a9f38efSVaclav Hapla 
9730a9f38efSVaclav Hapla   PetscFunctionBegin;
9740a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
9750a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
9760a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
9770a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
9780a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
9790a9f38efSVaclav Hapla   PetscFunctionReturn(0);
9800a9f38efSVaclav Hapla }
9810a9f38efSVaclav Hapla 
9820a9f38efSVaclav Hapla /*@
983e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
984ecce1506SVaclav Hapla 
985ecce1506SVaclav Hapla   Input Parameters:
986ecce1506SVaclav Hapla + viewer - The HDF5 viewer
987ecce1506SVaclav Hapla - obj    - The named object
988ecce1506SVaclav Hapla 
989ecce1506SVaclav Hapla   Output Parameter:
990ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
991ecce1506SVaclav Hapla 
992e3f143f7SVaclav Hapla   Notes:
993e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
994e3f143f7SVaclav Hapla 
995ecce1506SVaclav Hapla   Level: advanced
996ecce1506SVaclav Hapla 
997e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
998ecce1506SVaclav Hapla @*/
999ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1000ecce1506SVaclav Hapla {
100156cc0592SVaclav Hapla   H5O_type_t type;
10026c132bc1SVaclav Hapla   char *path;
1003ecce1506SVaclav Hapla   PetscErrorCode ierr;
1004ecce1506SVaclav Hapla 
1005ecce1506SVaclav Hapla   PetscFunctionBegin;
1006c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1007c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1008c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1009ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1010e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
10116c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
10126c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
101356cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
10146c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1015ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1016ecce1506SVaclav Hapla }
1017ecce1506SVaclav Hapla 
1018df863907SAlex Fikl /*@C
1019ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1020e7bdbf8eSMatthew G. Knepley 
1021e7bdbf8eSMatthew G. Knepley   Input Parameters:
1022e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
102310e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1024e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1025e7bdbf8eSMatthew G. Knepley 
1026e7bdbf8eSMatthew G. Knepley   Output Parameter:
1027e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1028e7bdbf8eSMatthew G. Knepley 
1029e7bdbf8eSMatthew G. Knepley   Level: advanced
1030e7bdbf8eSMatthew G. Knepley 
1031577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1032e7bdbf8eSMatthew G. Knepley @*/
103310e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1034e7bdbf8eSMatthew G. Knepley {
103510e69e8fSVaclav Hapla   char           *parent;
1036e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1037e7bdbf8eSMatthew G. Knepley 
1038e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
10395cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
104010e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1041c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1042c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
104310e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1044bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
104510e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
104610e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
104706db490cSVaclav Hapla   PetscFunctionReturn(0);
104806db490cSVaclav Hapla }
104906db490cSVaclav Hapla 
1050577e0e04SVaclav Hapla /*@C
1051577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1052577e0e04SVaclav Hapla 
1053577e0e04SVaclav Hapla   Input Parameters:
1054577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1055577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1056577e0e04SVaclav Hapla - name   - The attribute name
1057577e0e04SVaclav Hapla 
1058577e0e04SVaclav Hapla   Output Parameter:
1059577e0e04SVaclav Hapla . has    - Flag for attribute existence
1060577e0e04SVaclav Hapla 
1061577e0e04SVaclav Hapla   Notes:
1062577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1063577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1064577e0e04SVaclav Hapla 
1065577e0e04SVaclav Hapla   Level: advanced
1066577e0e04SVaclav Hapla 
1067577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1068577e0e04SVaclav Hapla @*/
1069577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1070577e0e04SVaclav Hapla {
1071577e0e04SVaclav Hapla   PetscErrorCode ierr;
1072577e0e04SVaclav Hapla 
1073577e0e04SVaclav Hapla   PetscFunctionBegin;
1074577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1075577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1076577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1077577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1078577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1079577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1080577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1081577e0e04SVaclav Hapla }
1082577e0e04SVaclav Hapla 
108306db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
108406db490cSVaclav Hapla {
108506db490cSVaclav Hapla   hid_t          h5;
108606db490cSVaclav Hapla   htri_t         hhas;
108706db490cSVaclav Hapla   PetscErrorCode ierr;
108806db490cSVaclav Hapla 
108906db490cSVaclav Hapla   PetscFunctionBegin;
109006db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
10912f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
10925cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1093e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1094e7bdbf8eSMatthew G. Knepley }
1095e7bdbf8eSMatthew G. Knepley 
1096a75e6a4aSMatthew G. Knepley /*
1097a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1098a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1099a75e6a4aSMatthew G. Knepley */
1100d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1101a75e6a4aSMatthew G. Knepley 
1102a75e6a4aSMatthew G. Knepley /*@C
1103a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1104a75e6a4aSMatthew G. Knepley 
1105a75e6a4aSMatthew G. Knepley   Collective on MPI_Comm
1106a75e6a4aSMatthew G. Knepley 
1107a75e6a4aSMatthew G. Knepley   Input Parameter:
1108a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1109a75e6a4aSMatthew G. Knepley 
1110a75e6a4aSMatthew G. Knepley   Level: intermediate
1111a75e6a4aSMatthew G. Knepley 
1112a75e6a4aSMatthew G. Knepley   Options Database Keys:
1113a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1114a75e6a4aSMatthew G. Knepley 
1115a75e6a4aSMatthew G. Knepley   Environmental variables:
1116a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1117a75e6a4aSMatthew G. Knepley 
1118a75e6a4aSMatthew G. Knepley   Notes:
1119a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1120a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1121a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1122a75e6a4aSMatthew G. Knepley 
1123a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1124a75e6a4aSMatthew G. Knepley @*/
1125a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1126a75e6a4aSMatthew G. Knepley {
1127a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1128a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1129a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1130a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1131a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1132a75e6a4aSMatthew G. Knepley 
1133a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1134a75e6a4aSMatthew 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);}
1135a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
113612801b39SBarry Smith     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0);
1137a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1138a75e6a4aSMatthew G. Knepley   }
113947435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1140a75e6a4aSMatthew G. Knepley   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1141a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1142a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
1143a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1144a75e6a4aSMatthew G. Knepley     if (!flg) {
1145a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
1146a75e6a4aSMatthew G. Knepley       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1147a75e6a4aSMatthew G. Knepley     }
1148a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
1149a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
1150a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
1151a75e6a4aSMatthew G. Knepley     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);}
115247435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
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   }
1155a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
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   PetscFunctionReturn(viewer);
1158a75e6a4aSMatthew G. Knepley }
1159