xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 7e4fd57318ebf5ced89996365070bb2d8eebba11)
1bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h>
27d964842SVaclav Hapla #include <petsc/private/viewerhdf5impl.h>
3d70abbfaSBarry Smith #include <petscviewerhdf5.h>    /*I   "petscviewerhdf5.h"   I*/
45c6c1daeSBarry Smith 
5bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
706db490cSVaclav Hapla 
86c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath)
96c132bc1SVaclav Hapla {
106c132bc1SVaclav Hapla   const char *group;
116c132bc1SVaclav Hapla   char buf[PETSC_MAX_PATH_LEN]="";
126c132bc1SVaclav Hapla   PetscErrorCode ierr;
136c132bc1SVaclav Hapla 
146c132bc1SVaclav Hapla   PetscFunctionBegin;
156c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
166c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, group);CHKERRQ(ierr);
176c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
186c132bc1SVaclav Hapla   ierr = PetscStrcat(buf, objname);CHKERRQ(ierr);
196c132bc1SVaclav Hapla   ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr);
206c132bc1SVaclav Hapla   PetscFunctionReturn(0);
216c132bc1SVaclav Hapla }
226c132bc1SVaclav Hapla 
23577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
24577e0e04SVaclav Hapla {
25577e0e04SVaclav Hapla   PetscBool has;
26577e0e04SVaclav Hapla   const char *group;
27577e0e04SVaclav Hapla   PetscErrorCode ierr;
28577e0e04SVaclav Hapla 
29577e0e04SVaclav Hapla   PetscFunctionBegin;
30577e0e04SVaclav Hapla   if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
31577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
32577e0e04SVaclav Hapla   if (!has) {
33577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
34577e0e04SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group);
35577e0e04SVaclav Hapla   }
36577e0e04SVaclav Hapla   PetscFunctionReturn(0);
37577e0e04SVaclav Hapla }
38577e0e04SVaclav Hapla 
394416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
4082ea9b62SBarry Smith {
4182ea9b62SBarry Smith   PetscErrorCode   ierr;
42ee8b9732SVaclav Hapla   PetscBool        flg = PETSC_FALSE, set;
4382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
4482ea9b62SBarry Smith 
4582ea9b62SBarry Smith   PetscFunctionBegin;
4682ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
4782ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
489a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
49ee8b9732SVaclav Hapla   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
50ee8b9732SVaclav Hapla   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
5182ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5282ea9b62SBarry Smith   PetscFunctionReturn(0);
5382ea9b62SBarry Smith }
5482ea9b62SBarry Smith 
551b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
561b793a25SVaclav Hapla {
571b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
5803fa8834SVaclav Hapla   PetscBool         flg;
591b793a25SVaclav Hapla   PetscErrorCode    ierr;
601b793a25SVaclav Hapla 
611b793a25SVaclav Hapla   PetscFunctionBegin;
621b793a25SVaclav Hapla   if (hdf5->filename) {
631b793a25SVaclav Hapla     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
641b793a25SVaclav Hapla   }
651b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
661b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
6703fa8834SVaclav Hapla   ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr);
6803fa8834SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr);
691b793a25SVaclav Hapla   PetscFunctionReturn(0);
701b793a25SVaclav Hapla }
711b793a25SVaclav Hapla 
725c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
735c6c1daeSBarry Smith {
745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
755c6c1daeSBarry Smith   PetscErrorCode   ierr;
765c6c1daeSBarry Smith 
775c6c1daeSBarry Smith   PetscFunctionBegin;
785c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
79729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
805c6c1daeSBarry Smith   PetscFunctionReturn(0);
815c6c1daeSBarry Smith }
825c6c1daeSBarry Smith 
839b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
845c6c1daeSBarry Smith {
855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
865c6c1daeSBarry Smith   PetscErrorCode   ierr;
875c6c1daeSBarry Smith 
885c6c1daeSBarry Smith   PetscFunctionBegin;
899c5072fbSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
905c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
915c6c1daeSBarry Smith   while (hdf5->groups) {
927d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
935c6c1daeSBarry Smith 
945c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
955c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
965c6c1daeSBarry Smith     hdf5->groups = tmp;
975c6c1daeSBarry Smith   }
985c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
990b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
100d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1010b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
102058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
103058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1045c6c1daeSBarry Smith   PetscFunctionReturn(0);
1055c6c1daeSBarry Smith }
1065c6c1daeSBarry Smith 
1079b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1085c6c1daeSBarry Smith {
1095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1105c6c1daeSBarry Smith 
1115c6c1daeSBarry Smith   PetscFunctionBegin;
1125c6c1daeSBarry Smith   hdf5->btype = type;
1135c6c1daeSBarry Smith   PetscFunctionReturn(0);
1145c6c1daeSBarry Smith }
1155c6c1daeSBarry Smith 
1160b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1170b62783dSVaclav Hapla {
1180b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1190b62783dSVaclav Hapla 
1200b62783dSVaclav Hapla   PetscFunctionBegin;
1210b62783dSVaclav Hapla   *type = hdf5->btype;
1220b62783dSVaclav Hapla   PetscFunctionReturn(0);
1230b62783dSVaclav Hapla }
1240b62783dSVaclav Hapla 
1259b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
12682ea9b62SBarry Smith {
12782ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
12882ea9b62SBarry Smith 
12982ea9b62SBarry Smith   PetscFunctionBegin;
13082ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13182ea9b62SBarry Smith   PetscFunctionReturn(0);
13282ea9b62SBarry Smith }
13382ea9b62SBarry Smith 
134df863907SAlex Fikl /*@
13582ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13682ea9b62SBarry Smith        dimension of 2.
13782ea9b62SBarry Smith 
13882ea9b62SBarry Smith     Logically Collective on PetscViewer
13982ea9b62SBarry Smith 
14082ea9b62SBarry Smith   Input Parameters:
14182ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
14282ea9b62SBarry 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
14382ea9b62SBarry Smith 
14482ea9b62SBarry Smith   Options Database:
14582ea9b62SBarry 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
14682ea9b62SBarry Smith 
14782ea9b62SBarry Smith 
14895452b02SPatrick Sanan   Notes:
14995452b02SPatrick 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
15082ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith   Level: intermediate
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
15582ea9b62SBarry Smith 
15682ea9b62SBarry Smith @*/
15782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
15882ea9b62SBarry Smith {
15982ea9b62SBarry Smith   PetscErrorCode ierr;
16082ea9b62SBarry Smith 
16182ea9b62SBarry Smith   PetscFunctionBegin;
16282ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
16382ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
16482ea9b62SBarry Smith   PetscFunctionReturn(0);
16582ea9b62SBarry Smith }
16682ea9b62SBarry Smith 
167df863907SAlex Fikl /*@
16882ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
16982ea9b62SBarry Smith        dimension of 2.
17082ea9b62SBarry Smith 
17182ea9b62SBarry Smith     Logically Collective on PetscViewer
17282ea9b62SBarry Smith 
17382ea9b62SBarry Smith   Input Parameter:
17482ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
17582ea9b62SBarry Smith 
17682ea9b62SBarry Smith   Output Parameter:
17782ea9b62SBarry 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
17882ea9b62SBarry Smith 
17995452b02SPatrick Sanan   Notes:
18095452b02SPatrick 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
18182ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
18282ea9b62SBarry Smith 
18382ea9b62SBarry Smith   Level: intermediate
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
18682ea9b62SBarry Smith 
18782ea9b62SBarry Smith @*/
18882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
18982ea9b62SBarry Smith {
19082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
19182ea9b62SBarry Smith 
19282ea9b62SBarry Smith   PetscFunctionBegin;
19382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
19482ea9b62SBarry Smith   *flg = hdf5->basedimension2;
19582ea9b62SBarry Smith   PetscFunctionReturn(0);
19682ea9b62SBarry Smith }
19782ea9b62SBarry Smith 
1989b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
1999a0502c6SHåkon Strandenes {
2009a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2019a0502c6SHåkon Strandenes 
2029a0502c6SHåkon Strandenes   PetscFunctionBegin;
2039a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2049a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2059a0502c6SHåkon Strandenes }
2069a0502c6SHåkon Strandenes 
207df863907SAlex Fikl /*@
2089a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2099a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2109a0502c6SHåkon Strandenes 
2119a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2129a0502c6SHåkon Strandenes 
2139a0502c6SHåkon Strandenes   Input Parameters:
2149a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2159a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2169a0502c6SHåkon Strandenes 
2179a0502c6SHåkon Strandenes   Options Database:
2189a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2199a0502c6SHåkon Strandenes 
2209a0502c6SHåkon Strandenes 
22195452b02SPatrick Sanan   Notes:
22295452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2239a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2249a0502c6SHåkon Strandenes 
2259a0502c6SHåkon Strandenes   Level: intermediate
2269a0502c6SHåkon Strandenes 
2279a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2289a0502c6SHåkon Strandenes           PetscReal
2299a0502c6SHåkon Strandenes 
2309a0502c6SHåkon Strandenes @*/
2319a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2329a0502c6SHåkon Strandenes {
2339a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2349a0502c6SHåkon Strandenes 
2359a0502c6SHåkon Strandenes   PetscFunctionBegin;
2369a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2379a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2389a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2399a0502c6SHåkon Strandenes }
2409a0502c6SHåkon Strandenes 
241df863907SAlex Fikl /*@
2429a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2439a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2449a0502c6SHåkon Strandenes 
2459a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2469a0502c6SHåkon Strandenes 
2479a0502c6SHåkon Strandenes   Input Parameter:
2489a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2499a0502c6SHåkon Strandenes 
2509a0502c6SHåkon Strandenes   Output Parameter:
2519a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2529a0502c6SHåkon Strandenes 
25395452b02SPatrick Sanan   Notes:
25495452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2559a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2569a0502c6SHåkon Strandenes 
2579a0502c6SHåkon Strandenes   Level: intermediate
2589a0502c6SHåkon Strandenes 
2599a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2609a0502c6SHåkon Strandenes           PetscReal
2619a0502c6SHåkon Strandenes 
2629a0502c6SHåkon Strandenes @*/
2639a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2649a0502c6SHåkon Strandenes {
2659a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2669a0502c6SHåkon Strandenes 
2679a0502c6SHåkon Strandenes   PetscFunctionBegin;
2689a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2699a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2709a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2719a0502c6SHåkon Strandenes }
2729a0502c6SHåkon Strandenes 
273ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
274ee8b9732SVaclav Hapla {
275ee8b9732SVaclav Hapla   PetscFunctionBegin;
276ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
277ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
278c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
279ee8b9732SVaclav Hapla   {
280ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
281ee8b9732SVaclav Hapla     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
282ee8b9732SVaclav Hapla   }
283ee8b9732SVaclav Hapla #else
284ee8b9732SVaclav Hapla   if (flg) {
285ee8b9732SVaclav Hapla     PetscErrorCode ierr;
286c796909eSBarry Smith     ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n");CHKERRQ(ierr);
287ee8b9732SVaclav Hapla   }
288ee8b9732SVaclav Hapla #endif
289ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
290ee8b9732SVaclav Hapla }
291ee8b9732SVaclav Hapla 
292ee8b9732SVaclav Hapla /*@
293ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
294ee8b9732SVaclav Hapla 
295ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
296ee8b9732SVaclav Hapla 
297ee8b9732SVaclav Hapla   Input Parameters:
298ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
299ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
300ee8b9732SVaclav Hapla 
301ee8b9732SVaclav Hapla   Options Database:
302ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
303ee8b9732SVaclav Hapla 
304ee8b9732SVaclav Hapla   Notes:
305ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
30653effed7SBarry Smith   However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions.
307ee8b9732SVaclav Hapla 
308ee8b9732SVaclav Hapla   Developer notes:
309ee8b9732SVaclav Hapla   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
310ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
311ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
312ee8b9732SVaclav Hapla 
313ee8b9732SVaclav Hapla   Level: intermediate
314ee8b9732SVaclav Hapla 
315ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
316ee8b9732SVaclav Hapla 
317ee8b9732SVaclav Hapla @*/
318ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
319ee8b9732SVaclav Hapla {
320ee8b9732SVaclav Hapla   PetscErrorCode ierr;
321ee8b9732SVaclav Hapla 
322ee8b9732SVaclav Hapla   PetscFunctionBegin;
323ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
324ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer,flg,2);
325ee8b9732SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
326ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
327ee8b9732SVaclav Hapla }
328ee8b9732SVaclav Hapla 
329ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
330ee8b9732SVaclav Hapla {
331c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
332ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
333ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
334c796909eSBarry Smith #endif
335ee8b9732SVaclav Hapla 
336ee8b9732SVaclav Hapla   PetscFunctionBegin;
337c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
338c796909eSBarry Smith   *flg = PETSC_FALSE;
339c796909eSBarry Smith #else
340ee8b9732SVaclav Hapla   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
341ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
342c796909eSBarry Smith #endif
343ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
344ee8b9732SVaclav Hapla }
345ee8b9732SVaclav Hapla 
346ee8b9732SVaclav Hapla /*@
347ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
348ee8b9732SVaclav Hapla 
349ee8b9732SVaclav Hapla   Not Collective
350ee8b9732SVaclav Hapla 
351ee8b9732SVaclav Hapla   Input Parameters:
352ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
353ee8b9732SVaclav Hapla 
354ee8b9732SVaclav Hapla   Output Parameters:
355ee8b9732SVaclav Hapla . flg - the flag
356ee8b9732SVaclav Hapla 
357ee8b9732SVaclav Hapla   Level: intermediate
358ee8b9732SVaclav Hapla 
359ee8b9732SVaclav Hapla   Notes:
360c796909eSBarry Smith   This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned.
361ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
362ee8b9732SVaclav Hapla 
363ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
364ee8b9732SVaclav Hapla 
365ee8b9732SVaclav Hapla @*/
366ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
367ee8b9732SVaclav Hapla {
368ee8b9732SVaclav Hapla   PetscErrorCode ierr;
369ee8b9732SVaclav Hapla 
370ee8b9732SVaclav Hapla   PetscFunctionBegin;
371ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
372534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
373ee8b9732SVaclav Hapla 
374ee8b9732SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
375ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
376ee8b9732SVaclav Hapla }
377ee8b9732SVaclav Hapla 
3789b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3795c6c1daeSBarry Smith {
3805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3815c6c1daeSBarry Smith   hid_t             plist_id;
3825c6c1daeSBarry Smith   PetscErrorCode    ierr;
3835c6c1daeSBarry Smith 
3845c6c1daeSBarry Smith   PetscFunctionBegin;
385f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
386f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
3875c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3885c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
389729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
390c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
391d28bfa55SVaclav Hapla   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
392c796909eSBarry Smith #endif
3935c6c1daeSBarry Smith   /* Create or open the file collectively */
3945c6c1daeSBarry Smith   switch (hdf5->btype) {
3955c6c1daeSBarry Smith   case FILE_MODE_READ:
396729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
3975c6c1daeSBarry Smith     break;
3985c6c1daeSBarry Smith   case FILE_MODE_APPEND:
399*7e4fd573SVaclav Hapla   case FILE_MODE_UPDATE:
400*7e4fd573SVaclav Hapla   {
401*7e4fd573SVaclav Hapla     PetscBool flg;
402*7e4fd573SVaclav Hapla     ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr);
403*7e4fd573SVaclav Hapla     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
404*7e4fd573SVaclav Hapla     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4055c6c1daeSBarry Smith     break;
406*7e4fd573SVaclav Hapla   }
4075c6c1daeSBarry Smith   case FILE_MODE_WRITE:
408729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
4095c6c1daeSBarry Smith     break;
410*7e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED:
411*7e4fd573SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
4125c6c1daeSBarry Smith   default:
413*7e4fd573SVaclav Hapla     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
4145c6c1daeSBarry Smith   }
4155c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
416729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
4175c6c1daeSBarry Smith   PetscFunctionReturn(0);
4185c6c1daeSBarry Smith }
4195c6c1daeSBarry Smith 
420d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
421d1232d7fSToby Isaac {
422d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
423d1232d7fSToby Isaac 
424d1232d7fSToby Isaac   PetscFunctionBegin;
425d1232d7fSToby Isaac   *name = vhdf5->filename;
426d1232d7fSToby Isaac   PetscFunctionReturn(0);
427d1232d7fSToby Isaac }
428d1232d7fSToby Isaac 
429b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
430b723ab35SVaclav Hapla {
4315dc64a97SVaclav Hapla   /*
432b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
433b723ab35SVaclav Hapla   PetscErrorCode   ierr;
4345dc64a97SVaclav Hapla   */
435b723ab35SVaclav Hapla 
436b723ab35SVaclav Hapla   PetscFunctionBegin;
437b723ab35SVaclav Hapla   PetscFunctionReturn(0);
438b723ab35SVaclav Hapla }
439b723ab35SVaclav Hapla 
4408556b5ebSBarry Smith /*MC
4418556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4428556b5ebSBarry Smith 
4438556b5ebSBarry Smith 
4448556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4458556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4468556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4478556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4488556b5ebSBarry Smith 
4491b266c99SBarry Smith   Level: beginner
4508556b5ebSBarry Smith M*/
451d1232d7fSToby Isaac 
4528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4535c6c1daeSBarry Smith {
4545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4555c6c1daeSBarry Smith   PetscErrorCode   ierr;
4565c6c1daeSBarry Smith 
4575c6c1daeSBarry Smith   PetscFunctionBegin;
458b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4595c6c1daeSBarry Smith 
4605c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4615c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
46282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
463b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4641b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
465908793a3SLisandro Dalcin   v->ops->flush          = NULL;
466*7e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
467908793a3SLisandro Dalcin   hdf5->filename         = NULL;
4685c6c1daeSBarry Smith   hdf5->timestep         = -1;
4690298fd71SBarry Smith   hdf5->groups           = NULL;
4705c6c1daeSBarry Smith 
4719c5072fbSVaclav Hapla   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
4729c5072fbSVaclav Hapla 
4730b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
474d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4750b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
4760b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
47782ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4789a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
479ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
480ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
4815c6c1daeSBarry Smith   PetscFunctionReturn(0);
4825c6c1daeSBarry Smith }
4835c6c1daeSBarry Smith 
4845c6c1daeSBarry Smith /*@C
4855c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4865c6c1daeSBarry Smith 
487d083f849SBarry Smith    Collective
4885c6c1daeSBarry Smith 
4895c6c1daeSBarry Smith    Input Parameters:
4905c6c1daeSBarry Smith +  comm - MPI communicator
4915c6c1daeSBarry Smith .  name - name of file
4925c6c1daeSBarry Smith -  type - type of file
4935c6c1daeSBarry Smith 
4945c6c1daeSBarry Smith    Output Parameter:
4955c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
4965c6c1daeSBarry Smith 
49782ea9b62SBarry Smith   Options Database:
498a2b725a8SWilliam Gropp +  -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
499a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
50082ea9b62SBarry Smith 
5015c6c1daeSBarry Smith    Level: beginner
5025c6c1daeSBarry Smith 
503*7e4fd573SVaclav Hapla    Notes:
504*7e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
505*7e4fd573SVaclav Hapla +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
506*7e4fd573SVaclav Hapla .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
507*7e4fd573SVaclav Hapla .  FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
508*7e4fd573SVaclav Hapla -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
509*7e4fd573SVaclav Hapla 
510*7e4fd573SVaclav Hapla    In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified.
511*7e4fd573SVaclav Hapla 
5125c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5135c6c1daeSBarry Smith 
5145c6c1daeSBarry Smith 
5156a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5169a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
517a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5185c6c1daeSBarry Smith @*/
5195c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5205c6c1daeSBarry Smith {
5215c6c1daeSBarry Smith   PetscErrorCode ierr;
5225c6c1daeSBarry Smith 
5235c6c1daeSBarry Smith   PetscFunctionBegin;
5245c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5255c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5265c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5275c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
528b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
5295c6c1daeSBarry Smith   PetscFunctionReturn(0);
5305c6c1daeSBarry Smith }
5315c6c1daeSBarry Smith 
5325c6c1daeSBarry Smith /*@C
5335c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5345c6c1daeSBarry Smith 
5355c6c1daeSBarry Smith   Not collective
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith   Input Parameter:
5385c6c1daeSBarry Smith . viewer - the PetscViewer
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith   Output Parameter:
5415c6c1daeSBarry Smith . file_id - The file id
5425c6c1daeSBarry Smith 
5435c6c1daeSBarry Smith   Level: intermediate
5445c6c1daeSBarry Smith 
5455c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5465c6c1daeSBarry Smith @*/
5475c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5485c6c1daeSBarry Smith {
5495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5505c6c1daeSBarry Smith 
5515c6c1daeSBarry Smith   PetscFunctionBegin;
5525c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5535c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5545c6c1daeSBarry Smith   PetscFunctionReturn(0);
5555c6c1daeSBarry Smith }
5565c6c1daeSBarry Smith 
5575c6c1daeSBarry Smith /*@C
5585c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5595c6c1daeSBarry Smith 
5605c6c1daeSBarry Smith   Not collective
5615c6c1daeSBarry Smith 
5625c6c1daeSBarry Smith   Input Parameters:
5635c6c1daeSBarry Smith + viewer - the PetscViewer
5645c6c1daeSBarry Smith - name - The group name
5655c6c1daeSBarry Smith 
5665c6c1daeSBarry Smith   Level: intermediate
5675c6c1daeSBarry Smith 
568820fc2d1SVaclav Hapla   Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/".
569820fc2d1SVaclav Hapla 
570874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5715c6c1daeSBarry Smith @*/
572be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
5735c6c1daeSBarry Smith {
5745c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5757d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
5765c6c1daeSBarry Smith   PetscErrorCode   ierr;
5775c6c1daeSBarry Smith 
5785c6c1daeSBarry Smith   PetscFunctionBegin;
5795c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
580820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name,2);
581820fc2d1SVaclav Hapla   if (name && name[0]) {
582820fc2d1SVaclav Hapla      size_t i,len;
583820fc2d1SVaclav Hapla      ierr = PetscStrlen(name, &len);CHKERRQ(ierr);
584820fc2d1SVaclav Hapla      for (i=0; i<len; i++) if (name[i] != '/') break;
585820fc2d1SVaclav Hapla      if (i == len) name = NULL;
586820fc2d1SVaclav Hapla   } else name = NULL;
58795dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
5885c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr);
5895c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
5905c6c1daeSBarry Smith   hdf5->groups    = groupNode;
5915c6c1daeSBarry Smith   PetscFunctionReturn(0);
5925c6c1daeSBarry Smith }
5935c6c1daeSBarry Smith 
5943ef9c667SSatish Balay /*@
5955c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
5965c6c1daeSBarry Smith 
5975c6c1daeSBarry Smith   Not collective
5985c6c1daeSBarry Smith 
5995c6c1daeSBarry Smith   Input Parameter:
6005c6c1daeSBarry Smith . viewer - the PetscViewer
6015c6c1daeSBarry Smith 
6025c6c1daeSBarry Smith   Level: intermediate
6035c6c1daeSBarry Smith 
604874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
6055c6c1daeSBarry Smith @*/
6065c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
6075c6c1daeSBarry Smith {
6085c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6097d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6105c6c1daeSBarry Smith   PetscErrorCode   ierr;
6115c6c1daeSBarry Smith 
6125c6c1daeSBarry Smith   PetscFunctionBegin;
6135c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
61482f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6155c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6165c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6175c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6185c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6195c6c1daeSBarry Smith   PetscFunctionReturn(0);
6205c6c1daeSBarry Smith }
6215c6c1daeSBarry Smith 
6225c6c1daeSBarry Smith /*@C
623874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
624874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6255c6c1daeSBarry Smith 
6265c6c1daeSBarry Smith   Not collective
6275c6c1daeSBarry Smith 
6285c6c1daeSBarry Smith   Input Parameter:
6295c6c1daeSBarry Smith . viewer - the PetscViewer
6305c6c1daeSBarry Smith 
6315c6c1daeSBarry Smith   Output Parameter:
6325c6c1daeSBarry Smith . name - The group name
6335c6c1daeSBarry Smith 
6345c6c1daeSBarry Smith   Level: intermediate
6355c6c1daeSBarry Smith 
636874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6375c6c1daeSBarry Smith @*/
638be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
6395c6c1daeSBarry Smith {
6405c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6415c6c1daeSBarry Smith 
6425c6c1daeSBarry Smith   PetscFunctionBegin;
6435c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
644c959eef4SJed Brown   PetscValidPointer(name,2);
645a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6460298fd71SBarry Smith   else *name = NULL;
6475c6c1daeSBarry Smith   PetscFunctionReturn(0);
6485c6c1daeSBarry Smith }
6495c6c1daeSBarry Smith 
6508c557b5aSVaclav Hapla /*@
651874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
652874270d9SVaclav Hapla   and return this group's ID and file ID.
653874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
654874270d9SVaclav Hapla 
655874270d9SVaclav Hapla   Not collective
656874270d9SVaclav Hapla 
657874270d9SVaclav Hapla   Input Parameter:
658874270d9SVaclav Hapla . viewer - the PetscViewer
659874270d9SVaclav Hapla 
660874270d9SVaclav Hapla   Output Parameter:
661874270d9SVaclav Hapla + fileId - The HDF5 file ID
662874270d9SVaclav Hapla - groupId - The HDF5 group ID
663874270d9SVaclav Hapla 
664e74616beSVaclav Hapla   Notes:
665e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
666e74616beSVaclav Hapla 
667874270d9SVaclav Hapla   Level: intermediate
668874270d9SVaclav Hapla 
669874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
670874270d9SVaclav Hapla @*/
67154dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
67254dbf706SVaclav Hapla {
67390e03537SVaclav Hapla   hid_t          file_id;
67490e03537SVaclav Hapla   H5O_type_t     type;
67554dbf706SVaclav Hapla   const char     *groupName = NULL;
676e74616beSVaclav Hapla   PetscBool      create;
67754dbf706SVaclav Hapla   PetscErrorCode ierr;
67854dbf706SVaclav Hapla 
67954dbf706SVaclav Hapla   PetscFunctionBegin;
680e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
68154dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
68254dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
683e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
68490e03537SVaclav 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);
68590e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
68654dbf706SVaclav Hapla   *fileId  = file_id;
68754dbf706SVaclav Hapla   PetscFunctionReturn(0);
68854dbf706SVaclav Hapla }
68954dbf706SVaclav Hapla 
6903ef9c667SSatish Balay /*@
6915c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
6925c6c1daeSBarry Smith 
6935c6c1daeSBarry Smith   Not collective
6945c6c1daeSBarry Smith 
6955c6c1daeSBarry Smith   Input Parameter:
6965c6c1daeSBarry Smith . viewer - the PetscViewer
6975c6c1daeSBarry Smith 
6985c6c1daeSBarry Smith   Level: intermediate
6995c6c1daeSBarry Smith 
7005c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
7015c6c1daeSBarry Smith @*/
7025c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
7035c6c1daeSBarry Smith {
7045c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7055c6c1daeSBarry Smith 
7065c6c1daeSBarry Smith   PetscFunctionBegin;
7075c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7085c6c1daeSBarry Smith   ++hdf5->timestep;
7095c6c1daeSBarry Smith   PetscFunctionReturn(0);
7105c6c1daeSBarry Smith }
7115c6c1daeSBarry Smith 
7123ef9c667SSatish Balay /*@
7135c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
7145c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
7155c6c1daeSBarry Smith 
7165c6c1daeSBarry Smith   Not collective
7175c6c1daeSBarry Smith 
7185c6c1daeSBarry Smith   Input Parameters:
7195c6c1daeSBarry Smith + viewer - the PetscViewer
7205c6c1daeSBarry Smith - timestep - The timestep number
7215c6c1daeSBarry Smith 
7225c6c1daeSBarry Smith   Level: intermediate
7235c6c1daeSBarry Smith 
7245c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
7255c6c1daeSBarry Smith @*/
7265c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
7275c6c1daeSBarry Smith {
7285c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith   PetscFunctionBegin;
7315c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7325c6c1daeSBarry Smith   hdf5->timestep = timestep;
7335c6c1daeSBarry Smith   PetscFunctionReturn(0);
7345c6c1daeSBarry Smith }
7355c6c1daeSBarry Smith 
7363ef9c667SSatish Balay /*@
7375c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7385c6c1daeSBarry Smith 
7395c6c1daeSBarry Smith   Not collective
7405c6c1daeSBarry Smith 
7415c6c1daeSBarry Smith   Input Parameter:
7425c6c1daeSBarry Smith . viewer - the PetscViewer
7435c6c1daeSBarry Smith 
7445c6c1daeSBarry Smith   Output Parameter:
7455c6c1daeSBarry Smith . timestep - The timestep number
7465c6c1daeSBarry Smith 
7475c6c1daeSBarry Smith   Level: intermediate
7485c6c1daeSBarry Smith 
7495c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7505c6c1daeSBarry Smith @*/
7515c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7525c6c1daeSBarry Smith {
7535c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7545c6c1daeSBarry Smith 
7555c6c1daeSBarry Smith   PetscFunctionBegin;
7565c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7575c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7585c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7595c6c1daeSBarry Smith   PetscFunctionReturn(0);
7605c6c1daeSBarry Smith }
7615c6c1daeSBarry Smith 
76236bce6e8SMatthew G. Knepley /*@C
76336bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
76436bce6e8SMatthew G. Knepley 
76536bce6e8SMatthew G. Knepley   Not collective
76636bce6e8SMatthew G. Knepley 
76736bce6e8SMatthew G. Knepley   Input Parameter:
76836bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
76936bce6e8SMatthew G. Knepley 
77036bce6e8SMatthew G. Knepley   Output Parameter:
77136bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
77236bce6e8SMatthew G. Knepley 
77336bce6e8SMatthew G. Knepley   Level: advanced
77436bce6e8SMatthew G. Knepley 
77536bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
77636bce6e8SMatthew G. Knepley @*/
77736bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
77836bce6e8SMatthew G. Knepley {
77936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
78036bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
78136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
78236bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
78336bce6e8SMatthew G. Knepley #else
78436bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
78536bce6e8SMatthew G. Knepley #endif
78636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
78736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
78836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
789c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
790de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
79136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
79236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
79336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
7947e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
79536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
79636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
79736bce6e8SMatthew G. Knepley }
79836bce6e8SMatthew G. Knepley 
79936bce6e8SMatthew G. Knepley /*@C
80036bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
80136bce6e8SMatthew G. Knepley 
80236bce6e8SMatthew G. Knepley   Not collective
80336bce6e8SMatthew G. Knepley 
80436bce6e8SMatthew G. Knepley   Input Parameter:
80536bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
80636bce6e8SMatthew G. Knepley 
80736bce6e8SMatthew G. Knepley   Output Parameter:
80836bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
80936bce6e8SMatthew G. Knepley 
81036bce6e8SMatthew G. Knepley   Level: advanced
81136bce6e8SMatthew G. Knepley 
81236bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
81336bce6e8SMatthew G. Knepley @*/
81436bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
81536bce6e8SMatthew G. Knepley {
81636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
81736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
81836bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
81936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
82036bce6e8SMatthew G. Knepley #else
82136bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
82236bce6e8SMatthew G. Knepley #endif
82336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
82436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
82536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
82636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
82736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
82836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8297e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
83036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
83136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
83236bce6e8SMatthew G. Knepley }
83336bce6e8SMatthew G. Knepley 
834df863907SAlex Fikl /*@C
835b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
83636bce6e8SMatthew G. Knepley 
83736bce6e8SMatthew G. Knepley   Input Parameters:
83836bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
83957d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
84036bce6e8SMatthew G. Knepley . name   - The attribute name
84136bce6e8SMatthew G. Knepley . datatype - The attribute type
84236bce6e8SMatthew G. Knepley - value    - The attribute value
84336bce6e8SMatthew G. Knepley 
84436bce6e8SMatthew G. Knepley   Level: advanced
84536bce6e8SMatthew G. Knepley 
846577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
84736bce6e8SMatthew G. Knepley @*/
84857d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value)
84936bce6e8SMatthew G. Knepley {
85057d22f7dSVaclav Hapla   char           *parent;
85160568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
852080f144cSVaclav Hapla   PetscBool      has;
85336bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
85436bce6e8SMatthew G. Knepley 
85536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8565cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
85757d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
858c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
859b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
86057d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
861bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
862b67695ecSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);
86336bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
8647e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
8657e97c476SMatthew G. Knepley     size_t len;
8667e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
867729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
8687e97c476SMatthew G. Knepley   }
86936bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
870729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
87160568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
872080f144cSVaclav Hapla   if (has) {
873080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
874080f144cSVaclav Hapla   } else {
87560568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
876080f144cSVaclav Hapla   }
877729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
878729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
879729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
88060568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
881729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
88257d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
88336bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
88436bce6e8SMatthew G. Knepley }
88536bce6e8SMatthew G. Knepley 
886df863907SAlex Fikl /*@C
887577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
888577e0e04SVaclav Hapla 
889577e0e04SVaclav Hapla   Input Parameters:
890577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
891577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
892577e0e04SVaclav Hapla . name     - The attribute name
893577e0e04SVaclav Hapla . datatype - The attribute type
894577e0e04SVaclav Hapla - value    - The attribute value
895577e0e04SVaclav Hapla 
896577e0e04SVaclav Hapla   Notes:
897577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
898577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
899577e0e04SVaclav Hapla 
900577e0e04SVaclav Hapla   Level: advanced
901577e0e04SVaclav Hapla 
902577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
903577e0e04SVaclav Hapla @*/
904577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
905577e0e04SVaclav Hapla {
906577e0e04SVaclav Hapla   PetscErrorCode ierr;
907577e0e04SVaclav Hapla 
908577e0e04SVaclav Hapla   PetscFunctionBegin;
909577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
910577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
911577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
912b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
913577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
914577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
915577e0e04SVaclav Hapla   PetscFunctionReturn(0);
916577e0e04SVaclav Hapla }
917577e0e04SVaclav Hapla 
918577e0e04SVaclav Hapla /*@C
919b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
920ad0c61feSMatthew G. Knepley 
921ad0c61feSMatthew G. Knepley   Input Parameters:
922ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
92357d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
924ad0c61feSMatthew G. Knepley . name   - The attribute name
925ad0c61feSMatthew G. Knepley - datatype - The attribute type
926ad0c61feSMatthew G. Knepley 
927ad0c61feSMatthew G. Knepley   Output Parameter:
928ad0c61feSMatthew G. Knepley . value    - The attribute value
929ad0c61feSMatthew G. Knepley 
930708d4cb3SBarry Smith   Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
931708d4cb3SBarry Smith 
932ad0c61feSMatthew G. Knepley   Level: advanced
933ad0c61feSMatthew G. Knepley 
934577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
935ad0c61feSMatthew G. Knepley @*/
93657d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value)
937ad0c61feSMatthew G. Knepley {
93857d22f7dSVaclav Hapla   char           *parent;
939f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
940969748fdSVaclav Hapla   PetscBool      has;
941ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
942ad0c61feSMatthew G. Knepley 
943ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9445cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
94557d22f7dSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset, 2);
946c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
947b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
94857d22f7dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
949969748fdSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
950969748fdSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);}
951969748fdSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name);
952ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
953ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
95460568592SMatthew G. Knepley   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT));
95560568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
956f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
957f0b58479SMatthew G. Knepley     size_t len;
958e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
95915b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
960708d4cb3SBarry Smith     ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr);
961708d4cb3SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
962708d4cb3SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
963708d4cb3SBarry Smith   } else {
96470efba33SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
965708d4cb3SBarry Smith   }
966729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
967e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
968e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
96957d22f7dSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
970ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
971ad0c61feSMatthew G. Knepley }
972ad0c61feSMatthew G. Knepley 
973577e0e04SVaclav Hapla /*@C
974577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
975577e0e04SVaclav Hapla 
976577e0e04SVaclav Hapla   Input Parameters:
977577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
978577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
979577e0e04SVaclav Hapla . name     - The attribute name
980577e0e04SVaclav Hapla - datatype - The attribute type
981577e0e04SVaclav Hapla 
982577e0e04SVaclav Hapla   Output Parameter:
983577e0e04SVaclav Hapla . value    - The attribute value
984577e0e04SVaclav Hapla 
985577e0e04SVaclav Hapla   Notes:
986577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
987577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
988577e0e04SVaclav Hapla 
989577e0e04SVaclav Hapla   Level: advanced
990577e0e04SVaclav Hapla 
991577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
992577e0e04SVaclav Hapla @*/
993577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
994577e0e04SVaclav Hapla {
995577e0e04SVaclav Hapla   PetscErrorCode ierr;
996577e0e04SVaclav Hapla 
997577e0e04SVaclav Hapla   PetscFunctionBegin;
998577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
999577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1000577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1001b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
1002577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1003577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
1004577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1005577e0e04SVaclav Hapla }
1006577e0e04SVaclav Hapla 
1007a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1008a75c3fd4SVaclav Hapla {
1009a75c3fd4SVaclav Hapla   htri_t exists;
1010a75c3fd4SVaclav Hapla   hid_t group;
1011a75c3fd4SVaclav Hapla 
1012a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1013a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1014a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1015a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1016a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1017a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
1018a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1019a75c3fd4SVaclav Hapla   }
1020a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1021a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1022a75c3fd4SVaclav Hapla }
1023a75c3fd4SVaclav Hapla 
1024bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
10255cdeb108SMatthew G. Knepley {
102690e03537SVaclav Hapla   const char     rootGroupName[] = "/";
10275cdeb108SMatthew G. Knepley   hid_t          h5;
1028e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
102915b861d2SVaclav Hapla   PetscInt       i;
103015b861d2SVaclav Hapla   int            n;
103185688ae2SVaclav Hapla   char           **hierarchy;
103285688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10335cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10345cdeb108SMatthew G. Knepley 
10355cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10365cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
103790e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
103890e03537SVaclav Hapla   else name = rootGroupName;
1039ccd66a83SVaclav Hapla   if (has) {
104056cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10415cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1042ccd66a83SVaclav Hapla   }
1043ccd66a83SVaclav Hapla   if (otype) {
1044ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
104556cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1046ccd66a83SVaclav Hapla   }
1047c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
104885688ae2SVaclav Hapla 
104985688ae2SVaclav Hapla   /*
105085688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
105185688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
105285688ae2SVaclav Hapla      1) whether it's a valid link
105385688ae2SVaclav Hapla      2) whether this link resolves to an object
105485688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
105585688ae2SVaclav Hapla   */
105685688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
105785688ae2SVaclav Hapla   if (!n) {
105885688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1059ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1060ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
106185688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
106285688ae2SVaclav Hapla     PetscFunctionReturn(0);
106385688ae2SVaclav Hapla   }
106485688ae2SVaclav Hapla   for (i=0; i<n; i++) {
106585688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
106685688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1067a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1068a75c3fd4SVaclav Hapla     if (!exists) break;
106985688ae2SVaclav Hapla   }
107085688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
107185688ae2SVaclav Hapla 
107285688ae2SVaclav Hapla   /* If the object exists, get its type */
1073ccd66a83SVaclav Hapla   if (exists && otype) {
10745cdeb108SMatthew G. Knepley     H5O_info_t info;
10755cdeb108SMatthew G. Knepley 
107676276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
107704633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
107856cc0592SVaclav Hapla     *otype = info.type;
10795cdeb108SMatthew G. Knepley   }
1080ccd66a83SVaclav Hapla   if (has) *has = exists;
10815cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
10825cdeb108SMatthew G. Knepley }
10835cdeb108SMatthew G. Knepley 
1084ecce1506SVaclav Hapla /*@
10850a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
10860a9f38efSVaclav Hapla 
10870a9f38efSVaclav Hapla   Input Parameters:
10880a9f38efSVaclav Hapla . viewer - The HDF5 viewer
10890a9f38efSVaclav Hapla 
10900a9f38efSVaclav Hapla   Output Parameter:
10910a9f38efSVaclav Hapla . has    - Flag for group existence
10920a9f38efSVaclav Hapla 
10930a9f38efSVaclav Hapla   Notes:
10940a9f38efSVaclav Hapla   If the path exists but is not a group, this returns PETSC_FALSE as well.
10950a9f38efSVaclav Hapla 
10960a9f38efSVaclav Hapla   Level: advanced
10970a9f38efSVaclav Hapla 
10980a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup()
10990a9f38efSVaclav Hapla @*/
11000a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has)
11010a9f38efSVaclav Hapla {
11020a9f38efSVaclav Hapla   H5O_type_t type;
11030a9f38efSVaclav Hapla   const char *name;
11040a9f38efSVaclav Hapla   PetscErrorCode ierr;
11050a9f38efSVaclav Hapla 
11060a9f38efSVaclav Hapla   PetscFunctionBegin;
11070a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11080a9f38efSVaclav Hapla   PetscValidIntPointer(has,2);
11090a9f38efSVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr);
11100a9f38efSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr);
11110a9f38efSVaclav Hapla   *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE;
11120a9f38efSVaclav Hapla   PetscFunctionReturn(0);
11130a9f38efSVaclav Hapla }
11140a9f38efSVaclav Hapla 
11150a9f38efSVaclav Hapla /*@
1116e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1117ecce1506SVaclav Hapla 
1118ecce1506SVaclav Hapla   Input Parameters:
1119ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1120ecce1506SVaclav Hapla - obj    - The named object
1121ecce1506SVaclav Hapla 
1122ecce1506SVaclav Hapla   Output Parameter:
1123ecce1506SVaclav Hapla . has    - Flag for dataset existence; PETSC_FALSE for unnamed object
1124ecce1506SVaclav Hapla 
1125e3f143f7SVaclav Hapla   Notes:
1126e3f143f7SVaclav Hapla   If the path exists but is not a dataset, this returns PETSC_FALSE as well.
1127e3f143f7SVaclav Hapla 
1128ecce1506SVaclav Hapla   Level: advanced
1129ecce1506SVaclav Hapla 
1130e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1131ecce1506SVaclav Hapla @*/
1132ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1133ecce1506SVaclav Hapla {
113456cc0592SVaclav Hapla   H5O_type_t type;
11356c132bc1SVaclav Hapla   char *path;
1136ecce1506SVaclav Hapla   PetscErrorCode ierr;
1137ecce1506SVaclav Hapla 
1138ecce1506SVaclav Hapla   PetscFunctionBegin;
1139c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1140c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
1141c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,3);
1142ecce1506SVaclav Hapla   *has = PETSC_FALSE;
1143e3f143f7SVaclav Hapla   if (!obj->name) PetscFunctionReturn(0);
11446c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr);
11456c132bc1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr);
114656cc0592SVaclav Hapla   *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE;
11476c132bc1SVaclav Hapla   ierr = PetscFree(path);CHKERRQ(ierr);
1148ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1149ecce1506SVaclav Hapla }
1150ecce1506SVaclav Hapla 
1151df863907SAlex Fikl /*@C
1152ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1153e7bdbf8eSMatthew G. Knepley 
1154e7bdbf8eSMatthew G. Knepley   Input Parameters:
1155e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
115610e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute.
1157e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1158e7bdbf8eSMatthew G. Knepley 
1159e7bdbf8eSMatthew G. Knepley   Output Parameter:
1160e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1161e7bdbf8eSMatthew G. Knepley 
1162e7bdbf8eSMatthew G. Knepley   Level: advanced
1163e7bdbf8eSMatthew G. Knepley 
1164577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1165e7bdbf8eSMatthew G. Knepley @*/
116610e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has)
1167e7bdbf8eSMatthew G. Knepley {
116810e69e8fSVaclav Hapla   char           *parent;
1169e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1170e7bdbf8eSMatthew G. Knepley 
1171e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
11725cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
117310e69e8fSVaclav Hapla   if (dataset) PetscValidCharPointer(dataset,2);
1174c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
1175c57a1a9aSVaclav Hapla   PetscValidIntPointer(has,4);
117610e69e8fSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr);
1177bb286ee1SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
117810e69e8fSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);}
117910e69e8fSVaclav Hapla   ierr = PetscFree(parent);CHKERRQ(ierr);
118006db490cSVaclav Hapla   PetscFunctionReturn(0);
118106db490cSVaclav Hapla }
118206db490cSVaclav Hapla 
1183577e0e04SVaclav Hapla /*@C
1184577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1185577e0e04SVaclav Hapla 
1186577e0e04SVaclav Hapla   Input Parameters:
1187577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1188577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1189577e0e04SVaclav Hapla - name   - The attribute name
1190577e0e04SVaclav Hapla 
1191577e0e04SVaclav Hapla   Output Parameter:
1192577e0e04SVaclav Hapla . has    - Flag for attribute existence
1193577e0e04SVaclav Hapla 
1194577e0e04SVaclav Hapla   Notes:
1195577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1196577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1197577e0e04SVaclav Hapla 
1198577e0e04SVaclav Hapla   Level: advanced
1199577e0e04SVaclav Hapla 
1200577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1201577e0e04SVaclav Hapla @*/
1202577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1203577e0e04SVaclav Hapla {
1204577e0e04SVaclav Hapla   PetscErrorCode ierr;
1205577e0e04SVaclav Hapla 
1206577e0e04SVaclav Hapla   PetscFunctionBegin;
1207577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1208577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1209577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1210577e0e04SVaclav Hapla   PetscValidIntPointer(has,4);
1211577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1212577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1213577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1214577e0e04SVaclav Hapla }
1215577e0e04SVaclav Hapla 
121606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
121706db490cSVaclav Hapla {
121806db490cSVaclav Hapla   hid_t          h5;
121906db490cSVaclav Hapla   htri_t         hhas;
122006db490cSVaclav Hapla   PetscErrorCode ierr;
122106db490cSVaclav Hapla 
122206db490cSVaclav Hapla   PetscFunctionBegin;
122306db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
12242f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
12255cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1226e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1227e7bdbf8eSMatthew G. Knepley }
1228e7bdbf8eSMatthew G. Knepley 
1229a75e6a4aSMatthew G. Knepley /*
1230a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1231a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1232a75e6a4aSMatthew G. Knepley */
1233d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1234a75e6a4aSMatthew G. Knepley 
1235a75e6a4aSMatthew G. Knepley /*@C
1236a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1237a75e6a4aSMatthew G. Knepley 
1238d083f849SBarry Smith   Collective
1239a75e6a4aSMatthew G. Knepley 
1240a75e6a4aSMatthew G. Knepley   Input Parameter:
1241a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1242a75e6a4aSMatthew G. Knepley 
1243a75e6a4aSMatthew G. Knepley   Level: intermediate
1244a75e6a4aSMatthew G. Knepley 
1245a75e6a4aSMatthew G. Knepley   Options Database Keys:
1246a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1247a75e6a4aSMatthew G. Knepley 
1248a75e6a4aSMatthew G. Knepley   Environmental variables:
1249a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1250a75e6a4aSMatthew G. Knepley 
1251a75e6a4aSMatthew G. Knepley   Notes:
1252a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1253a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1254a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1255a75e6a4aSMatthew G. Knepley 
1256a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1257a75e6a4aSMatthew G. Knepley @*/
1258a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1259a75e6a4aSMatthew G. Knepley {
1260a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1261a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1262a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1263a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1264a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1265a75e6a4aSMatthew G. Knepley 
1266a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1267908793a3SLisandro Dalcin   ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1268a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1269908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1270908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1271a75e6a4aSMatthew G. Knepley   }
127247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1273908793a3SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1274a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1275a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
12762cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1277a75e6a4aSMatthew G. Knepley     if (!flg) {
1278a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
12792cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1280a75e6a4aSMatthew G. Knepley     }
1281a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
12822cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1283a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
12842cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
128547435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1286908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1287a75e6a4aSMatthew G. Knepley   }
1288a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
12892cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1290a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1291a75e6a4aSMatthew G. Knepley }
1292