xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 792fecdfe9134cce4d631112660ddd34f063bc17)
16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/
25c6c1daeSBarry Smith 
3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*);
406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*);
506db490cSVaclav Hapla 
64302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[])
76c132bc1SVaclav Hapla {
84302210dSVaclav Hapla   PetscBool      relative = PETSC_FALSE;
96c132bc1SVaclav Hapla   const char     *group;
106c132bc1SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN] = "";
116c132bc1SVaclav Hapla 
126c132bc1SVaclav Hapla   PetscFunctionBegin;
139566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetGroup(viewer, &group));
149566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative));
154302210dSVaclav Hapla   if (relative) {
169566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(buf, group));
179566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
189566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, path));
199566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(buf, abspath));
204302210dSVaclav Hapla   } else {
219566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(path, abspath));
224302210dSVaclav Hapla   }
236c132bc1SVaclav Hapla   PetscFunctionReturn(0);
246c132bc1SVaclav Hapla }
256c132bc1SVaclav Hapla 
26577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
27577e0e04SVaclav Hapla {
28577e0e04SVaclav Hapla   PetscBool      has;
29577e0e04SVaclav Hapla 
30577e0e04SVaclav Hapla   PetscFunctionBegin;
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
32577e0e04SVaclav Hapla   if (!has) {
3389e0ef10SVaclav Hapla     const char *group;
349566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetGroup(viewer, &group));
3598921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
36577e0e04SVaclav Hapla   }
37577e0e04SVaclav Hapla   PetscFunctionReturn(0);
38577e0e04SVaclav Hapla }
39577e0e04SVaclav Hapla 
404416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v)
4182ea9b62SBarry Smith {
42ee8b9732SVaclav Hapla   PetscBool        flg = PETSC_FALSE, set;
4382ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data;
4482ea9b62SBarry Smith 
4582ea9b62SBarry Smith   PetscFunctionBegin;
46d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject,"HDF5 PetscViewer Options");
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL));
489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL));
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set));
509566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetCollective(v,flg));
5119a20e4cSMatthew G. Knepley   flg  = PETSC_FALSE;
529566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping","Set default timestepping state","PetscViewerHDF5SetDefaultTimestepping",flg,&flg,&set));
539566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v,flg));
54d0609cedSBarry Smith   PetscOptionsHeadEnd();
5582ea9b62SBarry Smith   PetscFunctionReturn(0);
5682ea9b62SBarry Smith }
5782ea9b62SBarry Smith 
581b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)
591b793a25SVaclav Hapla {
601b793a25SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*)v->data;
6103fa8834SVaclav Hapla   PetscBool         flg;
621b793a25SVaclav Hapla 
631b793a25SVaclav Hapla   PetscFunctionBegin;
641b793a25SVaclav Hapla   if (hdf5->filename) {
659566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename));
661b793a25SVaclav Hapla   }
679566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]));
689566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]));
699566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetCollective(v,&flg));
709566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent"));
719566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer,"Default timestepping: %s\n",PetscBools[hdf5->defTimestepping]));
721b793a25SVaclav Hapla   PetscFunctionReturn(0);
731b793a25SVaclav Hapla }
741b793a25SVaclav Hapla 
755c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
765c6c1daeSBarry Smith {
775c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
785c6c1daeSBarry Smith 
795c6c1daeSBarry Smith   PetscFunctionBegin;
809566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5->filename));
81*792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose,(hdf5->file_id));
825c6c1daeSBarry Smith   PetscFunctionReturn(0);
835c6c1daeSBarry Smith }
845c6c1daeSBarry Smith 
856226335aSVaclav Hapla static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
866226335aSVaclav Hapla {
876226335aSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data;
886226335aSVaclav Hapla 
896226335aSVaclav Hapla   PetscFunctionBegin;
90*792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fflush,(hdf5->file_id, H5F_SCOPE_LOCAL));
916226335aSVaclav Hapla   PetscFunctionReturn(0);
926226335aSVaclav Hapla }
936226335aSVaclav Hapla 
949b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
955c6c1daeSBarry Smith {
965c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
975c6c1daeSBarry Smith 
985c6c1daeSBarry Smith   PetscFunctionBegin;
99*792fecdfSBarry Smith   PetscCallHDF5(H5Pclose,(hdf5->dxpl_id));
1009566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_HDF5(viewer));
1015c6c1daeSBarry Smith   while (hdf5->groups) {
1027d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
1035c6c1daeSBarry Smith 
1049566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups->name));
1059566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups));
1065c6c1daeSBarry Smith     hdf5->groups = tmp;
1075c6c1daeSBarry Smith   }
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL));
1109566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL));
1119566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL));
1122e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL));
1139566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL));
1149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL));
1152e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetCollective_C",NULL));
1162e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetCollective_C",NULL));
1179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetDefaultTimestepping_C",NULL));
1189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetDefaultTimestepping_C",NULL));
1195c6c1daeSBarry Smith   PetscFunctionReturn(0);
1205c6c1daeSBarry Smith }
1215c6c1daeSBarry Smith 
1229b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
1235c6c1daeSBarry Smith {
1245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1255c6c1daeSBarry Smith 
1265c6c1daeSBarry Smith   PetscFunctionBegin;
1275c6c1daeSBarry Smith   hdf5->btype = type;
1285c6c1daeSBarry Smith   PetscFunctionReturn(0);
1295c6c1daeSBarry Smith }
1305c6c1daeSBarry Smith 
1310b62783dSVaclav Hapla static PetscErrorCode  PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
1320b62783dSVaclav Hapla {
1330b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
1340b62783dSVaclav Hapla 
1350b62783dSVaclav Hapla   PetscFunctionBegin;
1360b62783dSVaclav Hapla   *type = hdf5->btype;
1370b62783dSVaclav Hapla   PetscFunctionReturn(0);
1380b62783dSVaclav Hapla }
1390b62783dSVaclav Hapla 
1409b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
14182ea9b62SBarry Smith {
14282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
14382ea9b62SBarry Smith 
14482ea9b62SBarry Smith   PetscFunctionBegin;
14582ea9b62SBarry Smith   hdf5->basedimension2 = flg;
14682ea9b62SBarry Smith   PetscFunctionReturn(0);
14782ea9b62SBarry Smith }
14882ea9b62SBarry Smith 
149df863907SAlex Fikl /*@
15082ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
15182ea9b62SBarry Smith        dimension of 2.
15282ea9b62SBarry Smith 
15382ea9b62SBarry Smith     Logically Collective on PetscViewer
15482ea9b62SBarry Smith 
15582ea9b62SBarry Smith   Input Parameters:
15682ea9b62SBarry Smith +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
15782ea9b62SBarry 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
15882ea9b62SBarry Smith 
15982ea9b62SBarry Smith   Options Database:
16082ea9b62SBarry 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
16182ea9b62SBarry Smith 
16295452b02SPatrick Sanan   Notes:
16395452b02SPatrick 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
16482ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
16582ea9b62SBarry Smith 
16682ea9b62SBarry Smith   Level: intermediate
16782ea9b62SBarry Smith 
168db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
16982ea9b62SBarry Smith 
17082ea9b62SBarry Smith @*/
17182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)
17282ea9b62SBarry Smith {
17382ea9b62SBarry Smith   PetscFunctionBegin;
17482ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
175cac4c232SBarry Smith   PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));
17682ea9b62SBarry Smith   PetscFunctionReturn(0);
17782ea9b62SBarry Smith }
17882ea9b62SBarry Smith 
179df863907SAlex Fikl /*@
18082ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
18182ea9b62SBarry Smith        dimension of 2.
18282ea9b62SBarry Smith 
18382ea9b62SBarry Smith     Logically Collective on PetscViewer
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith   Input Parameter:
18682ea9b62SBarry Smith .  viewer - the PetscViewer, must be of type HDF5
18782ea9b62SBarry Smith 
18882ea9b62SBarry Smith   Output Parameter:
18982ea9b62SBarry 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
19082ea9b62SBarry Smith 
19195452b02SPatrick Sanan   Notes:
19295452b02SPatrick 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
19382ea9b62SBarry Smith          of one when the dimension is lower. Others think the option is crazy.
19482ea9b62SBarry Smith 
19582ea9b62SBarry Smith   Level: intermediate
19682ea9b62SBarry Smith 
197db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
19882ea9b62SBarry Smith 
19982ea9b62SBarry Smith @*/
20082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg)
20182ea9b62SBarry Smith {
20282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
20382ea9b62SBarry Smith 
20482ea9b62SBarry Smith   PetscFunctionBegin;
20582ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
20682ea9b62SBarry Smith   *flg = hdf5->basedimension2;
20782ea9b62SBarry Smith   PetscFunctionReturn(0);
20882ea9b62SBarry Smith }
20982ea9b62SBarry Smith 
2109b1403beSVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
2119a0502c6SHåkon Strandenes {
2129a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2139a0502c6SHåkon Strandenes 
2149a0502c6SHåkon Strandenes   PetscFunctionBegin;
2159a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2169a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2179a0502c6SHåkon Strandenes }
2189a0502c6SHåkon Strandenes 
219df863907SAlex Fikl /*@
2209a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
2219a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2229a0502c6SHåkon Strandenes 
2239a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2249a0502c6SHåkon Strandenes 
2259a0502c6SHåkon Strandenes   Input Parameters:
2269a0502c6SHåkon Strandenes +  viewer - the PetscViewer; if it is not hdf5 then this command is ignored
2279a0502c6SHåkon Strandenes -  flg - if PETSC_TRUE the data will be written to disk with single precision
2289a0502c6SHåkon Strandenes 
2299a0502c6SHåkon Strandenes   Options Database:
2309a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2319a0502c6SHåkon Strandenes 
23295452b02SPatrick Sanan   Notes:
23395452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2349a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2359a0502c6SHåkon Strandenes 
2369a0502c6SHåkon Strandenes   Level: intermediate
2379a0502c6SHåkon Strandenes 
238db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
239db781477SPatrick Sanan           `PetscReal`
2409a0502c6SHåkon Strandenes 
2419a0502c6SHåkon Strandenes @*/
2429a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)
2439a0502c6SHåkon Strandenes {
2449a0502c6SHåkon Strandenes   PetscFunctionBegin;
2459a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
246cac4c232SBarry Smith   PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));
2479a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2489a0502c6SHåkon Strandenes }
2499a0502c6SHåkon Strandenes 
250df863907SAlex Fikl /*@
2519a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
2529a0502c6SHåkon Strandenes        compiled with double precision PetscReal.
2539a0502c6SHåkon Strandenes 
2549a0502c6SHåkon Strandenes     Logically Collective on PetscViewer
2559a0502c6SHåkon Strandenes 
2569a0502c6SHåkon Strandenes   Input Parameter:
2579a0502c6SHåkon Strandenes .  viewer - the PetscViewer, must be of type HDF5
2589a0502c6SHåkon Strandenes 
2599a0502c6SHåkon Strandenes   Output Parameter:
2609a0502c6SHåkon Strandenes .  flg - if PETSC_TRUE the data will be written to disk with single precision
2619a0502c6SHåkon Strandenes 
26295452b02SPatrick Sanan   Notes:
26395452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2649a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2659a0502c6SHåkon Strandenes 
2669a0502c6SHåkon Strandenes   Level: intermediate
2679a0502c6SHåkon Strandenes 
268db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
269db781477SPatrick Sanan           `PetscReal`
2709a0502c6SHåkon Strandenes 
2719a0502c6SHåkon Strandenes @*/
2729a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg)
2739a0502c6SHåkon Strandenes {
2749a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
2759a0502c6SHåkon Strandenes 
2769a0502c6SHåkon Strandenes   PetscFunctionBegin;
2779a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
2789a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2799a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2809a0502c6SHåkon Strandenes }
2819a0502c6SHåkon Strandenes 
282ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
283ee8b9732SVaclav Hapla {
284ee8b9732SVaclav Hapla   PetscFunctionBegin;
285ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
286ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
287c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL)
288ee8b9732SVaclav Hapla   {
289ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
290*792fecdfSBarry Smith     PetscCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
291ee8b9732SVaclav Hapla   }
292ee8b9732SVaclav Hapla #else
2939566063dSJacob Faibussowitsch   if (flg) PetscCall(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"));
294ee8b9732SVaclav Hapla #endif
295ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
296ee8b9732SVaclav Hapla }
297ee8b9732SVaclav Hapla 
298ee8b9732SVaclav Hapla /*@
299ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
300ee8b9732SVaclav Hapla 
301ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
302ee8b9732SVaclav Hapla 
303ee8b9732SVaclav Hapla   Input Parameters:
304ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
305ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default)
306ee8b9732SVaclav Hapla 
307ee8b9732SVaclav Hapla   Options Database:
308ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
309ee8b9732SVaclav Hapla 
310ee8b9732SVaclav Hapla   Notes:
311ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
31253effed7SBarry 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.
313ee8b9732SVaclav Hapla 
314ee8b9732SVaclav Hapla   Developer notes:
315ee8b9732SVaclav Hapla   In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively.
316ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
317ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
318ee8b9732SVaclav Hapla 
319ee8b9732SVaclav Hapla   Level: intermediate
320ee8b9732SVaclav Hapla 
321db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
322ee8b9732SVaclav Hapla 
323ee8b9732SVaclav Hapla @*/
324ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)
325ee8b9732SVaclav Hapla {
326ee8b9732SVaclav Hapla   PetscFunctionBegin;
327ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
328ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer,flg,2);
329cac4c232SBarry Smith   PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));
330ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
331ee8b9732SVaclav Hapla }
332ee8b9732SVaclav Hapla 
333ee8b9732SVaclav Hapla static PetscErrorCode  PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
334ee8b9732SVaclav Hapla {
335c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
336ee8b9732SVaclav Hapla   PetscViewer_HDF5  *hdf5 = (PetscViewer_HDF5*) viewer->data;
337ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
338c796909eSBarry Smith #endif
339ee8b9732SVaclav Hapla 
340ee8b9732SVaclav Hapla   PetscFunctionBegin;
341c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
342c796909eSBarry Smith   *flg = PETSC_FALSE;
343c796909eSBarry Smith #else
344*792fecdfSBarry Smith   PetscCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode));
345ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
346c796909eSBarry Smith #endif
347ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
348ee8b9732SVaclav Hapla }
349ee8b9732SVaclav Hapla 
350ee8b9732SVaclav Hapla /*@
351ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
352ee8b9732SVaclav Hapla 
353ee8b9732SVaclav Hapla   Not Collective
354ee8b9732SVaclav Hapla 
355ee8b9732SVaclav Hapla   Input Parameters:
356ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer
357ee8b9732SVaclav Hapla 
358ee8b9732SVaclav Hapla   Output Parameters:
359ee8b9732SVaclav Hapla . flg - the flag
360ee8b9732SVaclav Hapla 
361ee8b9732SVaclav Hapla   Level: intermediate
362ee8b9732SVaclav Hapla 
363ee8b9732SVaclav Hapla   Notes:
364c796909eSBarry 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.
365ee8b9732SVaclav Hapla   For more details, see PetscViewerHDF5SetCollective().
366ee8b9732SVaclav Hapla 
367db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
368ee8b9732SVaclav Hapla 
369ee8b9732SVaclav Hapla @*/
370ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg)
371ee8b9732SVaclav Hapla {
372ee8b9732SVaclav Hapla   PetscFunctionBegin;
373ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
374534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
375ee8b9732SVaclav Hapla 
376cac4c232SBarry Smith   PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));
377ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
378ee8b9732SVaclav Hapla }
379ee8b9732SVaclav Hapla 
3809b1403beSVaclav Hapla static PetscErrorCode  PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
3815c6c1daeSBarry Smith {
3825c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
3835c6c1daeSBarry Smith   hid_t             plist_id;
3845c6c1daeSBarry Smith 
3855c6c1daeSBarry Smith   PetscFunctionBegin;
386*792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose,(hdf5->file_id));
3879566063dSJacob Faibussowitsch   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
3889566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &hdf5->filename));
3895c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
390*792fecdfSBarry Smith   PetscCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS));
391c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
392*792fecdfSBarry Smith   PetscCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
393c796909eSBarry Smith #endif
3945c6c1daeSBarry Smith   /* Create or open the file collectively */
3955c6c1daeSBarry Smith   switch (hdf5->btype) {
3965c6c1daeSBarry Smith   case FILE_MODE_READ:
3978a2871f6SBarry Smith     if (PetscDefined(USE_DEBUG)) {
3988a2871f6SBarry Smith       PetscMPIInt rank;
3998a2871f6SBarry Smith       PetscBool   flg;
4008a2871f6SBarry Smith 
4019566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank));
4028a2871f6SBarry Smith       if (rank == 0) {
4039566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
4045f80ce2aSJacob Faibussowitsch         PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"File %s requested for reading does not exist",hdf5->filename);
4058a2871f6SBarry Smith       }
4069566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
4078a2871f6SBarry Smith     }
408*792fecdfSBarry Smith     PetscCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id));
4095c6c1daeSBarry Smith     break;
4105c6c1daeSBarry Smith   case FILE_MODE_APPEND:
4117e4fd573SVaclav Hapla   case FILE_MODE_UPDATE:
4127e4fd573SVaclav Hapla   {
4137e4fd573SVaclav Hapla     PetscBool flg;
4149566063dSJacob Faibussowitsch     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
415*792fecdfSBarry Smith     if (flg) PetscCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id));
416*792fecdfSBarry Smith     else     PetscCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4175c6c1daeSBarry Smith     break;
4187e4fd573SVaclav Hapla   }
4195c6c1daeSBarry Smith   case FILE_MODE_WRITE:
420*792fecdfSBarry Smith     PetscCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
4215c6c1daeSBarry Smith     break;
4227e4fd573SVaclav Hapla   case FILE_MODE_UNDEFINED:
4237e4fd573SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
4245c6c1daeSBarry Smith   default:
42598921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]);
4265c6c1daeSBarry Smith   }
4275f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->file_id >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
428*792fecdfSBarry Smith   PetscCallHDF5(H5Pclose,(plist_id));
4295c6c1daeSBarry Smith   PetscFunctionReturn(0);
4305c6c1daeSBarry Smith }
4315c6c1daeSBarry Smith 
432d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name)
433d1232d7fSToby Isaac {
434d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data;
435d1232d7fSToby Isaac 
436d1232d7fSToby Isaac   PetscFunctionBegin;
437d1232d7fSToby Isaac   *name = vhdf5->filename;
438d1232d7fSToby Isaac   PetscFunctionReturn(0);
439d1232d7fSToby Isaac }
440d1232d7fSToby Isaac 
441b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
442b723ab35SVaclav Hapla {
4435dc64a97SVaclav Hapla   /*
444b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4455dc64a97SVaclav Hapla   */
446b723ab35SVaclav Hapla 
447b723ab35SVaclav Hapla   PetscFunctionBegin;
448b723ab35SVaclav Hapla   PetscFunctionReturn(0);
449b723ab35SVaclav Hapla }
450b723ab35SVaclav Hapla 
45119a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
45219a20e4cSMatthew G. Knepley {
45319a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
45419a20e4cSMatthew G. Knepley 
45519a20e4cSMatthew G. Knepley   PetscFunctionBegin;
45619a20e4cSMatthew G. Knepley   hdf5->defTimestepping = flg;
45719a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
45819a20e4cSMatthew G. Knepley }
45919a20e4cSMatthew G. Knepley 
46019a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
46119a20e4cSMatthew G. Knepley {
46219a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
46319a20e4cSMatthew G. Knepley 
46419a20e4cSMatthew G. Knepley   PetscFunctionBegin;
46519a20e4cSMatthew G. Knepley   *flg = hdf5->defTimestepping;
46619a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
46719a20e4cSMatthew G. Knepley }
46819a20e4cSMatthew G. Knepley 
46919a20e4cSMatthew G. Knepley /*@
47019a20e4cSMatthew G. Knepley   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
47119a20e4cSMatthew G. Knepley 
47219a20e4cSMatthew G. Knepley   Logically Collective on PetscViewer
47319a20e4cSMatthew G. Knepley 
47419a20e4cSMatthew G. Knepley   Input Parameters:
47519a20e4cSMatthew G. Knepley + viewer - the PetscViewer; if it is not hdf5 then this command is ignored
47619a20e4cSMatthew G. Knepley - flg    - if PETSC_TRUE we will assume that timestepping is on
47719a20e4cSMatthew G. Knepley 
47819a20e4cSMatthew G. Knepley   Options Database:
47919a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
48019a20e4cSMatthew G. Knepley 
48119a20e4cSMatthew G. Knepley   Notes:
48219a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
48319a20e4cSMatthew G. Knepley 
48419a20e4cSMatthew G. Knepley   Level: intermediate
48519a20e4cSMatthew G. Knepley 
486db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
48719a20e4cSMatthew G. Knepley @*/
48819a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
48919a20e4cSMatthew G. Knepley {
49019a20e4cSMatthew G. Knepley   PetscFunctionBegin;
49119a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
492cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg));
49319a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
49419a20e4cSMatthew G. Knepley }
49519a20e4cSMatthew G. Knepley 
49619a20e4cSMatthew G. Knepley /*@
49719a20e4cSMatthew G. Knepley   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
49819a20e4cSMatthew G. Knepley 
49919a20e4cSMatthew G. Knepley   Not collective
50019a20e4cSMatthew G. Knepley 
50119a20e4cSMatthew G. Knepley   Input Parameter:
50219a20e4cSMatthew G. Knepley . viewer - the PetscViewer
50319a20e4cSMatthew G. Knepley 
50419a20e4cSMatthew G. Knepley   Output Parameter:
50519a20e4cSMatthew G. Knepley . flg    - if PETSC_TRUE we will assume that timestepping is on
50619a20e4cSMatthew G. Knepley 
50719a20e4cSMatthew G. Knepley   Notes:
50819a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
50919a20e4cSMatthew G. Knepley 
51019a20e4cSMatthew G. Knepley   Level: intermediate
51119a20e4cSMatthew G. Knepley 
512db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
51319a20e4cSMatthew G. Knepley @*/
51419a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
51519a20e4cSMatthew G. Knepley {
51619a20e4cSMatthew G. Knepley   PetscFunctionBegin;
51719a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
518cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg));
51919a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
52019a20e4cSMatthew G. Knepley }
52119a20e4cSMatthew G. Knepley 
5228556b5ebSBarry Smith /*MC
5238556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
5248556b5ebSBarry Smith 
525db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
526db781477SPatrick Sanan           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
527db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
528db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
5298556b5ebSBarry Smith 
5301b266c99SBarry Smith   Level: beginner
5318556b5ebSBarry Smith M*/
532d1232d7fSToby Isaac 
5338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
5345c6c1daeSBarry Smith {
5355c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith   PetscFunctionBegin;
53899335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
53999335c34SVaclav Hapla   {
54099335c34SVaclav Hapla     PetscMPIInt size;
5419566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
5425f80ce2aSJacob Faibussowitsch     PetscCheck(size <= 1,PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)");
54399335c34SVaclav Hapla   }
54499335c34SVaclav Hapla #endif
54599335c34SVaclav Hapla 
5469566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(v,&hdf5));
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith   v->data                = (void*) hdf5;
5495c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
55082ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
551b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
5521b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
5536226335aSVaclav Hapla   v->ops->flush          = PetscViewerFlush_HDF5;
5547e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
555908793a3SLisandro Dalcin   hdf5->filename         = NULL;
5565c6c1daeSBarry Smith   hdf5->timestep         = -1;
5570298fd71SBarry Smith   hdf5->groups           = NULL;
5585c6c1daeSBarry Smith 
559*792fecdfSBarry Smith   PetscCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER));
5609c5072fbSVaclav Hapla 
5619566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5));
5629566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5));
5639566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5));
5649566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5));
5659566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5));
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5));
5679566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5));
5689566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5));
5699566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5));
5709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5));
5715c6c1daeSBarry Smith   PetscFunctionReturn(0);
5725c6c1daeSBarry Smith }
5735c6c1daeSBarry Smith 
5745c6c1daeSBarry Smith /*@C
5755c6c1daeSBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output.
5765c6c1daeSBarry Smith 
577d083f849SBarry Smith    Collective
5785c6c1daeSBarry Smith 
5795c6c1daeSBarry Smith    Input Parameters:
5805c6c1daeSBarry Smith +  comm - MPI communicator
5815c6c1daeSBarry Smith .  name - name of file
5825c6c1daeSBarry Smith -  type - type of file
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith    Output Parameter:
5855c6c1daeSBarry Smith .  hdf5v - PetscViewer for HDF5 input/output to use with the specified file
5865c6c1daeSBarry Smith 
58782ea9b62SBarry Smith   Options Database:
588a2b725a8SWilliam 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
589a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
59082ea9b62SBarry Smith 
5915c6c1daeSBarry Smith    Level: beginner
5925c6c1daeSBarry Smith 
5937e4fd573SVaclav Hapla    Notes:
5947e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
5957e4fd573SVaclav Hapla +  FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
5967e4fd573SVaclav Hapla .  FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
5977e4fd573SVaclav 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]
5987e4fd573SVaclav Hapla -  FILE_MODE_UPDATE - same as FILE_MODE_APPEND
5997e4fd573SVaclav Hapla 
6007e4fd573SVaclav 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.
6017e4fd573SVaclav Hapla 
6025c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
6035c6c1daeSBarry Smith 
604db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
605db781477SPatrick Sanan           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
606db781477SPatrick Sanan           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
6075c6c1daeSBarry Smith @*/
6085c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
6095c6c1daeSBarry Smith {
6105c6c1daeSBarry Smith   PetscFunctionBegin;
6119566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, hdf5v));
6129566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
6139566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
6149566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*hdf5v, name));
6159566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*hdf5v));
6165c6c1daeSBarry Smith   PetscFunctionReturn(0);
6175c6c1daeSBarry Smith }
6185c6c1daeSBarry Smith 
6195c6c1daeSBarry Smith /*@C
6205c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
6215c6c1daeSBarry Smith 
6225c6c1daeSBarry Smith   Not collective
6235c6c1daeSBarry Smith 
6245c6c1daeSBarry Smith   Input Parameter:
6255c6c1daeSBarry Smith . viewer - the PetscViewer
6265c6c1daeSBarry Smith 
6275c6c1daeSBarry Smith   Output Parameter:
6285c6c1daeSBarry Smith . file_id - The file id
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith   Level: intermediate
6315c6c1daeSBarry Smith 
632db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`
6335c6c1daeSBarry Smith @*/
6345c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
6355c6c1daeSBarry Smith {
6365c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
6375c6c1daeSBarry Smith 
6385c6c1daeSBarry Smith   PetscFunctionBegin;
6395c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
6405c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
6415c6c1daeSBarry Smith   PetscFunctionReturn(0);
6425c6c1daeSBarry Smith }
6435c6c1daeSBarry Smith 
6445c6c1daeSBarry Smith /*@C
6455c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
6465c6c1daeSBarry Smith 
6475c6c1daeSBarry Smith   Not collective
6485c6c1daeSBarry Smith 
6495c6c1daeSBarry Smith   Input Parameters:
6505c6c1daeSBarry Smith + viewer - the PetscViewer
6515c6c1daeSBarry Smith - name - The group name
6525c6c1daeSBarry Smith 
6535c6c1daeSBarry Smith   Level: intermediate
6545c6c1daeSBarry Smith 
6552e361344SVaclav Hapla   Notes:
6562e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
6572e361344SVaclav 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.
6582e361344SVaclav Hapla   . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
6592e361344SVaclav Hapla   - "." means the current group is pushed again.
6602e361344SVaclav Hapla 
6612e361344SVaclav Hapla   Example:
6622e361344SVaclav Hapla   Suppose the current group is "/a".
6632e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
6642e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
6652e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
6662e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
6672e361344SVaclav Hapla 
6682e361344SVaclav Hapla   Developer Notes:
6692e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
670820fc2d1SVaclav Hapla 
671c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
6725c6c1daeSBarry Smith @*/
673be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
6745c6c1daeSBarry Smith {
6755c6c1daeSBarry Smith   PetscViewer_HDF5          *hdf5 = (PetscViewer_HDF5*) viewer->data;
6767d964842SVaclav Hapla   PetscViewerHDF5GroupList  *groupNode;
6772e361344SVaclav Hapla   size_t                    i,len;
6782e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
6792e361344SVaclav Hapla   const char                *gname;
6805c6c1daeSBarry Smith 
6815c6c1daeSBarry Smith   PetscFunctionBegin;
6825c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
683820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name,2);
6849566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
6852e361344SVaclav Hapla   gname = NULL;
6862e361344SVaclav Hapla   if (len) {
6872e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
6882e361344SVaclav Hapla       /* use current name */
6892e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
6902e361344SVaclav Hapla     } else if (name[0] == '/') {
6912e361344SVaclav Hapla       /* absolute */
6922e361344SVaclav Hapla       for (i=1; i<len; i++) {
6932e361344SVaclav Hapla         if (name[i] != '/') {
6942e361344SVaclav Hapla           gname = name;
6952e361344SVaclav Hapla           break;
6962e361344SVaclav Hapla         }
6972e361344SVaclav Hapla       }
6982e361344SVaclav Hapla     } else {
6992e361344SVaclav Hapla       /* relative */
7002e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
7019566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
7022e361344SVaclav Hapla       gname = buf;
7032e361344SVaclav Hapla     }
7042e361344SVaclav Hapla   }
7059566063dSJacob Faibussowitsch   PetscCall(PetscNew(&groupNode));
7069566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(gname, (char**) &groupNode->name));
7075c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
7085c6c1daeSBarry Smith   hdf5->groups    = groupNode;
7095c6c1daeSBarry Smith   PetscFunctionReturn(0);
7105c6c1daeSBarry Smith }
7115c6c1daeSBarry Smith 
7123ef9c667SSatish Balay /*@
7135c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
7145c6c1daeSBarry Smith 
7155c6c1daeSBarry Smith   Not collective
7165c6c1daeSBarry Smith 
7175c6c1daeSBarry Smith   Input Parameter:
7185c6c1daeSBarry Smith . viewer - the PetscViewer
7195c6c1daeSBarry Smith 
7205c6c1daeSBarry Smith   Level: intermediate
7215c6c1daeSBarry Smith 
722c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
7235c6c1daeSBarry Smith @*/
7245c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5PopGroup(PetscViewer viewer)
7255c6c1daeSBarry Smith {
7265c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
7277d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
7285c6c1daeSBarry Smith 
7295c6c1daeSBarry Smith   PetscFunctionBegin;
7305c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
7315f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
7325c6c1daeSBarry Smith   groupNode    = hdf5->groups;
7335c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
7349566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
7359566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
7365c6c1daeSBarry Smith   PetscFunctionReturn(0);
7375c6c1daeSBarry Smith }
7385c6c1daeSBarry Smith 
7395c6c1daeSBarry Smith /*@C
740874270d9SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup().
741874270d9SVaclav Hapla   If none has been assigned, returns NULL.
7425c6c1daeSBarry Smith 
7435c6c1daeSBarry Smith   Not collective
7445c6c1daeSBarry Smith 
7455c6c1daeSBarry Smith   Input Parameter:
7465c6c1daeSBarry Smith . viewer - the PetscViewer
7475c6c1daeSBarry Smith 
7485c6c1daeSBarry Smith   Output Parameter:
7495c6c1daeSBarry Smith . name - The group name
7505c6c1daeSBarry Smith 
7515c6c1daeSBarry Smith   Level: intermediate
7525c6c1daeSBarry Smith 
753c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`
7545c6c1daeSBarry Smith @*/
755be7aa430SVaclav Hapla PetscErrorCode  PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[])
7565c6c1daeSBarry Smith {
7575c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data;
7585c6c1daeSBarry Smith 
7595c6c1daeSBarry Smith   PetscFunctionBegin;
7605c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
761c959eef4SJed Brown   PetscValidPointer(name,2);
762a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
7630298fd71SBarry Smith   else *name = NULL;
7645c6c1daeSBarry Smith   PetscFunctionReturn(0);
7655c6c1daeSBarry Smith }
7665c6c1daeSBarry Smith 
7678c557b5aSVaclav Hapla /*@
768874270d9SVaclav Hapla   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(),
769874270d9SVaclav Hapla   and return this group's ID and file ID.
770874270d9SVaclav Hapla   If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID.
771874270d9SVaclav Hapla 
772874270d9SVaclav Hapla   Not collective
773874270d9SVaclav Hapla 
774874270d9SVaclav Hapla   Input Parameter:
775874270d9SVaclav Hapla . viewer - the PetscViewer
776874270d9SVaclav Hapla 
777d8d19677SJose E. Roman   Output Parameters:
778874270d9SVaclav Hapla + fileId - The HDF5 file ID
779874270d9SVaclav Hapla - groupId - The HDF5 group ID
780874270d9SVaclav Hapla 
781e74616beSVaclav Hapla   Notes:
782e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
783e74616beSVaclav Hapla 
784874270d9SVaclav Hapla   Level: intermediate
785874270d9SVaclav Hapla 
786c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
787874270d9SVaclav Hapla @*/
78854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId)
78954dbf706SVaclav Hapla {
79090e03537SVaclav Hapla   hid_t          file_id;
79190e03537SVaclav Hapla   H5O_type_t     type;
79276d59af2SVaclav Hapla   const char     *groupName = NULL, *fileName = NULL;
79376d59af2SVaclav Hapla   PetscBool      writable, has;
79454dbf706SVaclav Hapla 
79554dbf706SVaclav Hapla   PetscFunctionBegin;
7969566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
7979566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
7989566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
7999566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName));
8009566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
80176d59af2SVaclav Hapla   if (!has) {
8025f80ce2aSJacob Faibussowitsch     PetscCheck(writable,PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName);
803f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
80476d59af2SVaclav Hapla   }
8055f80ce2aSJacob Faibussowitsch   PetscCheck(type == H5O_TYPE_GROUP,PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName);
806*792fecdfSBarry Smith   PetscCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT));
80754dbf706SVaclav Hapla   *fileId  = file_id;
80854dbf706SVaclav Hapla   PetscFunctionReturn(0);
80954dbf706SVaclav Hapla }
81054dbf706SVaclav Hapla 
8113ef9c667SSatish Balay /*@
812d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
8135c6c1daeSBarry Smith 
8145c6c1daeSBarry Smith   Not collective
8155c6c1daeSBarry Smith 
8165c6c1daeSBarry Smith   Input Parameter:
8175c6c1daeSBarry Smith . viewer - the PetscViewer
8185c6c1daeSBarry Smith 
8195c6c1daeSBarry Smith   Level: intermediate
8205c6c1daeSBarry Smith 
821d7dd068bSVaclav Hapla   Notes:
822d7dd068bSVaclav Hapla   On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0.
823d7dd068bSVaclav Hapla   Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep().
824d7dd068bSVaclav Hapla   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()].
825d7dd068bSVaclav Hapla   Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
826d7dd068bSVaclav Hapla   Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping().
827d7dd068bSVaclav Hapla 
828d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
829d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
830d7dd068bSVaclav Hapla 
831d7dd068bSVaclav Hapla   Developer notes:
832d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
833d7dd068bSVaclav Hapla 
834db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
835d7dd068bSVaclav Hapla @*/
836d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5PushTimestepping(PetscViewer viewer)
837d7dd068bSVaclav Hapla {
838d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
839d7dd068bSVaclav Hapla 
840d7dd068bSVaclav Hapla   PetscFunctionBegin;
841d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
8425f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
843d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
844d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
845d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
846d7dd068bSVaclav Hapla }
847d7dd068bSVaclav Hapla 
848d7dd068bSVaclav Hapla /*@
849d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
850d7dd068bSVaclav Hapla 
851d7dd068bSVaclav Hapla   Not collective
852d7dd068bSVaclav Hapla 
853d7dd068bSVaclav Hapla   Input Parameter:
854d7dd068bSVaclav Hapla . viewer - the PetscViewer
855d7dd068bSVaclav Hapla 
856d7dd068bSVaclav Hapla   Level: intermediate
857d7dd068bSVaclav Hapla 
858d7dd068bSVaclav Hapla   Notes:
859d7dd068bSVaclav Hapla   See PetscViewerHDF5PushTimestepping() for details.
860d7dd068bSVaclav Hapla 
861db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
862d7dd068bSVaclav Hapla @*/
863d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5PopTimestepping(PetscViewer viewer)
864d7dd068bSVaclav Hapla {
865d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
866d7dd068bSVaclav Hapla 
867d7dd068bSVaclav Hapla   PetscFunctionBegin;
868d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
8695f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
870d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
871d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
872d7dd068bSVaclav Hapla }
873d7dd068bSVaclav Hapla 
874d7dd068bSVaclav Hapla /*@
875d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
876d7dd068bSVaclav Hapla 
877d7dd068bSVaclav Hapla   Not collective
878d7dd068bSVaclav Hapla 
879d7dd068bSVaclav Hapla   Input Parameter:
880d7dd068bSVaclav Hapla . viewer - the PetscViewer
881d7dd068bSVaclav Hapla 
882d7dd068bSVaclav Hapla   Output Parameter:
883d7dd068bSVaclav Hapla . flg - is timestepping active?
884d7dd068bSVaclav Hapla 
885d7dd068bSVaclav Hapla   Level: intermediate
886d7dd068bSVaclav Hapla 
887d7dd068bSVaclav Hapla   Notes:
888d7dd068bSVaclav Hapla   See PetscViewerHDF5PushTimestepping() for details.
889d7dd068bSVaclav Hapla 
890db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
891d7dd068bSVaclav Hapla @*/
892d7dd068bSVaclav Hapla PetscErrorCode  PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
893d7dd068bSVaclav Hapla {
894d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
895d7dd068bSVaclav Hapla 
896d7dd068bSVaclav Hapla   PetscFunctionBegin;
897d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
898d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
899d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
900d7dd068bSVaclav Hapla }
901d7dd068bSVaclav Hapla 
902d7dd068bSVaclav Hapla /*@
903d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
904d7dd068bSVaclav Hapla 
905d7dd068bSVaclav Hapla   Not collective
906d7dd068bSVaclav Hapla 
907d7dd068bSVaclav Hapla   Input Parameter:
908d7dd068bSVaclav Hapla . viewer - the PetscViewer
909d7dd068bSVaclav Hapla 
910d7dd068bSVaclav Hapla   Level: intermediate
911d7dd068bSVaclav Hapla 
912d7dd068bSVaclav Hapla   Notes:
913d7dd068bSVaclav Hapla   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
914d7dd068bSVaclav Hapla 
915db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
9165c6c1daeSBarry Smith @*/
9175c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
9185c6c1daeSBarry Smith {
9195c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
9205c6c1daeSBarry Smith 
9215c6c1daeSBarry Smith   PetscFunctionBegin;
9225c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
9235f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9245c6c1daeSBarry Smith   ++hdf5->timestep;
9255c6c1daeSBarry Smith   PetscFunctionReturn(0);
9265c6c1daeSBarry Smith }
9275c6c1daeSBarry Smith 
9283ef9c667SSatish Balay /*@
929d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
9305c6c1daeSBarry Smith 
931d7dd068bSVaclav Hapla   Logically collective
9325c6c1daeSBarry Smith 
9335c6c1daeSBarry Smith   Input Parameters:
9345c6c1daeSBarry Smith + viewer - the PetscViewer
935d7dd068bSVaclav Hapla - timestep - The timestep
9365c6c1daeSBarry Smith 
9375c6c1daeSBarry Smith   Level: intermediate
9385c6c1daeSBarry Smith 
939d7dd068bSVaclav Hapla   Notes:
940d7dd068bSVaclav Hapla   This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
941d7dd068bSVaclav Hapla 
942db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
9435c6c1daeSBarry Smith @*/
9445c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
9455c6c1daeSBarry Smith {
9465c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
9475c6c1daeSBarry Smith 
9485c6c1daeSBarry Smith   PetscFunctionBegin;
9495c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
950d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
9515f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
9525f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9535c6c1daeSBarry Smith   hdf5->timestep = timestep;
9545c6c1daeSBarry Smith   PetscFunctionReturn(0);
9555c6c1daeSBarry Smith }
9565c6c1daeSBarry Smith 
9573ef9c667SSatish Balay /*@
9585c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
9595c6c1daeSBarry Smith 
9605c6c1daeSBarry Smith   Not collective
9615c6c1daeSBarry Smith 
9625c6c1daeSBarry Smith   Input Parameter:
9635c6c1daeSBarry Smith . viewer - the PetscViewer
9645c6c1daeSBarry Smith 
9655c6c1daeSBarry Smith   Output Parameter:
966d7dd068bSVaclav Hapla . timestep - The timestep
9675c6c1daeSBarry Smith 
9685c6c1daeSBarry Smith   Level: intermediate
9695c6c1daeSBarry Smith 
970d7dd068bSVaclav Hapla   Notes:
971d7dd068bSVaclav Hapla   This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details.
972d7dd068bSVaclav Hapla 
973db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
9745c6c1daeSBarry Smith @*/
9755c6c1daeSBarry Smith PetscErrorCode  PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
9765c6c1daeSBarry Smith {
9775c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
9785c6c1daeSBarry Smith 
9795c6c1daeSBarry Smith   PetscFunctionBegin;
9805c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
981d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep,2);
9825f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9835c6c1daeSBarry Smith   *timestep = hdf5->timestep;
9845c6c1daeSBarry Smith   PetscFunctionReturn(0);
9855c6c1daeSBarry Smith }
9865c6c1daeSBarry Smith 
98736bce6e8SMatthew G. Knepley /*@C
98836bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
98936bce6e8SMatthew G. Knepley 
99036bce6e8SMatthew G. Knepley   Not collective
99136bce6e8SMatthew G. Knepley 
99236bce6e8SMatthew G. Knepley   Input Parameter:
99336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
99436bce6e8SMatthew G. Knepley 
99536bce6e8SMatthew G. Knepley   Output Parameter:
99636bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...)
99736bce6e8SMatthew G. Knepley 
99836bce6e8SMatthew G. Knepley   Level: advanced
99936bce6e8SMatthew G. Knepley 
1000db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
100136bce6e8SMatthew G. Knepley @*/
100236bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
100336bce6e8SMatthew G. Knepley {
100436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
100536bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
100636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
100736bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_LLONG;
100836bce6e8SMatthew G. Knepley #else
100936bce6e8SMatthew G. Knepley                                        *htype = H5T_NATIVE_INT;
101036bce6e8SMatthew G. Knepley #endif
101136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE)      *htype = H5T_NATIVE_DOUBLE;
101236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG)        *htype = H5T_NATIVE_LONG;
101336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT)       *htype = H5T_NATIVE_SHORT;
1014c9450970SBarry Smith   else if (ptype == PETSC_ENUM)        *htype = H5T_NATIVE_INT;
1015de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL)        *htype = H5T_NATIVE_INT;
101636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT)       *htype = H5T_NATIVE_FLOAT;
101736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR)        *htype = H5T_NATIVE_CHAR;
101836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
10197e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING)      *htype = H5Tcopy(H5T_C_S1);
102036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
102136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
102236bce6e8SMatthew G. Knepley }
102336bce6e8SMatthew G. Knepley 
102436bce6e8SMatthew G. Knepley /*@C
102536bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
102636bce6e8SMatthew G. Knepley 
102736bce6e8SMatthew G. Knepley   Not collective
102836bce6e8SMatthew G. Knepley 
102936bce6e8SMatthew G. Knepley   Input Parameter:
103036bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...)
103136bce6e8SMatthew G. Knepley 
103236bce6e8SMatthew G. Knepley   Output Parameter:
103336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE)
103436bce6e8SMatthew G. Knepley 
103536bce6e8SMatthew G. Knepley   Level: advanced
103636bce6e8SMatthew G. Knepley 
1037db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
103836bce6e8SMatthew G. Knepley @*/
103936bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
104036bce6e8SMatthew G. Knepley {
104136bce6e8SMatthew G. Knepley   PetscFunctionBegin;
104236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
104336bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_LONG;
104436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG)  *ptype = PETSC_INT;
104536bce6e8SMatthew G. Knepley #else
104636bce6e8SMatthew G. Knepley   if      (htype == H5T_NATIVE_INT)    *ptype = PETSC_INT;
104736bce6e8SMatthew G. Knepley #endif
104836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
104936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG)   *ptype = PETSC_LONG;
105036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT)  *ptype = PETSC_SHORT;
105136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT)  *ptype = PETSC_FLOAT;
105236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR)   *ptype = PETSC_CHAR;
105336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR)  *ptype = PETSC_CHAR;
10547e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1)          *ptype = PETSC_STRING;
105536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
105636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
105736bce6e8SMatthew G. Knepley }
105836bce6e8SMatthew G. Knepley 
1059df863907SAlex Fikl /*@C
1060b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
106136bce6e8SMatthew G. Knepley 
1062343a20b2SVaclav Hapla   Collective
1063343a20b2SVaclav Hapla 
106436bce6e8SMatthew G. Knepley   Input Parameters:
106536bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer
10664302210dSVaclav Hapla . parent - The parent dataset/group name
106736bce6e8SMatthew G. Knepley . name   - The attribute name
106836bce6e8SMatthew G. Knepley . datatype - The attribute type
106936bce6e8SMatthew G. Knepley - value    - The attribute value
107036bce6e8SMatthew G. Knepley 
107136bce6e8SMatthew G. Knepley   Level: advanced
107236bce6e8SMatthew G. Knepley 
10734302210dSVaclav Hapla   Notes:
1074343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
10754302210dSVaclav Hapla 
1076c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
107736bce6e8SMatthew G. Knepley @*/
10784302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
107936bce6e8SMatthew G. Knepley {
10804302210dSVaclav Hapla   char           *parentAbsPath;
108160568592SMatthew G. Knepley   hid_t          h5, dataspace, obj, attribute, dtype;
1082080f144cSVaclav Hapla   PetscBool      has;
108336bce6e8SMatthew G. Knepley 
108436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
10855cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
10864302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1087c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
10884302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer,datatype,4);
1089b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
10909566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
10919566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
10929566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
10939566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
10947e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
10957e97c476SMatthew G. Knepley     size_t len;
10969566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *) value, &len));
1097*792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size,(dtype, len+1));
10987e97c476SMatthew G. Knepley   }
10999566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1100*792fecdfSBarry Smith   PetscCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR));
1101*792fecdfSBarry Smith   PetscCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1102080f144cSVaclav Hapla   if (has) {
1103*792fecdfSBarry Smith     PetscCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1104080f144cSVaclav Hapla   } else {
1105*792fecdfSBarry Smith     PetscCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1106080f144cSVaclav Hapla   }
1107*792fecdfSBarry Smith   PetscCallHDF5(H5Awrite,(attribute, dtype, value));
1108*792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose,(dtype));
1109*792fecdfSBarry Smith   PetscCallHDF5(H5Aclose,(attribute));
1110*792fecdfSBarry Smith   PetscCallHDF5(H5Oclose,(obj));
1111*792fecdfSBarry Smith   PetscCallHDF5(H5Sclose,(dataspace));
11129566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
111336bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
111436bce6e8SMatthew G. Knepley }
111536bce6e8SMatthew G. Knepley 
1116df863907SAlex Fikl /*@C
1117577e0e04SVaclav Hapla  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name
1118577e0e04SVaclav Hapla 
1119343a20b2SVaclav Hapla   Collective
1120343a20b2SVaclav Hapla 
1121577e0e04SVaclav Hapla   Input Parameters:
1122577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
1123577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1124577e0e04SVaclav Hapla . name     - The attribute name
1125577e0e04SVaclav Hapla . datatype - The attribute type
1126577e0e04SVaclav Hapla - value    - The attribute value
1127577e0e04SVaclav Hapla 
1128577e0e04SVaclav Hapla   Notes:
1129343a20b2SVaclav Hapla   This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1130577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1131577e0e04SVaclav Hapla 
1132577e0e04SVaclav Hapla   Level: advanced
1133577e0e04SVaclav Hapla 
1134c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1135577e0e04SVaclav Hapla @*/
1136577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1137577e0e04SVaclav Hapla {
1138577e0e04SVaclav Hapla   PetscFunctionBegin;
1139577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1140577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1141577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1142b7e81f0fSVaclav Hapla   PetscValidPointer(value,5);
11439566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
11449566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1145577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1146577e0e04SVaclav Hapla }
1147577e0e04SVaclav Hapla 
1148577e0e04SVaclav Hapla /*@C
1149b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1150ad0c61feSMatthew G. Knepley 
1151343a20b2SVaclav Hapla   Collective
1152343a20b2SVaclav Hapla 
1153ad0c61feSMatthew G. Knepley   Input Parameters:
1154ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer
11554302210dSVaclav Hapla . parent - The parent dataset/group name
1156ad0c61feSMatthew G. Knepley . name   - The attribute name
1157a2d6be1bSVaclav Hapla . datatype - The attribute type
1158a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1159ad0c61feSMatthew G. Knepley 
1160ad0c61feSMatthew G. Knepley   Output Parameter:
1161a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1162ad0c61feSMatthew G. Knepley 
1163a2d6be1bSVaclav Hapla   Notes:
1164a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1165a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1166a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1167a2d6be1bSVaclav Hapla $  flg = PETSC_FALSE;
11689566063dSJacob Faibussowitsch $  PetscCall(PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1169a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1170a2d6be1bSVaclav Hapla 
1171a2d6be1bSVaclav Hapla   If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed.
1172708d4cb3SBarry Smith 
1173343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
1174ad0c61feSMatthew G. Knepley 
1175343a20b2SVaclav Hapla   Level: advanced
11764302210dSVaclav Hapla 
1177c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1178ad0c61feSMatthew G. Knepley @*/
1179d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1180ad0c61feSMatthew G. Knepley {
11814302210dSVaclav Hapla   char           *parentAbsPath;
1182a2d6be1bSVaclav Hapla   hid_t          h5, obj, attribute, dtype;
1183969748fdSVaclav Hapla   PetscBool      has;
1184ad0c61feSMatthew G. Knepley 
1185ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
11865cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
11874302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1188c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1189a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1190a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
11919566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11929566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
11939566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
11949566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1195a2d6be1bSVaclav Hapla   if (!has) {
1196a2d6be1bSVaclav Hapla     if (defaultValue) {
1197a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1198a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
11999566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char**)defaultValue, (char**)value));
1200a2d6be1bSVaclav Hapla         } else {
1201a2d6be1bSVaclav Hapla           size_t len;
1202*792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len,H5Tget_size,(dtype));
12039566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1204a2d6be1bSVaclav Hapla         }
1205a2d6be1bSVaclav Hapla       }
12069566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
1207a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
120898921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1209a2d6be1bSVaclav Hapla   }
12109566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1211*792fecdfSBarry Smith   PetscCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT));
1212*792fecdfSBarry Smith   PetscCallHDF5Return(attribute,H5Aopen_name,(obj, name));
1213f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1214f0b58479SMatthew G. Knepley     size_t len;
1215a2d6be1bSVaclav Hapla     hid_t  atype;
1216*792fecdfSBarry Smith     PetscCallHDF5Return(atype,H5Aget_type,(attribute));
1217*792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len,H5Tget_size,(atype));
12189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len+1) * sizeof(char), value));
1219*792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size,(dtype, len+1));
1220*792fecdfSBarry Smith     PetscCallHDF5(H5Aread,(attribute, dtype, *(char**)value));
1221708d4cb3SBarry Smith   } else {
1222*792fecdfSBarry Smith     PetscCallHDF5(H5Aread,(attribute, dtype, value));
1223708d4cb3SBarry Smith   }
1224*792fecdfSBarry Smith   PetscCallHDF5(H5Aclose,(attribute));
1225e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1226*792fecdfSBarry Smith   PetscCallHDF5(H5Oclose,(obj));
12279566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
1228ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1229ad0c61feSMatthew G. Knepley }
1230ad0c61feSMatthew G. Knepley 
1231577e0e04SVaclav Hapla /*@C
1232577e0e04SVaclav Hapla  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name
1233577e0e04SVaclav Hapla 
1234343a20b2SVaclav Hapla   Collective
1235343a20b2SVaclav Hapla 
1236577e0e04SVaclav Hapla   Input Parameters:
1237577e0e04SVaclav Hapla + viewer   - The HDF5 viewer
1238577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1239577e0e04SVaclav Hapla . name     - The attribute name
1240577e0e04SVaclav Hapla - datatype - The attribute type
1241577e0e04SVaclav Hapla 
1242577e0e04SVaclav Hapla   Output Parameter:
1243577e0e04SVaclav Hapla . value    - The attribute value
1244577e0e04SVaclav Hapla 
1245577e0e04SVaclav Hapla   Notes:
1246577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1247577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1248577e0e04SVaclav Hapla 
1249577e0e04SVaclav Hapla   Level: advanced
1250577e0e04SVaclav Hapla 
1251c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1252577e0e04SVaclav Hapla @*/
1253a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1254577e0e04SVaclav Hapla {
1255577e0e04SVaclav Hapla   PetscFunctionBegin;
1256577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1257577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1258577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
1259064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
12609566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12619566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1262577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1263577e0e04SVaclav Hapla }
1264577e0e04SVaclav Hapla 
12659fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1266a75c3fd4SVaclav Hapla {
1267a75c3fd4SVaclav Hapla   htri_t exists;
1268a75c3fd4SVaclav Hapla   hid_t group;
1269a75c3fd4SVaclav Hapla 
1270a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1271*792fecdfSBarry Smith   PetscCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT));
1272*792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT));
1273a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1274*792fecdfSBarry Smith     PetscCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1275*792fecdfSBarry Smith     PetscCallHDF5(H5Gclose,(group));
1276a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1277a75c3fd4SVaclav Hapla   }
1278a75c3fd4SVaclav Hapla   *exists_ = (PetscBool) exists;
1279a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1280a75c3fd4SVaclav Hapla }
1281a75c3fd4SVaclav Hapla 
1282bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
12835cdeb108SMatthew G. Knepley {
128490e03537SVaclav Hapla   const char     rootGroupName[] = "/";
12855cdeb108SMatthew G. Knepley   hid_t          h5;
1286e5bf9ebcSVaclav Hapla   PetscBool      exists=PETSC_FALSE;
128715b861d2SVaclav Hapla   PetscInt       i;
128815b861d2SVaclav Hapla   int            n;
128985688ae2SVaclav Hapla   char           **hierarchy;
129085688ae2SVaclav Hapla   char           buf[PETSC_MAX_PATH_LEN]="";
12915cdeb108SMatthew G. Knepley 
12925cdeb108SMatthew G. Knepley   PetscFunctionBegin;
12935cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
129490e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
129590e03537SVaclav Hapla   else name = rootGroupName;
1296ccd66a83SVaclav Hapla   if (has) {
1297064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
12985cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1299ccd66a83SVaclav Hapla   }
1300ccd66a83SVaclav Hapla   if (otype) {
1301064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
130256cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1303ccd66a83SVaclav Hapla   }
13049566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
130585688ae2SVaclav Hapla 
130685688ae2SVaclav Hapla   /*
130785688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
130885688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
130985688ae2SVaclav Hapla      1) whether it's a valid link
131085688ae2SVaclav Hapla      2) whether this link resolves to an object
131185688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
131285688ae2SVaclav Hapla   */
13139566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name,'/',&n,&hierarchy));
131485688ae2SVaclav Hapla   if (!n) {
131585688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1316ccd66a83SVaclav Hapla     if (has)   *has   = PETSC_TRUE;
1317ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
13189566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n,hierarchy));
131985688ae2SVaclav Hapla     PetscFunctionReturn(0);
132085688ae2SVaclav Hapla   }
132185688ae2SVaclav Hapla   for (i=0; i<n; i++) {
13229566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf,"/"));
13239566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf,hierarchy[i]));
13249566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1325a75c3fd4SVaclav Hapla     if (!exists) break;
132685688ae2SVaclav Hapla   }
13279566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n,hierarchy));
132885688ae2SVaclav Hapla 
132985688ae2SVaclav Hapla   /* If the object exists, get its type */
1330ccd66a83SVaclav Hapla   if (exists && otype) {
13315cdeb108SMatthew G. Knepley     H5O_info_t info;
13325cdeb108SMatthew G. Knepley 
133376276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1334*792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT));
133556cc0592SVaclav Hapla     *otype = info.type;
13365cdeb108SMatthew G. Knepley   }
1337ccd66a83SVaclav Hapla   if (has) *has = exists;
13385cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
13395cdeb108SMatthew G. Knepley }
13405cdeb108SMatthew G. Knepley 
13414302210dSVaclav Hapla /*@C
13420a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
13430a9f38efSVaclav Hapla 
1344343a20b2SVaclav Hapla   Collective
1345343a20b2SVaclav Hapla 
13460a9f38efSVaclav Hapla   Input Parameters:
13474302210dSVaclav Hapla + viewer - The HDF5 viewer
13484302210dSVaclav Hapla - path - The group path
13490a9f38efSVaclav Hapla 
13500a9f38efSVaclav Hapla   Output Parameter:
13510a9f38efSVaclav Hapla . has - Flag for group existence
13520a9f38efSVaclav Hapla 
13530a9f38efSVaclav Hapla   Level: advanced
13540a9f38efSVaclav Hapla 
13554302210dSVaclav Hapla   Notes:
1356785c443eSVaclav Hapla   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
1357785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
1358343a20b2SVaclav Hapla   If path exists but is not a group, PETSC_FALSE is returned.
13594302210dSVaclav Hapla 
1360db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
13610a9f38efSVaclav Hapla @*/
13624302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
13630a9f38efSVaclav Hapla {
13640a9f38efSVaclav Hapla   H5O_type_t     type;
13654302210dSVaclav Hapla   char           *abspath;
13660a9f38efSVaclav Hapla 
13670a9f38efSVaclav Hapla   PetscFunctionBegin;
13680a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
13694302210dSVaclav Hapla   if (path) PetscValidCharPointer(path,2);
13704302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
13719566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
13729566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
13734302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
13749566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
13750a9f38efSVaclav Hapla   PetscFunctionReturn(0);
13760a9f38efSVaclav Hapla }
13770a9f38efSVaclav Hapla 
137889e0ef10SVaclav Hapla /*@C
137989e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
138089e0ef10SVaclav Hapla 
1381343a20b2SVaclav Hapla   Collective
1382343a20b2SVaclav Hapla 
138389e0ef10SVaclav Hapla   Input Parameters:
138489e0ef10SVaclav Hapla + viewer - The HDF5 viewer
138589e0ef10SVaclav Hapla - path - The dataset path
138689e0ef10SVaclav Hapla 
138789e0ef10SVaclav Hapla   Output Parameter:
138889e0ef10SVaclav Hapla . has - Flag whether dataset exists
138989e0ef10SVaclav Hapla 
139089e0ef10SVaclav Hapla   Level: advanced
139189e0ef10SVaclav Hapla 
139289e0ef10SVaclav Hapla   Notes:
1393343a20b2SVaclav Hapla   If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group.
139489e0ef10SVaclav Hapla   If path is NULL or empty, has is set to PETSC_FALSE.
1395343a20b2SVaclav Hapla   If path exists but is not a dataset, has is set to PETSC_FALSE as well.
139689e0ef10SVaclav Hapla 
1397db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
139889e0ef10SVaclav Hapla @*/
139989e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
140089e0ef10SVaclav Hapla {
140189e0ef10SVaclav Hapla   H5O_type_t     type;
140289e0ef10SVaclav Hapla   char           *abspath;
140389e0ef10SVaclav Hapla 
140489e0ef10SVaclav Hapla   PetscFunctionBegin;
140589e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
140689e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path,2);
140789e0ef10SVaclav Hapla   PetscValidBoolPointer(has,3);
14089566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
14099566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
141089e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
14119566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
141289e0ef10SVaclav Hapla   PetscFunctionReturn(0);
141389e0ef10SVaclav Hapla }
141489e0ef10SVaclav Hapla 
14150a9f38efSVaclav Hapla /*@
1416e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1417ecce1506SVaclav Hapla 
1418343a20b2SVaclav Hapla   Collective
1419343a20b2SVaclav Hapla 
1420ecce1506SVaclav Hapla   Input Parameters:
1421ecce1506SVaclav Hapla + viewer - The HDF5 viewer
1422ecce1506SVaclav Hapla - obj    - The named object
1423ecce1506SVaclav Hapla 
1424ecce1506SVaclav Hapla   Output Parameter:
142589e0ef10SVaclav Hapla . has    - Flag for dataset existence
1426ecce1506SVaclav Hapla 
1427e3f143f7SVaclav Hapla   Notes:
142889e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1429343a20b2SVaclav Hapla   If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well.
1430e3f143f7SVaclav Hapla 
1431ecce1506SVaclav Hapla   Level: advanced
1432ecce1506SVaclav Hapla 
1433db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1434ecce1506SVaclav Hapla @*/
1435ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1436ecce1506SVaclav Hapla {
143789e0ef10SVaclav Hapla   size_t         len;
1438ecce1506SVaclav Hapla 
1439ecce1506SVaclav Hapla   PetscFunctionBegin;
1440c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1441c57a1a9aSVaclav Hapla   PetscValidHeader(obj,2);
14424302210dSVaclav Hapla   PetscValidBoolPointer(has,3);
14439566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
14445f80ce2aSJacob Faibussowitsch   PetscCheck(len,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
14459566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1446ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1447ecce1506SVaclav Hapla }
1448ecce1506SVaclav Hapla 
1449df863907SAlex Fikl /*@C
1450ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1451e7bdbf8eSMatthew G. Knepley 
1452343a20b2SVaclav Hapla   Collective
1453343a20b2SVaclav Hapla 
1454e7bdbf8eSMatthew G. Knepley   Input Parameters:
1455e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer
14564302210dSVaclav Hapla . parent - The parent dataset/group name
1457e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1458e7bdbf8eSMatthew G. Knepley 
1459e7bdbf8eSMatthew G. Knepley   Output Parameter:
1460e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1461e7bdbf8eSMatthew G. Knepley 
1462e7bdbf8eSMatthew G. Knepley   Level: advanced
1463e7bdbf8eSMatthew G. Knepley 
14644302210dSVaclav Hapla   Notes:
1465343a20b2SVaclav Hapla   If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group.
14664302210dSVaclav Hapla 
1467c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1468e7bdbf8eSMatthew G. Knepley @*/
14694302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1470e7bdbf8eSMatthew G. Knepley {
14714302210dSVaclav Hapla   char           *parentAbsPath;
1472e7bdbf8eSMatthew G. Knepley 
1473e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
14745cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
14754302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent,2);
1476c57a1a9aSVaclav Hapla   PetscValidCharPointer(name,3);
14774302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
14789566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
14799566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
14809566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
14819566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
148206db490cSVaclav Hapla   PetscFunctionReturn(0);
148306db490cSVaclav Hapla }
148406db490cSVaclav Hapla 
1485577e0e04SVaclav Hapla /*@C
1486577e0e04SVaclav Hapla  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name
1487577e0e04SVaclav Hapla 
1488343a20b2SVaclav Hapla   Collective
1489343a20b2SVaclav Hapla 
1490577e0e04SVaclav Hapla   Input Parameters:
1491577e0e04SVaclav Hapla + viewer - The HDF5 viewer
1492577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1493577e0e04SVaclav Hapla - name   - The attribute name
1494577e0e04SVaclav Hapla 
1495577e0e04SVaclav Hapla   Output Parameter:
1496577e0e04SVaclav Hapla . has    - Flag for attribute existence
1497577e0e04SVaclav Hapla 
1498577e0e04SVaclav Hapla   Notes:
1499577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1500577e0e04SVaclav Hapla   You might want to check first if it does using PetscViewerHDF5HasObject().
1501577e0e04SVaclav Hapla 
1502577e0e04SVaclav Hapla   Level: advanced
1503577e0e04SVaclav Hapla 
1504c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1505577e0e04SVaclav Hapla @*/
1506577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1507577e0e04SVaclav Hapla {
1508577e0e04SVaclav Hapla   PetscFunctionBegin;
1509577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1);
1510577e0e04SVaclav Hapla   PetscValidHeader(obj,2);
1511577e0e04SVaclav Hapla   PetscValidCharPointer(name,3);
15124302210dSVaclav Hapla   PetscValidBoolPointer(has,4);
15139566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
15149566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1515577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1516577e0e04SVaclav Hapla }
1517577e0e04SVaclav Hapla 
151806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
151906db490cSVaclav Hapla {
152006db490cSVaclav Hapla   hid_t          h5;
152106db490cSVaclav Hapla   htri_t         hhas;
152206db490cSVaclav Hapla 
152306db490cSVaclav Hapla   PetscFunctionBegin;
15249566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1525*792fecdfSBarry Smith   PetscCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT));
15265cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1527e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1528e7bdbf8eSMatthew G. Knepley }
1529e7bdbf8eSMatthew G. Knepley 
1530a75e6a4aSMatthew G. Knepley /*
1531a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1532a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1533a75e6a4aSMatthew G. Knepley */
1534d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1535a75e6a4aSMatthew G. Knepley 
1536a75e6a4aSMatthew G. Knepley /*@C
1537a75e6a4aSMatthew G. Knepley   PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator.
1538a75e6a4aSMatthew G. Knepley 
1539d083f849SBarry Smith   Collective
1540a75e6a4aSMatthew G. Knepley 
1541a75e6a4aSMatthew G. Knepley   Input Parameter:
1542a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer
1543a75e6a4aSMatthew G. Knepley 
1544a75e6a4aSMatthew G. Knepley   Level: intermediate
1545a75e6a4aSMatthew G. Knepley 
1546a75e6a4aSMatthew G. Knepley   Options Database Keys:
154710699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1548a75e6a4aSMatthew G. Knepley 
1549a75e6a4aSMatthew G. Knepley   Environmental variables:
155010699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file
1551a75e6a4aSMatthew G. Knepley 
1552a75e6a4aSMatthew G. Knepley   Notes:
1553a75e6a4aSMatthew G. Knepley   Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return
1554a75e6a4aSMatthew G. Knepley   an error code.  The HDF5 PetscViewer is usually used in the form
1555a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1556a75e6a4aSMatthew G. Knepley 
1557db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1558a75e6a4aSMatthew G. Knepley @*/
1559a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1560a75e6a4aSMatthew G. Knepley {
1561a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1562a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1563a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1564a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1565a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1566a75e6a4aSMatthew G. Knepley 
1567a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
1568908793a3SLisandro 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);}
1569a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1570908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL);
1571908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1572a75e6a4aSMatthew G. Knepley   }
157347435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg);
1574908793a3SLisandro Dalcin   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1575a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1576a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg);
15772cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1578a75e6a4aSMatthew G. Knepley     if (!flg) {
1579a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname,"output.h5");
15802cb5e1ccSBarry Smith       if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1581a75e6a4aSMatthew G. Knepley     }
1582a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer);
15832cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1584a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15852cb5e1ccSBarry Smith     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
158647435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer);
1587908793a3SLisandro Dalcin     if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);}
1588a75e6a4aSMatthew G. Knepley   }
1589a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
15902cb5e1ccSBarry Smith   if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);}
1591a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1592a75e6a4aSMatthew G. Knepley }
1593