xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 89e0ef10bc51934e7d3e7d7ba7751b56bc779fa0)
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 
84302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[])
96c132bc1SVaclav Hapla {
104302210dSVaclav Hapla   PetscBool      relative = PETSC_FALSE;
116c132bc1SVaclav Hapla   const char     *group;
126c132bc1SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN] = "";
136c132bc1SVaclav Hapla   PetscErrorCode ierr;
146c132bc1SVaclav Hapla 
156c132bc1SVaclav Hapla   PetscFunctionBegin;
166c132bc1SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
174302210dSVaclav Hapla   ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr);
184302210dSVaclav Hapla   if (relative) {
194302210dSVaclav Hapla     ierr = PetscStrcpy(buf, group);CHKERRQ(ierr);
206c132bc1SVaclav Hapla     ierr = PetscStrcat(buf, "/");CHKERRQ(ierr);
214302210dSVaclav Hapla     ierr = PetscStrcat(buf, path);CHKERRQ(ierr);
224302210dSVaclav Hapla     ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr);
234302210dSVaclav Hapla   } else {
244302210dSVaclav Hapla     ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr);
254302210dSVaclav Hapla   }
266c132bc1SVaclav Hapla   PetscFunctionReturn(0);
276c132bc1SVaclav Hapla }
286c132bc1SVaclav Hapla 
29577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
30577e0e04SVaclav Hapla {
31577e0e04SVaclav Hapla   PetscBool      has;
32577e0e04SVaclav Hapla   PetscErrorCode ierr;
33577e0e04SVaclav Hapla 
34577e0e04SVaclav Hapla   PetscFunctionBegin;
35577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr);
36577e0e04SVaclav Hapla   if (!has) {
37*89e0ef10SVaclav Hapla     const char *group;
38577e0e04SVaclav Hapla     ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr);
39*89e0ef10SVaclav Hapla     SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
40577e0e04SVaclav Hapla   }
41577e0e04SVaclav Hapla   PetscFunctionReturn(0);
42577e0e04SVaclav Hapla }
43577e0e04SVaclav Hapla 
444416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
4582ea9b62SBarry Smith {
4682ea9b62SBarry Smith   PetscErrorCode   ierr;
47ee8b9732SVaclav Hapla   PetscBool        flg = PETSC_FALSE, set;
4882ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
4982ea9b62SBarry Smith 
5082ea9b62SBarry Smith   PetscFunctionBegin;
5182ea9b62SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr);
5282ea9b62SBarry Smith   ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr);
539a0502c6SHåkon Strandenes   ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr);
54ee8b9732SVaclav Hapla   ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr);
55ee8b9732SVaclav Hapla   if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);}
5682ea9b62SBarry Smith   ierr = PetscOptionsTail();CHKERRQ(ierr);
5782ea9b62SBarry Smith   PetscFunctionReturn(0);
5882ea9b62SBarry Smith }
5982ea9b62SBarry Smith 
601b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
611b793a25SVaclav Hapla {
621b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
6303fa8834SVaclav Hapla   PetscBool         flg;
641b793a25SVaclav Hapla   PetscErrorCode    ierr;
651b793a25SVaclav Hapla 
661b793a25SVaclav Hapla   PetscFunctionBegin;
671b793a25SVaclav Hapla   if (hdf5->filename) {
681b793a25SVaclav Hapla     ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr);
691b793a25SVaclav Hapla   }
701b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr);
711b793a25SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr);
7203fa8834SVaclav Hapla   ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr);
7303fa8834SVaclav Hapla   ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr);
741b793a25SVaclav Hapla   PetscFunctionReturn(0);
751b793a25SVaclav Hapla }
761b793a25SVaclav Hapla 
775c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
785c6c1daeSBarry Smith {
795c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
805c6c1daeSBarry Smith   PetscErrorCode   ierr;
815c6c1daeSBarry Smith 
825c6c1daeSBarry Smith   PetscFunctionBegin;
835c6c1daeSBarry Smith   ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);
84729ab6d9SBarry Smith   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
855c6c1daeSBarry Smith   PetscFunctionReturn(0);
865c6c1daeSBarry Smith }
875c6c1daeSBarry Smith 
889b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
895c6c1daeSBarry Smith {
905c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
915c6c1daeSBarry Smith   PetscErrorCode   ierr;
925c6c1daeSBarry Smith 
935c6c1daeSBarry Smith   PetscFunctionBegin;
949c5072fbSVaclav Hapla   PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id));
955c6c1daeSBarry Smith   ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr);
965c6c1daeSBarry Smith   while (hdf5->groups) {
977d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
985c6c1daeSBarry Smith 
995c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups->name);CHKERRQ(ierr);
1005c6c1daeSBarry Smith     ierr         = PetscFree(hdf5->groups);CHKERRQ(ierr);
1015c6c1daeSBarry Smith     hdf5->groups = tmp;
1025c6c1daeSBarry Smith   }
1035c6c1daeSBarry Smith   ierr = PetscFree(hdf5);CHKERRQ(ierr);
1040b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr);
105d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr);
1060b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr);
107058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr);
108058bd781SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr);
1095c6c1daeSBarry Smith   PetscFunctionReturn(0);
1105c6c1daeSBarry Smith }
1115c6c1daeSBarry Smith 
1129b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1135c6c1daeSBarry Smith {
1145c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1155c6c1daeSBarry Smith 
1165c6c1daeSBarry Smith   PetscFunctionBegin;
1175c6c1daeSBarry Smith   hdf5->btype = type;
1185c6c1daeSBarry Smith   PetscFunctionReturn(0);
1195c6c1daeSBarry Smith }
1205c6c1daeSBarry Smith 
1210b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1220b62783dSVaclav Hapla {
1230b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1240b62783dSVaclav Hapla 
1250b62783dSVaclav Hapla   PetscFunctionBegin;
1260b62783dSVaclav Hapla   *type = hdf5->btype;
1270b62783dSVaclav Hapla   PetscFunctionReturn(0);
1280b62783dSVaclav Hapla }
1290b62783dSVaclav Hapla 
1309b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
13182ea9b62SBarry Smith {
13282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
13382ea9b62SBarry Smith 
13482ea9b62SBarry Smith   PetscFunctionBegin;
13582ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13682ea9b62SBarry Smith   PetscFunctionReturn(0);
13782ea9b62SBarry Smith }
13882ea9b62SBarry Smith 
139df863907SAlex Fikl /*@
14082ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
14182ea9b62SBarry Smith        dimension of 2.
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith     Logically Collective on PetscViewer
14482ea9b62SBarry Smith 
14582ea9b62SBarry Smith   Input Parameters:
14682ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
14782ea9b62SBarry 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
14882ea9b62SBarry Smith 
14982ea9b62SBarry Smith   Options Database:
15082ea9b62SBarry 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
15182ea9b62SBarry Smith 
15282ea9b62SBarry Smith 
15395452b02SPatrick Sanan   Notes:
15495452b02SPatrick 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
15582ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
15682ea9b62SBarry Smith 
15782ea9b62SBarry Smith   Level: intermediate
15882ea9b62SBarry Smith 
15982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
16082ea9b62SBarry Smith 
16182ea9b62SBarry Smith @*/
16282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
16382ea9b62SBarry Smith {
16482ea9b62SBarry Smith   PetscErrorCode ierr;
16582ea9b62SBarry Smith 
16682ea9b62SBarry Smith   PetscFunctionBegin;
16782ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
16882ea9b62SBarry Smith   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
16982ea9b62SBarry Smith   PetscFunctionReturn(0);
17082ea9b62SBarry Smith }
17182ea9b62SBarry Smith 
172df863907SAlex Fikl /*@
17382ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
17482ea9b62SBarry Smith        dimension of 2.
17582ea9b62SBarry Smith 
17682ea9b62SBarry Smith     Logically Collective on PetscViewer
17782ea9b62SBarry Smith 
17882ea9b62SBarry Smith   Input Parameter:
17982ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
18082ea9b62SBarry Smith 
18182ea9b62SBarry Smith   Output Parameter:
18282ea9b62SBarry 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
18382ea9b62SBarry Smith 
18495452b02SPatrick Sanan   Notes:
18595452b02SPatrick 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
18682ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
18782ea9b62SBarry Smith 
18882ea9b62SBarry Smith   Level: intermediate
18982ea9b62SBarry Smith 
19082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen()
19182ea9b62SBarry Smith 
19282ea9b62SBarry Smith @*/
19382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
19482ea9b62SBarry Smith {
19582ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
19682ea9b62SBarry Smith 
19782ea9b62SBarry Smith   PetscFunctionBegin;
19882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
19982ea9b62SBarry Smith   *flg = hdf5->basedimension2;
20082ea9b62SBarry Smith   PetscFunctionReturn(0);
20182ea9b62SBarry Smith }
20282ea9b62SBarry Smith 
2039b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2049a0502c6SHåkon Strandenes {
2059a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2069a0502c6SHåkon Strandenes 
2079a0502c6SHåkon Strandenes   PetscFunctionBegin;
2089a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2099a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2109a0502c6SHåkon Strandenes }
2119a0502c6SHåkon Strandenes 
212df863907SAlex Fikl /*@
2139a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2149a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2159a0502c6SHåkon Strandenes 
2169a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2179a0502c6SHåkon Strandenes 
2189a0502c6SHåkon Strandenes   Input Parameters:
2199a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2209a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2219a0502c6SHåkon Strandenes 
2229a0502c6SHåkon Strandenes   Options Database:
2239a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2249a0502c6SHåkon Strandenes 
2259a0502c6SHåkon Strandenes 
22695452b02SPatrick Sanan   Notes:
22795452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2289a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2299a0502c6SHåkon Strandenes 
2309a0502c6SHåkon Strandenes   Level: intermediate
2319a0502c6SHåkon Strandenes 
2329a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2339a0502c6SHåkon Strandenes           PetscReal
2349a0502c6SHåkon Strandenes 
2359a0502c6SHåkon Strandenes @*/
2369a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2379a0502c6SHåkon Strandenes {
2389a0502c6SHåkon Strandenes   PetscErrorCode ierr;
2399a0502c6SHåkon Strandenes 
2409a0502c6SHåkon Strandenes   PetscFunctionBegin;
2419a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2429a0502c6SHåkon Strandenes   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
2439a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2449a0502c6SHåkon Strandenes }
2459a0502c6SHåkon Strandenes 
246df863907SAlex Fikl /*@
2479a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2489a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2499a0502c6SHåkon Strandenes 
2509a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2519a0502c6SHåkon Strandenes 
2529a0502c6SHåkon Strandenes   Input Parameter:
2539a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2549a0502c6SHåkon Strandenes 
2559a0502c6SHåkon Strandenes   Output Parameter:
2569a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2579a0502c6SHåkon Strandenes 
25895452b02SPatrick Sanan   Notes:
25995452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2609a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2619a0502c6SHåkon Strandenes 
2629a0502c6SHåkon Strandenes   Level: intermediate
2639a0502c6SHåkon Strandenes 
2649a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(),
2659a0502c6SHåkon Strandenes           PetscReal
2669a0502c6SHåkon Strandenes 
2679a0502c6SHåkon Strandenes @*/
2689a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2699a0502c6SHåkon Strandenes {
2709a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2719a0502c6SHåkon Strandenes 
2729a0502c6SHåkon Strandenes   PetscFunctionBegin;
2739a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2749a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2759a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2769a0502c6SHåkon Strandenes }
2779a0502c6SHåkon Strandenes 
278ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
279ee8b9732SVaclav Hapla {
280ee8b9732SVaclav Hapla   PetscFunctionBegin;
281ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
282ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
283c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
284ee8b9732SVaclav Hapla   {
285ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
286ee8b9732SVaclav Hapla     PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
287ee8b9732SVaclav Hapla   }
288ee8b9732SVaclav Hapla #else
289ee8b9732SVaclav Hapla   if (flg) {
290ee8b9732SVaclav Hapla     PetscErrorCode ierr;
291c796909eSBarry 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);
292ee8b9732SVaclav Hapla   }
293ee8b9732SVaclav Hapla #endif
294ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
295ee8b9732SVaclav Hapla }
296ee8b9732SVaclav Hapla 
297ee8b9732SVaclav Hapla /*@
298ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
299ee8b9732SVaclav Hapla 
300ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
301ee8b9732SVaclav Hapla 
302ee8b9732SVaclav Hapla   Input Parameters:
303ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
304ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
305ee8b9732SVaclav Hapla 
306ee8b9732SVaclav Hapla   Options Database:
307ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
308ee8b9732SVaclav Hapla 
309ee8b9732SVaclav Hapla   Notes:
310ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
31153effed7SBarry 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.
312ee8b9732SVaclav Hapla 
313ee8b9732SVaclav Hapla   Developer notes:
314ee8b9732SVaclav Hapla   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
315ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
316ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
317ee8b9732SVaclav Hapla 
318ee8b9732SVaclav Hapla   Level: intermediate
319ee8b9732SVaclav Hapla 
320ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
321ee8b9732SVaclav Hapla 
322ee8b9732SVaclav Hapla @*/
323ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
324ee8b9732SVaclav Hapla {
325ee8b9732SVaclav Hapla   PetscErrorCode ierr;
326ee8b9732SVaclav Hapla 
327ee8b9732SVaclav Hapla   PetscFunctionBegin;
328ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
329ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer,flg,2);
330ee8b9732SVaclav Hapla   ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr);
331ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
332ee8b9732SVaclav Hapla }
333ee8b9732SVaclav Hapla 
334ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
335ee8b9732SVaclav Hapla {
336c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
337ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
338ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
339c796909eSBarry Smith #endif
340ee8b9732SVaclav Hapla 
341ee8b9732SVaclav Hapla   PetscFunctionBegin;
342c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
343c796909eSBarry Smith   *flg = PETSC_FALSE;
344c796909eSBarry Smith #else
345ee8b9732SVaclav Hapla   PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
346ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
347c796909eSBarry Smith #endif
348ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
349ee8b9732SVaclav Hapla }
350ee8b9732SVaclav Hapla 
351ee8b9732SVaclav Hapla /*@
352ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
353ee8b9732SVaclav Hapla 
354ee8b9732SVaclav Hapla   Not Collective
355ee8b9732SVaclav Hapla 
356ee8b9732SVaclav Hapla   Input Parameters:
357ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
358ee8b9732SVaclav Hapla 
359ee8b9732SVaclav Hapla   Output Parameters:
360ee8b9732SVaclav Hapla . flg - the flag
361ee8b9732SVaclav Hapla 
362ee8b9732SVaclav Hapla   Level: intermediate
363ee8b9732SVaclav Hapla 
364ee8b9732SVaclav Hapla   Notes:
365c796909eSBarry 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.
366ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
367ee8b9732SVaclav Hapla 
368ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open()
369ee8b9732SVaclav Hapla 
370ee8b9732SVaclav Hapla @*/
371ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
372ee8b9732SVaclav Hapla {
373ee8b9732SVaclav Hapla   PetscErrorCode ierr;
374ee8b9732SVaclav Hapla 
375ee8b9732SVaclav Hapla   PetscFunctionBegin;
376ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
377534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
378ee8b9732SVaclav Hapla 
379ee8b9732SVaclav Hapla   ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr);
380ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
381ee8b9732SVaclav Hapla }
382ee8b9732SVaclav Hapla 
3839b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3845c6c1daeSBarry Smith {
3855c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3865c6c1daeSBarry Smith   hid_t             plist_id;
3875c6c1daeSBarry Smith   PetscErrorCode    ierr;
3885c6c1daeSBarry Smith 
3895c6c1daeSBarry Smith   PetscFunctionBegin;
390f683cc58SToby Isaac   if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id));
391f683cc58SToby Isaac   if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);}
3925c6c1daeSBarry Smith   ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr);
3935c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
394729ab6d9SBarry Smith   PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
395c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
396d28bfa55SVaclav Hapla   PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
397c796909eSBarry Smith #endif
3985c6c1daeSBarry Smith   /* Create or open the file collectively */
3995c6c1daeSBarry Smith   switch (hdf5->btype) {
4005c6c1daeSBarry Smith   case FILE_MODE_READ:
401729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
4025c6c1daeSBarry Smith     break;
4035c6c1daeSBarry Smith   case FILE_MODE_APPEND:
4047e4fd573SVaclav Hapla   case FILE_MODE_UPDATE:
4057e4fd573SVaclav Hapla   {
4067e4fd573SVaclav Hapla     PetscBool flg;
4077e4fd573SVaclav Hapla     ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr);
4087e4fd573SVaclav Hapla     if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
4097e4fd573SVaclav Hapla     else     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4105c6c1daeSBarry Smith     break;
4117e4fd573SVaclav Hapla   }
4125c6c1daeSBarry Smith   case FILE_MODE_WRITE:
413729ab6d9SBarry Smith     PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
4145c6c1daeSBarry Smith     break;
4157e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED:
4167e4fd573SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
4175c6c1daeSBarry Smith   default:
4187e4fd573SVaclav Hapla     SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
4195c6c1daeSBarry Smith   }
4205c6c1daeSBarry Smith   if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
421729ab6d9SBarry Smith   PetscStackCallHDF5(H5Pclose,(plist_id));
4225c6c1daeSBarry Smith   PetscFunctionReturn(0);
4235c6c1daeSBarry Smith }
4245c6c1daeSBarry Smith 
425d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
426d1232d7fSToby Isaac {
427d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
428d1232d7fSToby Isaac 
429d1232d7fSToby Isaac   PetscFunctionBegin;
430d1232d7fSToby Isaac   *name = vhdf5->filename;
431d1232d7fSToby Isaac   PetscFunctionReturn(0);
432d1232d7fSToby Isaac }
433d1232d7fSToby Isaac 
434b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
435b723ab35SVaclav Hapla {
4365dc64a97SVaclav Hapla   /*
437b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
438b723ab35SVaclav Hapla   PetscErrorCode   ierr;
4395dc64a97SVaclav Hapla   */
440b723ab35SVaclav Hapla 
441b723ab35SVaclav Hapla   PetscFunctionBegin;
442b723ab35SVaclav Hapla   PetscFunctionReturn(0);
443b723ab35SVaclav Hapla }
444b723ab35SVaclav Hapla 
4458556b5ebSBarry Smith /*MC
4468556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4478556b5ebSBarry Smith 
4488556b5ebSBarry Smith 
4498556b5ebSBarry Smith .seealso:  PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
4508556b5ebSBarry Smith            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
4518556b5ebSBarry Smith            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
4528556b5ebSBarry Smith            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
4538556b5ebSBarry Smith 
4541b266c99SBarry Smith   Level: beginner
4558556b5ebSBarry Smith M*/
456d1232d7fSToby Isaac 
4578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
4585c6c1daeSBarry Smith {
4595c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4605c6c1daeSBarry Smith   PetscErrorCode   ierr;
4615c6c1daeSBarry Smith 
4625c6c1daeSBarry Smith   PetscFunctionBegin;
463b00a9115SJed Brown   ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr);
4645c6c1daeSBarry Smith 
4655c6c1daeSBarry Smith   v->data                = (void*) hdf5;
4665c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
46782ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
468b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
4691b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
470908793a3SLisandro Dalcin   v->ops->flush          = NULL;
4717e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
472908793a3SLisandro Dalcin   hdf5->filename         = NULL;
4735c6c1daeSBarry Smith   hdf5->timestep         = -1;
4740298fd71SBarry Smith   hdf5->groups           = NULL;
4755c6c1daeSBarry Smith 
4769c5072fbSVaclav Hapla   PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
4779c5072fbSVaclav Hapla 
4780b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr);
479d1232d7fSToby Isaac   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr);
4800b062f91SJed Brown   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr);
4810b62783dSVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr);
48282ea9b62SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr);
4839a0502c6SHåkon Strandenes   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr);
484ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr);
485ee8b9732SVaclav Hapla   ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr);
4865c6c1daeSBarry Smith   PetscFunctionReturn(0);
4875c6c1daeSBarry Smith }
4885c6c1daeSBarry Smith 
4895c6c1daeSBarry Smith /*@C
4905c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
4915c6c1daeSBarry Smith 
492d083f849SBarry Smith    Collective
4935c6c1daeSBarry Smith 
4945c6c1daeSBarry Smith    Input Parameters:
4955c6c1daeSBarry Smith +  comm - MPI communicator
4965c6c1daeSBarry Smith .  name - name of file
4975c6c1daeSBarry Smith -  type - type of file
4985c6c1daeSBarry Smith 
4995c6c1daeSBarry Smith    Output Parameter:
5005c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
5015c6c1daeSBarry Smith 
50282ea9b62SBarry Smith   Options Database:
503a2b725a8SWilliam 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
504a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
50582ea9b62SBarry Smith 
5065c6c1daeSBarry Smith    Level: beginner
5075c6c1daeSBarry Smith 
5087e4fd573SVaclav Hapla    Notes:
5097e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
5107e4fd573SVaclav Hapla +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
5117e4fd573SVaclav Hapla .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
5127e4fd573SVaclav 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]
5137e4fd573SVaclav Hapla -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
5147e4fd573SVaclav Hapla 
5157e4fd573SVaclav 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.
5167e4fd573SVaclav Hapla 
5175c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5185c6c1daeSBarry Smith 
5195c6c1daeSBarry Smith 
5206a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(),
5219a0502c6SHåkon Strandenes           PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(),
522a56f64adSBarry Smith           MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName()
5235c6c1daeSBarry Smith @*/
5245c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
5255c6c1daeSBarry Smith {
5265c6c1daeSBarry Smith   PetscErrorCode ierr;
5275c6c1daeSBarry Smith 
5285c6c1daeSBarry Smith   PetscFunctionBegin;
5295c6c1daeSBarry Smith   ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr);
5305c6c1daeSBarry Smith   ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr);
5315c6c1daeSBarry Smith   ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr);
5325c6c1daeSBarry Smith   ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr);
533b663d2d2SVaclav Hapla   ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr);
5345c6c1daeSBarry Smith   PetscFunctionReturn(0);
5355c6c1daeSBarry Smith }
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith /*@C
5385c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5395c6c1daeSBarry Smith 
5405c6c1daeSBarry Smith   Not collective
5415c6c1daeSBarry Smith 
5425c6c1daeSBarry Smith   Input Parameter:
5435c6c1daeSBarry Smith . viewer - the PetscViewer
5445c6c1daeSBarry Smith 
5455c6c1daeSBarry Smith   Output Parameter:
5465c6c1daeSBarry Smith . file_id - The file id
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith   Level: intermediate
5495c6c1daeSBarry Smith 
5505c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open()
5515c6c1daeSBarry Smith @*/
5525c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
5535c6c1daeSBarry Smith {
5545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
5555c6c1daeSBarry Smith 
5565c6c1daeSBarry Smith   PetscFunctionBegin;
5575c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
5585c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5595c6c1daeSBarry Smith   PetscFunctionReturn(0);
5605c6c1daeSBarry Smith }
5615c6c1daeSBarry Smith 
5625c6c1daeSBarry Smith /*@C
5635c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5645c6c1daeSBarry Smith 
5655c6c1daeSBarry Smith   Not collective
5665c6c1daeSBarry Smith 
5675c6c1daeSBarry Smith   Input Parameters:
5685c6c1daeSBarry Smith + viewer - the PetscViewer
5695c6c1daeSBarry Smith - name - The group name
5705c6c1daeSBarry Smith 
5715c6c1daeSBarry Smith   Level: intermediate
5725c6c1daeSBarry Smith 
5732e361344SVaclav Hapla   Notes:
5742e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
5752e361344SVaclav Hapla   + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group.
5762e361344SVaclav Hapla   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
5772e361344SVaclav Hapla   - "." means the current group is pushed again.
5782e361344SVaclav Hapla 
5792e361344SVaclav Hapla   Example:
5802e361344SVaclav Hapla   Suppose the current group is "/a".
5812e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
5822e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
5832e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
5842e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
5852e361344SVaclav Hapla 
5862e361344SVaclav Hapla   Developer Notes:
5872e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
588820fc2d1SVaclav Hapla 
589874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
5905c6c1daeSBarry Smith @*/
591be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
5925c6c1daeSBarry Smith {
5935c6c1daeSBarry Smith   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
5947d964842SVaclav Hapla   PetscViewerHDF5GroupList  *groupNode;
5952e361344SVaclav Hapla   size_t                    i,len;
5962e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
5972e361344SVaclav Hapla   const char                *gname;
5985c6c1daeSBarry Smith   PetscErrorCode            ierr;
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith   PetscFunctionBegin;
6015c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
602820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name,2);
603820fc2d1SVaclav Hapla   ierr = PetscStrlen(name, &len);CHKERRQ(ierr);
6042e361344SVaclav Hapla   gname = NULL;
6052e361344SVaclav Hapla   if (len) {
6062e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
6072e361344SVaclav Hapla       /* use current name */
6082e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
6092e361344SVaclav Hapla     } else if (name[0] == '/') {
6102e361344SVaclav Hapla       /* absolute */
6112e361344SVaclav Hapla       for (i=1; i<len; i++) {
6122e361344SVaclav Hapla         if (name[i] != '/') {
6132e361344SVaclav Hapla           gname = name;
6142e361344SVaclav Hapla           break;
6152e361344SVaclav Hapla         }
6162e361344SVaclav Hapla       }
6172e361344SVaclav Hapla     } else {
6182e361344SVaclav Hapla       /* relative */
6192e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
6202e361344SVaclav Hapla       ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr);
6212e361344SVaclav Hapla       gname = buf;
6222e361344SVaclav Hapla     }
6232e361344SVaclav Hapla   }
62495dccacaSBarry Smith   ierr = PetscNew(&groupNode);CHKERRQ(ierr);
6252e361344SVaclav Hapla   ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr);
6265c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
6275c6c1daeSBarry Smith   hdf5->groups    = groupNode;
6285c6c1daeSBarry Smith   PetscFunctionReturn(0);
6295c6c1daeSBarry Smith }
6305c6c1daeSBarry Smith 
6313ef9c667SSatish Balay /*@
6325c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
6335c6c1daeSBarry Smith 
6345c6c1daeSBarry Smith   Not collective
6355c6c1daeSBarry Smith 
6365c6c1daeSBarry Smith   Input Parameter:
6375c6c1daeSBarry Smith . viewer - the PetscViewer
6385c6c1daeSBarry Smith 
6395c6c1daeSBarry Smith   Level: intermediate
6405c6c1daeSBarry Smith 
641874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup()
6425c6c1daeSBarry Smith @*/
6435c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
6445c6c1daeSBarry Smith {
6455c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6467d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6475c6c1daeSBarry Smith   PetscErrorCode   ierr;
6485c6c1daeSBarry Smith 
6495c6c1daeSBarry Smith   PetscFunctionBegin;
6505c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
65182f516ccSBarry Smith   if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6525c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6535c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6545c6c1daeSBarry Smith   ierr         = PetscFree(groupNode->name);CHKERRQ(ierr);
6555c6c1daeSBarry Smith   ierr         = PetscFree(groupNode);CHKERRQ(ierr);
6565c6c1daeSBarry Smith   PetscFunctionReturn(0);
6575c6c1daeSBarry Smith }
6585c6c1daeSBarry Smith 
6595c6c1daeSBarry Smith /*@C
660874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
661874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6625c6c1daeSBarry Smith 
6635c6c1daeSBarry Smith   Not collective
6645c6c1daeSBarry Smith 
6655c6c1daeSBarry Smith   Input Parameter:
6665c6c1daeSBarry Smith . viewer - the PetscViewer
6675c6c1daeSBarry Smith 
6685c6c1daeSBarry Smith   Output Parameter:
6695c6c1daeSBarry Smith . name - The group name
6705c6c1daeSBarry Smith 
6715c6c1daeSBarry Smith   Level: intermediate
6725c6c1daeSBarry Smith 
673874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup()
6745c6c1daeSBarry Smith @*/
675be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
6765c6c1daeSBarry Smith {
6775c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
6785c6c1daeSBarry Smith 
6795c6c1daeSBarry Smith   PetscFunctionBegin;
6805c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
681c959eef4SJed Brown   PetscValidPointer(name,2);
682a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
6830298fd71SBarry Smith   else *name = NULL;
6845c6c1daeSBarry Smith   PetscFunctionReturn(0);
6855c6c1daeSBarry Smith }
6865c6c1daeSBarry Smith 
6878c557b5aSVaclav Hapla /*@
688874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
689874270d9SVaclav Hapla   and return this group's ID and file ID.
690874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
691874270d9SVaclav Hapla 
692874270d9SVaclav Hapla   Not collective
693874270d9SVaclav Hapla 
694874270d9SVaclav Hapla   Input Parameter:
695874270d9SVaclav Hapla . viewer - the PetscViewer
696874270d9SVaclav Hapla 
697874270d9SVaclav Hapla   Output Parameter:
698874270d9SVaclav Hapla + fileId - The HDF5 file ID
699874270d9SVaclav Hapla - groupId - The HDF5 group ID
700874270d9SVaclav Hapla 
701e74616beSVaclav Hapla   Notes:
702e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
703e74616beSVaclav Hapla 
704874270d9SVaclav Hapla   Level: intermediate
705874270d9SVaclav Hapla 
706874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
707874270d9SVaclav Hapla @*/
70854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
70954dbf706SVaclav Hapla {
71090e03537SVaclav Hapla   hid_t          file_id;
71190e03537SVaclav Hapla   H5O_type_t     type;
71254dbf706SVaclav Hapla   const char     *groupName = NULL;
713e74616beSVaclav Hapla   PetscBool      create;
71454dbf706SVaclav Hapla   PetscErrorCode ierr;
71554dbf706SVaclav Hapla 
71654dbf706SVaclav Hapla   PetscFunctionBegin;
717e74616beSVaclav Hapla   ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr);
71854dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr);
71954dbf706SVaclav Hapla   ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr);
720e74616beSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr);
72190e03537SVaclav 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);
72290e03537SVaclav Hapla   PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
72354dbf706SVaclav Hapla   *fileId  = file_id;
72454dbf706SVaclav Hapla   PetscFunctionReturn(0);
72554dbf706SVaclav Hapla }
72654dbf706SVaclav Hapla 
7273ef9c667SSatish Balay /*@
7285c6c1daeSBarry Smith   PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time.
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith   Not collective
7315c6c1daeSBarry Smith 
7325c6c1daeSBarry Smith   Input Parameter:
7335c6c1daeSBarry Smith . viewer - the PetscViewer
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   Level: intermediate
7365c6c1daeSBarry Smith 
7375c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep()
7385c6c1daeSBarry Smith @*/
7395c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
7405c6c1daeSBarry Smith {
7415c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7425c6c1daeSBarry Smith 
7435c6c1daeSBarry Smith   PetscFunctionBegin;
7445c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7455c6c1daeSBarry Smith   ++hdf5->timestep;
7465c6c1daeSBarry Smith   PetscFunctionReturn(0);
7475c6c1daeSBarry Smith }
7485c6c1daeSBarry Smith 
7493ef9c667SSatish Balay /*@
7505c6c1daeSBarry Smith   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep
7515c6c1daeSBarry Smith   of -1 disables blocking with timesteps.
7525c6c1daeSBarry Smith 
7535c6c1daeSBarry Smith   Not collective
7545c6c1daeSBarry Smith 
7555c6c1daeSBarry Smith   Input Parameters:
7565c6c1daeSBarry Smith + viewer - the PetscViewer
7575c6c1daeSBarry Smith - timestep - The timestep number
7585c6c1daeSBarry Smith 
7595c6c1daeSBarry Smith   Level: intermediate
7605c6c1daeSBarry Smith 
7615c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep()
7625c6c1daeSBarry Smith @*/
7635c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
7645c6c1daeSBarry Smith {
7655c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7665c6c1daeSBarry Smith 
7675c6c1daeSBarry Smith   PetscFunctionBegin;
7685c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7695c6c1daeSBarry Smith   hdf5->timestep = timestep;
7705c6c1daeSBarry Smith   PetscFunctionReturn(0);
7715c6c1daeSBarry Smith }
7725c6c1daeSBarry Smith 
7733ef9c667SSatish Balay /*@
7745c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
7755c6c1daeSBarry Smith 
7765c6c1daeSBarry Smith   Not collective
7775c6c1daeSBarry Smith 
7785c6c1daeSBarry Smith   Input Parameter:
7795c6c1daeSBarry Smith . viewer - the PetscViewer
7805c6c1daeSBarry Smith 
7815c6c1daeSBarry Smith   Output Parameter:
7825c6c1daeSBarry Smith . timestep - The timestep number
7835c6c1daeSBarry Smith 
7845c6c1daeSBarry Smith   Level: intermediate
7855c6c1daeSBarry Smith 
7865c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep()
7875c6c1daeSBarry Smith @*/
7885c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
7895c6c1daeSBarry Smith {
7905c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7915c6c1daeSBarry Smith 
7925c6c1daeSBarry Smith   PetscFunctionBegin;
7935c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7945c6c1daeSBarry Smith   PetscValidPointer(timestep,2);
7955c6c1daeSBarry Smith   *timestep = hdf5->timestep;
7965c6c1daeSBarry Smith   PetscFunctionReturn(0);
7975c6c1daeSBarry Smith }
7985c6c1daeSBarry Smith 
79936bce6e8SMatthew G. Knepley /*@C
80036bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
80136bce6e8SMatthew G. Knepley 
80236bce6e8SMatthew G. Knepley   Not collective
80336bce6e8SMatthew G. Knepley 
80436bce6e8SMatthew G. Knepley   Input Parameter:
80536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
80636bce6e8SMatthew G. Knepley 
80736bce6e8SMatthew G. Knepley   Output Parameter:
80836bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_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 PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
81536bce6e8SMatthew G. Knepley {
81636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
81736bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
81836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
81936bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
82036bce6e8SMatthew G. Knepley #else
82136bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
82236bce6e8SMatthew G. Knepley #endif
82336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
82436bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
82536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
826c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
827de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
82836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
82936bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
83036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
8317e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
83236bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
83336bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
83436bce6e8SMatthew G. Knepley }
83536bce6e8SMatthew G. Knepley 
83636bce6e8SMatthew G. Knepley /*@C
83736bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
83836bce6e8SMatthew G. Knepley 
83936bce6e8SMatthew G. Knepley   Not collective
84036bce6e8SMatthew G. Knepley 
84136bce6e8SMatthew G. Knepley   Input Parameter:
84236bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
84336bce6e8SMatthew G. Knepley 
84436bce6e8SMatthew G. Knepley   Output Parameter:
84536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
84636bce6e8SMatthew G. Knepley 
84736bce6e8SMatthew G. Knepley   Level: advanced
84836bce6e8SMatthew G. Knepley 
84936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType()
85036bce6e8SMatthew G. Knepley @*/
85136bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
85236bce6e8SMatthew G. Knepley {
85336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
85436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
85536bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
85636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
85736bce6e8SMatthew G. Knepley #else
85836bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
85936bce6e8SMatthew G. Knepley #endif
86036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
86136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
86236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
86336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
86436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
86536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
8667e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
86736bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
86836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
86936bce6e8SMatthew G. Knepley }
87036bce6e8SMatthew G. Knepley 
871df863907SAlex Fikl /*@C
872b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
87336bce6e8SMatthew G. Knepley 
87436bce6e8SMatthew G. Knepley   Input Parameters:
87536bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
8764302210dSVaclav Hapla . parent - The parent dataset/group name
87736bce6e8SMatthew G. Knepley . name   - The attribute name
87836bce6e8SMatthew G. Knepley . datatype - The attribute type
87936bce6e8SMatthew G. Knepley - value    - The attribute value
88036bce6e8SMatthew G. Knepley 
88136bce6e8SMatthew G. Knepley   Level: advanced
88236bce6e8SMatthew G. Knepley 
8834302210dSVaclav Hapla   Notes:
8844302210dSVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
8854302210dSVaclav Hapla 
886577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
88736bce6e8SMatthew G. Knepley @*/
8884302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
88936bce6e8SMatthew G. Knepley {
8904302210dSVaclav Hapla   char           *parentAbsPath;
89160568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
892080f144cSVaclav Hapla   PetscBool      has;
89336bce6e8SMatthew G. Knepley   PetscErrorCode ierr;
89436bce6e8SMatthew G. Knepley 
89536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
8965cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
8974302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
898c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
8994302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer,datatype,4);
900b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
9014302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
9024302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr);
9034302210dSVaclav Hapla   ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);
90436bce6e8SMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
9057e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
9067e97c476SMatthew G. Knepley     size_t len;
9077e97c476SMatthew G. Knepley     ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr);
908729ab6d9SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
9097e97c476SMatthew G. Knepley   }
91036bce6e8SMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
911729ab6d9SBarry Smith   PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
9124302210dSVaclav Hapla   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
913080f144cSVaclav Hapla   if (has) {
914080f144cSVaclav Hapla     PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
915080f144cSVaclav Hapla   } else {
91660568592SMatthew G. Knepley     PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
917080f144cSVaclav Hapla   }
918729ab6d9SBarry Smith   PetscStackCallHDF5(H5Awrite,(attribute, dtype, value));
919729ab6d9SBarry Smith   if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype));
920729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
92160568592SMatthew G. Knepley   PetscStackCallHDF5(H5Oclose,(obj));
922729ab6d9SBarry Smith   PetscStackCallHDF5(H5Sclose,(dataspace));
9234302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
92436bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
92536bce6e8SMatthew G. Knepley }
92636bce6e8SMatthew G. Knepley 
927df863907SAlex Fikl /*@C
928577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
929577e0e04SVaclav Hapla 
930577e0e04SVaclav Hapla   Input Parameters:
931577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
932577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
933577e0e04SVaclav Hapla . name     - The attribute name
934577e0e04SVaclav Hapla . datatype - The attribute type
935577e0e04SVaclav Hapla - value    - The attribute value
936577e0e04SVaclav Hapla 
937577e0e04SVaclav Hapla   Notes:
9384302210dSVaclav Hapla   This fails if currentGroup/objectName doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
939577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
940577e0e04SVaclav Hapla 
941577e0e04SVaclav Hapla   Level: advanced
942577e0e04SVaclav Hapla 
943577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
944577e0e04SVaclav Hapla @*/
945577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
946577e0e04SVaclav Hapla {
947577e0e04SVaclav Hapla   PetscErrorCode ierr;
948577e0e04SVaclav Hapla 
949577e0e04SVaclav Hapla   PetscFunctionBegin;
950577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
951577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
952577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
953b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
954577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
955577e0e04SVaclav Hapla   ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
956577e0e04SVaclav Hapla   PetscFunctionReturn(0);
957577e0e04SVaclav Hapla }
958577e0e04SVaclav Hapla 
959577e0e04SVaclav Hapla /*@C
960b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
961ad0c61feSMatthew G. Knepley 
962ad0c61feSMatthew G. Knepley   Input Parameters:
963ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
9644302210dSVaclav Hapla . parent - The parent dataset/group name
965ad0c61feSMatthew G. Knepley . name   - The attribute name
966ad0c61feSMatthew G. Knepley - datatype - The attribute type
967ad0c61feSMatthew G. Knepley 
968ad0c61feSMatthew G. Knepley   Output Parameter:
969ad0c61feSMatthew G. Knepley . value    - The attribute value
970ad0c61feSMatthew G. Knepley 
971708d4cb3SBarry Smith   Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed.
972708d4cb3SBarry Smith 
973ad0c61feSMatthew G. Knepley   Level: advanced
974ad0c61feSMatthew G. Knepley 
9754302210dSVaclav Hapla   Notes:
9764302210dSVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
9774302210dSVaclav Hapla 
978577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
979ad0c61feSMatthew G. Knepley @*/
9804302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value)
981ad0c61feSMatthew G. Knepley {
9824302210dSVaclav Hapla   char           *parentAbsPath;
983f0b58479SMatthew G. Knepley   hid_t          h5, obj, attribute, atype, dtype;
984969748fdSVaclav Hapla   PetscBool      has;
985ad0c61feSMatthew G. Knepley   PetscErrorCode ierr;
986ad0c61feSMatthew G. Knepley 
987ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
9885cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
9894302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
990c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
991b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
9924302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
9934302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr);
9944302210dSVaclav Hapla   if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);}
9954302210dSVaclav Hapla   if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parentAbsPath, name);
996ad0c61feSMatthew G. Knepley   ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr);
997ad0c61feSMatthew G. Knepley   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
9984302210dSVaclav Hapla   PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
99960568592SMatthew G. Knepley   PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1000f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1001f0b58479SMatthew G. Knepley     size_t len;
1002e90c5075SMarek Pecha     PetscStackCallHDF5Return(atype,H5Aget_type,(attribute));
100315b861d2SVaclav Hapla     PetscStackCall("H5Tget_size",len = H5Tget_size(atype));
1004708d4cb3SBarry Smith     ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr);
1005708d4cb3SBarry Smith     PetscStackCallHDF5(H5Tset_size,(dtype, len+1));
1006708d4cb3SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1007708d4cb3SBarry Smith   } else {
100870efba33SBarry Smith     PetscStackCallHDF5(H5Aread,(attribute, dtype, value));
1009708d4cb3SBarry Smith   }
1010729ab6d9SBarry Smith   PetscStackCallHDF5(H5Aclose,(attribute));
1011e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1012e90c5075SMarek Pecha   PetscStackCallHDF5(H5Oclose,(obj));
10134302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
1014ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1015ad0c61feSMatthew G. Knepley }
1016ad0c61feSMatthew G. Knepley 
1017577e0e04SVaclav Hapla /*@C
1018577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
1019577e0e04SVaclav Hapla 
1020577e0e04SVaclav Hapla   Input Parameters:
1021577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
1022577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1023577e0e04SVaclav Hapla . name     - The attribute name
1024577e0e04SVaclav Hapla - datatype - The attribute type
1025577e0e04SVaclav Hapla 
1026577e0e04SVaclav Hapla   Output Parameter:
1027577e0e04SVaclav Hapla . value    - The attribute value
1028577e0e04SVaclav Hapla 
1029577e0e04SVaclav Hapla   Notes:
1030577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1031577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1032577e0e04SVaclav Hapla 
1033577e0e04SVaclav Hapla   Level: advanced
1034577e0e04SVaclav Hapla 
1035577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1036577e0e04SVaclav Hapla @*/
1037577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value)
1038577e0e04SVaclav Hapla {
1039577e0e04SVaclav Hapla   PetscErrorCode ierr;
1040577e0e04SVaclav Hapla 
1041577e0e04SVaclav Hapla   PetscFunctionBegin;
1042577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1043577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1044577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1045b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
1046577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1047577e0e04SVaclav Hapla   ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr);
1048577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1049577e0e04SVaclav Hapla }
1050577e0e04SVaclav Hapla 
1051a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1052a75c3fd4SVaclav Hapla {
1053a75c3fd4SVaclav Hapla   htri_t exists;
1054a75c3fd4SVaclav Hapla   hid_t group;
1055a75c3fd4SVaclav Hapla 
1056a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1057a75c3fd4SVaclav Hapla   PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1058a75c3fd4SVaclav Hapla   if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1059a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1060a75c3fd4SVaclav Hapla     PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1061a75c3fd4SVaclav Hapla     PetscStackCallHDF5(H5Gclose,(group));
1062a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1063a75c3fd4SVaclav Hapla   }
1064a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1065a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1066a75c3fd4SVaclav Hapla }
1067a75c3fd4SVaclav Hapla 
1068bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
10695cdeb108SMatthew G. Knepley {
107090e03537SVaclav Hapla   const char     rootGroupName[] = "/";
10715cdeb108SMatthew G. Knepley   hid_t          h5;
1072e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
107315b861d2SVaclav Hapla   PetscInt       i;
107415b861d2SVaclav Hapla   int            n;
107585688ae2SVaclav Hapla   char           **hierarchy;
107685688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
10775cdeb108SMatthew G. Knepley   PetscErrorCode ierr;
10785cdeb108SMatthew G. Knepley 
10795cdeb108SMatthew G. Knepley   PetscFunctionBegin;
10805cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
108190e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
108290e03537SVaclav Hapla   else name = rootGroupName;
1083ccd66a83SVaclav Hapla   if (has) {
108456cc0592SVaclav Hapla     PetscValidIntPointer(has, 3);
10855cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1086ccd66a83SVaclav Hapla   }
1087ccd66a83SVaclav Hapla   if (otype) {
1088ccd66a83SVaclav Hapla     PetscValidIntPointer(otype, 4);
108956cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1090ccd66a83SVaclav Hapla   }
1091c0cd0301SJed Brown   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
109285688ae2SVaclav Hapla 
109385688ae2SVaclav Hapla   /*
109485688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
109585688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
109685688ae2SVaclav Hapla      1) whether it's a valid link
109785688ae2SVaclav Hapla      2) whether this link resolves to an object
109885688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
109985688ae2SVaclav Hapla   */
110085688ae2SVaclav Hapla   ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr);
110185688ae2SVaclav Hapla   if (!n) {
110285688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1103ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1104ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
110585688ae2SVaclav Hapla     ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
110685688ae2SVaclav Hapla     PetscFunctionReturn(0);
110785688ae2SVaclav Hapla   }
110885688ae2SVaclav Hapla   for (i=0; i<n; i++) {
110985688ae2SVaclav Hapla     ierr = PetscStrcat(buf,"/");CHKERRQ(ierr);
111085688ae2SVaclav Hapla     ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr);
1111a75c3fd4SVaclav Hapla     ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr);
1112a75c3fd4SVaclav Hapla     if (!exists) break;
111385688ae2SVaclav Hapla   }
111485688ae2SVaclav Hapla   ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr);
111585688ae2SVaclav Hapla 
111685688ae2SVaclav Hapla   /* If the object exists, get its type */
1117ccd66a83SVaclav Hapla   if (exists && otype) {
11185cdeb108SMatthew G. Knepley     H5O_info_t info;
11195cdeb108SMatthew G. Knepley 
112076276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
112104633f7fSVaclav Hapla     PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
112256cc0592SVaclav Hapla     *otype = info.type;
11235cdeb108SMatthew G. Knepley   }
1124ccd66a83SVaclav Hapla   if (has) *has = exists;
11255cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
11265cdeb108SMatthew G. Knepley }
11275cdeb108SMatthew G. Knepley 
11284302210dSVaclav Hapla /*@C
11290a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
11300a9f38efSVaclav Hapla 
11310a9f38efSVaclav Hapla   Input Parameters:
11324302210dSVaclav Hapla + viewer - The HDF5 viewer
11334302210dSVaclav Hapla - path - The group path
11340a9f38efSVaclav Hapla 
11350a9f38efSVaclav Hapla   Output Parameter:
11360a9f38efSVaclav Hapla . has - Flag for group existence
11370a9f38efSVaclav Hapla 
11380a9f38efSVaclav Hapla   Level: advanced
11390a9f38efSVaclav Hapla 
11404302210dSVaclav Hapla   Notes:
11414302210dSVaclav Hapla   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
11424302210dSVaclav Hapla   If the path exists but is not a group, PETSC_FALSE is returned.
11434302210dSVaclav Hapla 
1144*89e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup()
11450a9f38efSVaclav Hapla @*/
11464302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
11470a9f38efSVaclav Hapla {
11480a9f38efSVaclav Hapla   H5O_type_t     type;
11494302210dSVaclav Hapla   char           *abspath;
11500a9f38efSVaclav Hapla   PetscErrorCode ierr;
11510a9f38efSVaclav Hapla 
11520a9f38efSVaclav Hapla   PetscFunctionBegin;
11530a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11544302210dSVaclav Hapla   if (path) PetscValidCharPointer(path,2);
11554302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
11564302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
11574302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
11584302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
11594302210dSVaclav Hapla   ierr = PetscFree(abspath);CHKERRQ(ierr);
11600a9f38efSVaclav Hapla   PetscFunctionReturn(0);
11610a9f38efSVaclav Hapla }
11620a9f38efSVaclav Hapla 
1163*89e0ef10SVaclav Hapla /*@C
1164*89e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
1165*89e0ef10SVaclav Hapla 
1166*89e0ef10SVaclav Hapla   Input Parameters:
1167*89e0ef10SVaclav Hapla + viewer - The HDF5 viewer
1168*89e0ef10SVaclav Hapla - path - The dataset path
1169*89e0ef10SVaclav Hapla 
1170*89e0ef10SVaclav Hapla   Output Parameter:
1171*89e0ef10SVaclav Hapla . has - Flag whether dataset exists
1172*89e0ef10SVaclav Hapla 
1173*89e0ef10SVaclav Hapla   Level: advanced
1174*89e0ef10SVaclav Hapla 
1175*89e0ef10SVaclav Hapla   Notes:
1176*89e0ef10SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group.
1177*89e0ef10SVaclav Hapla   If path is NULL or empty, has is set to PETSC_FALSE.
1178*89e0ef10SVaclav Hapla   If the path exists but is not a dataset, has is set to PETSC_FALSE as well.
1179*89e0ef10SVaclav Hapla 
1180*89e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1181*89e0ef10SVaclav Hapla @*/
1182*89e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1183*89e0ef10SVaclav Hapla {
1184*89e0ef10SVaclav Hapla   H5O_type_t     type;
1185*89e0ef10SVaclav Hapla   char           *abspath;
1186*89e0ef10SVaclav Hapla   PetscErrorCode ierr;
1187*89e0ef10SVaclav Hapla 
1188*89e0ef10SVaclav Hapla   PetscFunctionBegin;
1189*89e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1190*89e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path,2);
1191*89e0ef10SVaclav Hapla   PetscValidBoolPointer(has,3);
1192*89e0ef10SVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr);
1193*89e0ef10SVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr);
1194*89e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
1195*89e0ef10SVaclav Hapla   ierr = PetscFree(abspath);CHKERRQ(ierr);
1196*89e0ef10SVaclav Hapla   PetscFunctionReturn(0);
1197*89e0ef10SVaclav Hapla }
1198*89e0ef10SVaclav Hapla 
11990a9f38efSVaclav Hapla /*@
1200e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1201ecce1506SVaclav Hapla 
1202ecce1506SVaclav Hapla   Input Parameters:
1203ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1204ecce1506SVaclav Hapla - obj    - The named object
1205ecce1506SVaclav Hapla 
1206ecce1506SVaclav Hapla   Output Parameter:
1207*89e0ef10SVaclav Hapla . has    - Flag for dataset existence
1208ecce1506SVaclav Hapla 
1209e3f143f7SVaclav Hapla   Notes:
1210*89e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
12114302210dSVaclav Hapla   If the path currentGroup/objectName exists but is not a dataset, has is set to PETSC_FALSE as well.
1212e3f143f7SVaclav Hapla 
1213ecce1506SVaclav Hapla   Level: advanced
1214ecce1506SVaclav Hapla 
1215*89e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup()
1216ecce1506SVaclav Hapla @*/
1217ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1218ecce1506SVaclav Hapla {
1219*89e0ef10SVaclav Hapla   size_t         len;
1220ecce1506SVaclav Hapla   PetscErrorCode ierr;
1221ecce1506SVaclav Hapla 
1222ecce1506SVaclav Hapla   PetscFunctionBegin;
1223c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1224c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
12254302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
1226*89e0ef10SVaclav Hapla   ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr);
1227*89e0ef10SVaclav Hapla   if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
1228*89e0ef10SVaclav Hapla   ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr);
1229ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1230ecce1506SVaclav Hapla }
1231ecce1506SVaclav Hapla 
1232df863907SAlex Fikl /*@C
1233ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1234e7bdbf8eSMatthew G. Knepley 
1235e7bdbf8eSMatthew G. Knepley   Input Parameters:
1236e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
12374302210dSVaclav Hapla . parent - The parent dataset/group name
1238e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1239e7bdbf8eSMatthew G. Knepley 
1240e7bdbf8eSMatthew G. Knepley   Output Parameter:
1241e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1242e7bdbf8eSMatthew G. Knepley 
1243e7bdbf8eSMatthew G. Knepley   Level: advanced
1244e7bdbf8eSMatthew G. Knepley 
12454302210dSVaclav Hapla   Notes:
12464302210dSVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else the dataset/group is relative to the current pushed group. NULL means the current pushed group.
12474302210dSVaclav Hapla 
1248577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1249e7bdbf8eSMatthew G. Knepley @*/
12504302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1251e7bdbf8eSMatthew G. Knepley {
12524302210dSVaclav Hapla   char           *parentAbsPath;
1253e7bdbf8eSMatthew G. Knepley   PetscErrorCode ierr;
1254e7bdbf8eSMatthew G. Knepley 
1255e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
12565cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
12574302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent,2);
1258c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
12594302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
12604302210dSVaclav Hapla   ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr);
12614302210dSVaclav Hapla   ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr);
12624302210dSVaclav Hapla   if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);}
12634302210dSVaclav Hapla   ierr = PetscFree(parentAbsPath);CHKERRQ(ierr);
126406db490cSVaclav Hapla   PetscFunctionReturn(0);
126506db490cSVaclav Hapla }
126606db490cSVaclav Hapla 
1267577e0e04SVaclav Hapla /*@C
1268577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1269577e0e04SVaclav Hapla 
1270577e0e04SVaclav Hapla   Input Parameters:
1271577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1272577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1273577e0e04SVaclav Hapla - name   - The attribute name
1274577e0e04SVaclav Hapla 
1275577e0e04SVaclav Hapla   Output Parameter:
1276577e0e04SVaclav Hapla . has    - Flag for attribute existence
1277577e0e04SVaclav Hapla 
1278577e0e04SVaclav Hapla   Notes:
1279577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1280577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1281577e0e04SVaclav Hapla 
1282577e0e04SVaclav Hapla   Level: advanced
1283577e0e04SVaclav Hapla 
1284577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup()
1285577e0e04SVaclav Hapla @*/
1286577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1287577e0e04SVaclav Hapla {
1288577e0e04SVaclav Hapla   PetscErrorCode ierr;
1289577e0e04SVaclav Hapla 
1290577e0e04SVaclav Hapla   PetscFunctionBegin;
1291577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1292577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1293577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
12944302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
1295577e0e04SVaclav Hapla   ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr);
1296577e0e04SVaclav Hapla   ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr);
1297577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1298577e0e04SVaclav Hapla }
1299577e0e04SVaclav Hapla 
130006db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
130106db490cSVaclav Hapla {
130206db490cSVaclav Hapla   hid_t          h5;
130306db490cSVaclav Hapla   htri_t         hhas;
130406db490cSVaclav Hapla   PetscErrorCode ierr;
130506db490cSVaclav Hapla 
130606db490cSVaclav Hapla   PetscFunctionBegin;
130706db490cSVaclav Hapla   ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr);
13082f936e54SVaclav Hapla   PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
13095cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1310e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1311e7bdbf8eSMatthew G. Knepley }
1312e7bdbf8eSMatthew G. Knepley 
1313a75e6a4aSMatthew G. Knepley /*
1314a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1315a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1316a75e6a4aSMatthew G. Knepley */
1317d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1318a75e6a4aSMatthew G. Knepley 
1319a75e6a4aSMatthew G. Knepley /*@C
1320a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1321a75e6a4aSMatthew G. Knepley 
1322d083f849SBarry Smith   Collective
1323a75e6a4aSMatthew G. Knepley 
1324a75e6a4aSMatthew G. Knepley   Input Parameter:
1325a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1326a75e6a4aSMatthew G. Knepley 
1327a75e6a4aSMatthew G. Knepley   Level: intermediate
1328a75e6a4aSMatthew G. Knepley 
1329a75e6a4aSMatthew G. Knepley   Options Database Keys:
1330a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name>
1331a75e6a4aSMatthew G. Knepley 
1332a75e6a4aSMatthew G. Knepley   Environmental variables:
1333a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME
1334a75e6a4aSMatthew G. Knepley 
1335a75e6a4aSMatthew G. Knepley   Notes:
1336a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1337a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1338a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1339a75e6a4aSMatthew G. Knepley 
1340a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy()
1341a75e6a4aSMatthew G. Knepley @*/
1342a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1343a75e6a4aSMatthew G. Knepley {
1344a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1345a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1346a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1347a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1348a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1349a75e6a4aSMatthew G. Knepley 
1350a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1351908793a3SLisandro 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);}
1352a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1353908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1354908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1355a75e6a4aSMatthew G. Knepley   }
135647435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1357908793a3SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1358a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1359a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
13602cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1361a75e6a4aSMatthew G. Knepley     if (!flg) {
1362a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
13632cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1364a75e6a4aSMatthew G. Knepley     }
1365a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
13662cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1367a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
13682cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
136947435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1370908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1371a75e6a4aSMatthew G. Knepley   }
1372a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
13732cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1374a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1375a75e6a4aSMatthew G. Knepley }
1376