xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
69371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) {
74302210dSVaclav Hapla   PetscBool   relative = PETSC_FALSE;
86c132bc1SVaclav Hapla   const char *group;
96c132bc1SVaclav Hapla   char        buf[PETSC_MAX_PATH_LEN] = "";
106c132bc1SVaclav Hapla 
116c132bc1SVaclav Hapla   PetscFunctionBegin;
129566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetGroup(viewer, &group));
139566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative));
144302210dSVaclav Hapla   if (relative) {
159566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(buf, group));
169566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
179566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, path));
189566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(buf, abspath));
194302210dSVaclav Hapla   } else {
209566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(path, abspath));
214302210dSVaclav Hapla   }
226c132bc1SVaclav Hapla   PetscFunctionReturn(0);
236c132bc1SVaclav Hapla }
246c132bc1SVaclav Hapla 
259371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) {
26577e0e04SVaclav Hapla   PetscBool has;
27577e0e04SVaclav Hapla 
28577e0e04SVaclav Hapla   PetscFunctionBegin;
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
30577e0e04SVaclav Hapla   if (!has) {
3189e0ef10SVaclav Hapla     const char *group;
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5GetGroup(viewer, &group));
3398921bdaSJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/");
34577e0e04SVaclav Hapla   }
35577e0e04SVaclav Hapla   PetscFunctionReturn(0);
36577e0e04SVaclav Hapla }
37577e0e04SVaclav Hapla 
389371c9d4SSatish Balay static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject) {
39ee8b9732SVaclav Hapla   PetscBool         flg  = PETSC_FALSE, set;
4082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
4182ea9b62SBarry Smith 
4282ea9b62SBarry Smith   PetscFunctionBegin;
43d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
479566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
4819a20e4cSMatthew G. Knepley   flg = PETSC_FALSE;
499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
509566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
51d0609cedSBarry Smith   PetscOptionsHeadEnd();
5282ea9b62SBarry Smith   PetscFunctionReturn(0);
5382ea9b62SBarry Smith }
5482ea9b62SBarry Smith 
559371c9d4SSatish Balay static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer) {
561b793a25SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
5703fa8834SVaclav Hapla   PetscBool         flg;
581b793a25SVaclav Hapla 
591b793a25SVaclav Hapla   PetscFunctionBegin;
6048a46eb9SPierre Jolivet   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
619566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
629566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
639566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
649566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
661b793a25SVaclav Hapla   PetscFunctionReturn(0);
671b793a25SVaclav Hapla }
681b793a25SVaclav Hapla 
699371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) {
705c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
715c6c1daeSBarry Smith 
725c6c1daeSBarry Smith   PetscFunctionBegin;
739566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5->filename));
74792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
755c6c1daeSBarry Smith   PetscFunctionReturn(0);
765c6c1daeSBarry Smith }
775c6c1daeSBarry Smith 
789371c9d4SSatish Balay static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) {
796226335aSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
806226335aSVaclav Hapla 
816226335aSVaclav Hapla   PetscFunctionBegin;
82792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
836226335aSVaclav Hapla   PetscFunctionReturn(0);
846226335aSVaclav Hapla }
856226335aSVaclav Hapla 
869371c9d4SSatish Balay static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) {
875c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
885c6c1daeSBarry Smith 
895c6c1daeSBarry Smith   PetscFunctionBegin;
90792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
919566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_HDF5(viewer));
925c6c1daeSBarry Smith   while (hdf5->groups) {
937d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
945c6c1daeSBarry Smith 
959566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups->name));
969566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups));
975c6c1daeSBarry Smith     hdf5->groups = tmp;
985c6c1daeSBarry Smith   }
999566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5));
1009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
1029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
1032e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
1049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
1059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
1062e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
1072e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
1089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
1099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
1105c6c1daeSBarry Smith   PetscFunctionReturn(0);
1115c6c1daeSBarry Smith }
1125c6c1daeSBarry Smith 
1139371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) {
1145c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1155c6c1daeSBarry Smith 
1165c6c1daeSBarry Smith   PetscFunctionBegin;
1175c6c1daeSBarry Smith   hdf5->btype = type;
1185c6c1daeSBarry Smith   PetscFunctionReturn(0);
1195c6c1daeSBarry Smith }
1205c6c1daeSBarry Smith 
1219371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) {
1220b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1230b62783dSVaclav Hapla 
1240b62783dSVaclav Hapla   PetscFunctionBegin;
1250b62783dSVaclav Hapla   *type = hdf5->btype;
1260b62783dSVaclav Hapla   PetscFunctionReturn(0);
1270b62783dSVaclav Hapla }
1280b62783dSVaclav Hapla 
1299371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) {
13082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
13182ea9b62SBarry Smith 
13282ea9b62SBarry Smith   PetscFunctionBegin;
13382ea9b62SBarry Smith   hdf5->basedimension2 = flg;
13482ea9b62SBarry Smith   PetscFunctionReturn(0);
13582ea9b62SBarry Smith }
13682ea9b62SBarry Smith 
137df863907SAlex Fikl /*@
13882ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
13982ea9b62SBarry Smith        dimension of 2.
14082ea9b62SBarry Smith 
141811af0c4SBarry Smith     Logically Collective on viewer
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith   Input Parameters:
144811af0c4SBarry Smith +  viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
145811af0c4SBarry 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
14682ea9b62SBarry Smith 
147811af0c4SBarry Smith   Options Database Key:
14882ea9b62SBarry 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
14982ea9b62SBarry Smith 
150811af0c4SBarry Smith   Note:
15195452b02SPatrick 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
15282ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
15382ea9b62SBarry Smith 
15482ea9b62SBarry Smith   Level: intermediate
15582ea9b62SBarry Smith 
156811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
15782ea9b62SBarry Smith @*/
1589371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg) {
15982ea9b62SBarry Smith   PetscFunctionBegin;
16082ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
161cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
16282ea9b62SBarry Smith   PetscFunctionReturn(0);
16382ea9b62SBarry Smith }
16482ea9b62SBarry Smith 
165df863907SAlex Fikl /*@
16682ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
16782ea9b62SBarry Smith        dimension of 2.
16882ea9b62SBarry Smith 
169811af0c4SBarry Smith     Logically Collective on viewer
17082ea9b62SBarry Smith 
17182ea9b62SBarry Smith   Input Parameter:
172811af0c4SBarry Smith .  viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`
17382ea9b62SBarry Smith 
17482ea9b62SBarry Smith   Output Parameter:
175811af0c4SBarry 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
17682ea9b62SBarry Smith 
177811af0c4SBarry Smith   Note:
17895452b02SPatrick 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
17982ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
18082ea9b62SBarry Smith 
18182ea9b62SBarry Smith   Level: intermediate
18282ea9b62SBarry Smith 
183811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
18482ea9b62SBarry Smith @*/
1859371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg) {
18682ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
18782ea9b62SBarry Smith 
18882ea9b62SBarry Smith   PetscFunctionBegin;
18982ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
19082ea9b62SBarry Smith   *flg = hdf5->basedimension2;
19182ea9b62SBarry Smith   PetscFunctionReturn(0);
19282ea9b62SBarry Smith }
19382ea9b62SBarry Smith 
1949371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) {
1959a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1969a0502c6SHåkon Strandenes 
1979a0502c6SHåkon Strandenes   PetscFunctionBegin;
1989a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
1999a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2009a0502c6SHåkon Strandenes }
2019a0502c6SHåkon Strandenes 
202df863907SAlex Fikl /*@
2039a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
204811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2059a0502c6SHåkon Strandenes 
206811af0c4SBarry Smith     Logically Collective on viewer
2079a0502c6SHåkon Strandenes 
2089a0502c6SHåkon Strandenes   Input Parameters:
209811af0c4SBarry Smith +  viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
210811af0c4SBarry Smith -  flg - if `PETSC_TRUE` the data will be written to disk with single precision
2119a0502c6SHåkon Strandenes 
212811af0c4SBarry Smith   Options Database Key:
2139a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2149a0502c6SHåkon Strandenes 
215811af0c4SBarry Smith   Note:
21695452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2179a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2189a0502c6SHåkon Strandenes 
2199a0502c6SHåkon Strandenes   Level: intermediate
2209a0502c6SHåkon Strandenes 
221811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
222811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
2239a0502c6SHåkon Strandenes @*/
2249371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg) {
2259a0502c6SHåkon Strandenes   PetscFunctionBegin;
2269a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
227cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
2289a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2299a0502c6SHåkon Strandenes }
2309a0502c6SHåkon Strandenes 
231df863907SAlex Fikl /*@
2329a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
233811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2349a0502c6SHåkon Strandenes 
235811af0c4SBarry Smith     Logically Collective on viewer
2369a0502c6SHåkon Strandenes 
2379a0502c6SHåkon Strandenes   Input Parameter:
238811af0c4SBarry Smith .  viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`
2399a0502c6SHåkon Strandenes 
2409a0502c6SHåkon Strandenes   Output Parameter:
241811af0c4SBarry Smith .  flg - if `PETSC_TRUE` the data will be written to disk with single precision
2429a0502c6SHåkon Strandenes 
2439a0502c6SHåkon Strandenes   Level: intermediate
2449a0502c6SHåkon Strandenes 
245db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
246811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
2479a0502c6SHåkon Strandenes @*/
2489371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg) {
2499a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
2509a0502c6SHåkon Strandenes 
2519a0502c6SHåkon Strandenes   PetscFunctionBegin;
2529a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
2539a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2549a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2559a0502c6SHåkon Strandenes }
2569a0502c6SHåkon Strandenes 
2579371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) {
258ee8b9732SVaclav Hapla   PetscFunctionBegin;
259ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
260ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
261c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
262ee8b9732SVaclav Hapla   {
263ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
264792fecdfSBarry Smith     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
265ee8b9732SVaclav Hapla   }
266ee8b9732SVaclav Hapla #else
2679566063dSJacob 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"));
268ee8b9732SVaclav Hapla #endif
269ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
270ee8b9732SVaclav Hapla }
271ee8b9732SVaclav Hapla 
272ee8b9732SVaclav Hapla /*@
273ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
274ee8b9732SVaclav Hapla 
275ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
276ee8b9732SVaclav Hapla 
277ee8b9732SVaclav Hapla   Input Parameters:
278811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
279811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)
280ee8b9732SVaclav Hapla 
281811af0c4SBarry Smith   Options Database Key:
282ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
283ee8b9732SVaclav Hapla 
284811af0c4SBarry Smith   Note:
285ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
28653effed7SBarry 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.
287ee8b9732SVaclav Hapla 
288811af0c4SBarry Smith   Developer Note:
289811af0c4SBarry Smith   In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively.
290ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
291ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
292ee8b9732SVaclav Hapla 
293ee8b9732SVaclav Hapla   Level: intermediate
294ee8b9732SVaclav Hapla 
295811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
296ee8b9732SVaclav Hapla @*/
2979371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg) {
298ee8b9732SVaclav Hapla   PetscFunctionBegin;
299ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
300ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer, flg, 2);
301cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
302ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
303ee8b9732SVaclav Hapla }
304ee8b9732SVaclav Hapla 
3059371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) {
306c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
307ee8b9732SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
308ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
309c796909eSBarry Smith #endif
310ee8b9732SVaclav Hapla 
311ee8b9732SVaclav Hapla   PetscFunctionBegin;
312c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
313c796909eSBarry Smith   *flg = PETSC_FALSE;
314c796909eSBarry Smith #else
315792fecdfSBarry Smith   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
316ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
317c796909eSBarry Smith #endif
318ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
319ee8b9732SVaclav Hapla }
320ee8b9732SVaclav Hapla 
321ee8b9732SVaclav Hapla /*@
322ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
323ee8b9732SVaclav Hapla 
324ee8b9732SVaclav Hapla   Not Collective
325ee8b9732SVaclav Hapla 
326ee8b9732SVaclav Hapla   Input Parameters:
327811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer`
328ee8b9732SVaclav Hapla 
329ee8b9732SVaclav Hapla   Output Parameters:
330ee8b9732SVaclav Hapla . flg - the flag
331ee8b9732SVaclav Hapla 
332ee8b9732SVaclav Hapla   Level: intermediate
333ee8b9732SVaclav Hapla 
334811af0c4SBarry Smith   Note:
335811af0c4SBarry 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.
336811af0c4SBarry Smith   For more details, see `PetscViewerHDF5SetCollective()`.
337ee8b9732SVaclav Hapla 
338811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
339ee8b9732SVaclav Hapla @*/
3409371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg) {
341ee8b9732SVaclav Hapla   PetscFunctionBegin;
342ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
343534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
344ee8b9732SVaclav Hapla 
345cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
346ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
347ee8b9732SVaclav Hapla }
348ee8b9732SVaclav Hapla 
3499371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) {
3505c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
3515c6c1daeSBarry Smith   hid_t             plist_id;
3525c6c1daeSBarry Smith 
3535c6c1daeSBarry Smith   PetscFunctionBegin;
354792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
3559566063dSJacob Faibussowitsch   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
3569566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &hdf5->filename));
3575c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
358792fecdfSBarry Smith   PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS));
359c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
360792fecdfSBarry Smith   PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
361c796909eSBarry Smith #endif
3625c6c1daeSBarry Smith   /* Create or open the file collectively */
3635c6c1daeSBarry Smith   switch (hdf5->btype) {
3645c6c1daeSBarry Smith   case FILE_MODE_READ:
3658a2871f6SBarry Smith     if (PetscDefined(USE_DEBUG)) {
3668a2871f6SBarry Smith       PetscMPIInt rank;
3678a2871f6SBarry Smith       PetscBool   flg;
3688a2871f6SBarry Smith 
3699566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
3708a2871f6SBarry Smith       if (rank == 0) {
3719566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
3725f80ce2aSJacob Faibussowitsch         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
3738a2871f6SBarry Smith       }
3749566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
3758a2871f6SBarry Smith     }
376792fecdfSBarry Smith     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id));
3775c6c1daeSBarry Smith     break;
3785c6c1daeSBarry Smith   case FILE_MODE_APPEND:
3799371c9d4SSatish Balay   case FILE_MODE_UPDATE: {
3807e4fd573SVaclav Hapla     PetscBool flg;
3819566063dSJacob Faibussowitsch     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
382792fecdfSBarry Smith     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id));
383792fecdfSBarry Smith     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
3845c6c1daeSBarry Smith     break;
3857e4fd573SVaclav Hapla   }
3869371c9d4SSatish Balay   case FILE_MODE_WRITE: PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); break;
3879371c9d4SSatish Balay   case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
3889371c9d4SSatish Balay   default: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
3895c6c1daeSBarry Smith   }
3905f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
391792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (plist_id));
3925c6c1daeSBarry Smith   PetscFunctionReturn(0);
3935c6c1daeSBarry Smith }
3945c6c1daeSBarry Smith 
3959371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name) {
396d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
397d1232d7fSToby Isaac 
398d1232d7fSToby Isaac   PetscFunctionBegin;
399d1232d7fSToby Isaac   *name = vhdf5->filename;
400d1232d7fSToby Isaac   PetscFunctionReturn(0);
401d1232d7fSToby Isaac }
402d1232d7fSToby Isaac 
4039371c9d4SSatish Balay static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) {
4045dc64a97SVaclav Hapla   /*
405b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4065dc64a97SVaclav Hapla   */
407b723ab35SVaclav Hapla 
408b723ab35SVaclav Hapla   PetscFunctionBegin;
409b723ab35SVaclav Hapla   PetscFunctionReturn(0);
410b723ab35SVaclav Hapla }
411b723ab35SVaclav Hapla 
4129371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) {
41319a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
41419a20e4cSMatthew G. Knepley 
41519a20e4cSMatthew G. Knepley   PetscFunctionBegin;
41619a20e4cSMatthew G. Knepley   hdf5->defTimestepping = flg;
41719a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
41819a20e4cSMatthew G. Knepley }
41919a20e4cSMatthew G. Knepley 
4209371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) {
42119a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
42219a20e4cSMatthew G. Knepley 
42319a20e4cSMatthew G. Knepley   PetscFunctionBegin;
42419a20e4cSMatthew G. Knepley   *flg = hdf5->defTimestepping;
42519a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
42619a20e4cSMatthew G. Knepley }
42719a20e4cSMatthew G. Knepley 
42819a20e4cSMatthew G. Knepley /*@
42919a20e4cSMatthew G. Knepley   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
43019a20e4cSMatthew G. Knepley 
431811af0c4SBarry Smith   Logically Collective on viewer
43219a20e4cSMatthew G. Knepley 
43319a20e4cSMatthew G. Knepley   Input Parameters:
434811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
435811af0c4SBarry Smith - flg    - if `PETSC_TRUE` we will assume that timestepping is on
43619a20e4cSMatthew G. Knepley 
437811af0c4SBarry Smith   Options Database Key:
43819a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
43919a20e4cSMatthew G. Knepley 
440811af0c4SBarry Smith   Note:
44119a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
44219a20e4cSMatthew G. Knepley 
44319a20e4cSMatthew G. Knepley   Level: intermediate
44419a20e4cSMatthew G. Knepley 
445811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
44619a20e4cSMatthew G. Knepley @*/
4479371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) {
44819a20e4cSMatthew G. Knepley   PetscFunctionBegin;
44919a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
450cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
45119a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
45219a20e4cSMatthew G. Knepley }
45319a20e4cSMatthew G. Knepley 
45419a20e4cSMatthew G. Knepley /*@
45519a20e4cSMatthew G. Knepley   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
45619a20e4cSMatthew G. Knepley 
45719a20e4cSMatthew G. Knepley   Not collective
45819a20e4cSMatthew G. Knepley 
45919a20e4cSMatthew G. Knepley   Input Parameter:
460811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
46119a20e4cSMatthew G. Knepley 
46219a20e4cSMatthew G. Knepley   Output Parameter:
463811af0c4SBarry Smith . flg    - if `PETSC_TRUE` we will assume that timestepping is on
46419a20e4cSMatthew G. Knepley 
46519a20e4cSMatthew G. Knepley   Level: intermediate
46619a20e4cSMatthew G. Knepley 
467811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
46819a20e4cSMatthew G. Knepley @*/
4699371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) {
47019a20e4cSMatthew G. Knepley   PetscFunctionBegin;
47119a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
472cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
47319a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
47419a20e4cSMatthew G. Knepley }
47519a20e4cSMatthew G. Knepley 
4768556b5ebSBarry Smith /*MC
4778556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
4788556b5ebSBarry Smith 
479811af0c4SBarry Smith   Level: beginner
480811af0c4SBarry Smith 
481db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
482db781477SPatrick Sanan           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
483db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
484db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
4858556b5ebSBarry Smith M*/
486d1232d7fSToby Isaac 
4879371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) {
4885c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
4895c6c1daeSBarry Smith 
4905c6c1daeSBarry Smith   PetscFunctionBegin;
49199335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
49299335c34SVaclav Hapla   {
49399335c34SVaclav Hapla     PetscMPIInt size;
4949566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
4955f80ce2aSJacob 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)");
49699335c34SVaclav Hapla   }
49799335c34SVaclav Hapla #endif
49899335c34SVaclav Hapla 
499*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hdf5));
5005c6c1daeSBarry Smith 
5015c6c1daeSBarry Smith   v->data                = (void *)hdf5;
5025c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
50382ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
504b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
5051b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
5066226335aSVaclav Hapla   v->ops->flush          = PetscViewerFlush_HDF5;
5077e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
508908793a3SLisandro Dalcin   hdf5->filename         = NULL;
5095c6c1daeSBarry Smith   hdf5->timestep         = -1;
5100298fd71SBarry Smith   hdf5->groups           = NULL;
5115c6c1daeSBarry Smith 
512792fecdfSBarry Smith   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
5139c5072fbSVaclav Hapla 
5149566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
5159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
5169566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
5179566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
5189566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
5199566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
5219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
5229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
5239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
5245c6c1daeSBarry Smith   PetscFunctionReturn(0);
5255c6c1daeSBarry Smith }
5265c6c1daeSBarry Smith 
5275c6c1daeSBarry Smith /*@C
528811af0c4SBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`
5295c6c1daeSBarry Smith 
530d083f849SBarry Smith    Collective
5315c6c1daeSBarry Smith 
5325c6c1daeSBarry Smith    Input Parameters:
5335c6c1daeSBarry Smith +  comm - MPI communicator
5345c6c1daeSBarry Smith .  name - name of file
5355c6c1daeSBarry Smith -  type - type of file
5365c6c1daeSBarry Smith 
5375c6c1daeSBarry Smith    Output Parameter:
538811af0c4SBarry Smith .  hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
5395c6c1daeSBarry Smith 
540811af0c4SBarry Smith   Options Database Keys:
541a2b725a8SWilliam 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
542a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
54382ea9b62SBarry Smith 
5445c6c1daeSBarry Smith    Level: beginner
5455c6c1daeSBarry Smith 
5467e4fd573SVaclav Hapla    Notes:
5477e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
548811af0c4SBarry Smith .vb
549811af0c4SBarry Smith   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
550811af0c4SBarry Smith   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
551811af0c4SBarry Smith   FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL]
552811af0c4SBarry Smith   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
553811af0c4SBarry Smith .ve
5547e4fd573SVaclav Hapla 
555811af0c4SBarry Smith    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.
5567e4fd573SVaclav Hapla 
5575c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
5585c6c1daeSBarry Smith 
559811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
560db781477SPatrick Sanan           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
561db781477SPatrick Sanan           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
5625c6c1daeSBarry Smith @*/
5639371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) {
5645c6c1daeSBarry Smith   PetscFunctionBegin;
5659566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, hdf5v));
5669566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
5679566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
5689566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*hdf5v, name));
5699566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*hdf5v));
5705c6c1daeSBarry Smith   PetscFunctionReturn(0);
5715c6c1daeSBarry Smith }
5725c6c1daeSBarry Smith 
5735c6c1daeSBarry Smith /*@C
5745c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
5755c6c1daeSBarry Smith 
5765c6c1daeSBarry Smith   Not collective
5775c6c1daeSBarry Smith 
5785c6c1daeSBarry Smith   Input Parameter:
579811af0c4SBarry Smith . viewer - the `PetscViewer`
5805c6c1daeSBarry Smith 
5815c6c1daeSBarry Smith   Output Parameter:
5825c6c1daeSBarry Smith . file_id - The file id
5835c6c1daeSBarry Smith 
5845c6c1daeSBarry Smith   Level: intermediate
5855c6c1daeSBarry Smith 
586811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
5875c6c1daeSBarry Smith @*/
5889371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) {
5895c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
5905c6c1daeSBarry Smith 
5915c6c1daeSBarry Smith   PetscFunctionBegin;
5925c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
5935c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
5945c6c1daeSBarry Smith   PetscFunctionReturn(0);
5955c6c1daeSBarry Smith }
5965c6c1daeSBarry Smith 
5975c6c1daeSBarry Smith /*@C
5985c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
5995c6c1daeSBarry Smith 
6005c6c1daeSBarry Smith   Not collective
6015c6c1daeSBarry Smith 
6025c6c1daeSBarry Smith   Input Parameters:
603811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
6045c6c1daeSBarry Smith - name - The group name
6055c6c1daeSBarry Smith 
6065c6c1daeSBarry Smith   Level: intermediate
6075c6c1daeSBarry Smith 
6082e361344SVaclav Hapla   Notes:
6092e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
610811af0c4SBarry Smith .vb
611811af0c4SBarry Smith   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.
612811af0c4SBarry Smith   NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
613811af0c4SBarry Smith   "." means the current group is pushed again.
614811af0c4SBarry Smith .ve
6152e361344SVaclav Hapla 
6162e361344SVaclav Hapla   Example:
6172e361344SVaclav Hapla   Suppose the current group is "/a".
6182e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
6192e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
6202e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
6212e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
6222e361344SVaclav Hapla 
623811af0c4SBarry Smith   Developer Note:
6242e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
625820fc2d1SVaclav Hapla 
626811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
6275c6c1daeSBarry Smith @*/
6289371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) {
6295c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6307d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6312e361344SVaclav Hapla   size_t                    i, len;
6322e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
6332e361344SVaclav Hapla   const char               *gname;
6345c6c1daeSBarry Smith 
6355c6c1daeSBarry Smith   PetscFunctionBegin;
6365c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
637820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
6389566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
6392e361344SVaclav Hapla   gname = NULL;
6402e361344SVaclav Hapla   if (len) {
6412e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
6422e361344SVaclav Hapla       /* use current name */
6432e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
6442e361344SVaclav Hapla     } else if (name[0] == '/') {
6452e361344SVaclav Hapla       /* absolute */
6462e361344SVaclav Hapla       for (i = 1; i < len; i++) {
6472e361344SVaclav Hapla         if (name[i] != '/') {
6482e361344SVaclav Hapla           gname = name;
6492e361344SVaclav Hapla           break;
6502e361344SVaclav Hapla         }
6512e361344SVaclav Hapla       }
6522e361344SVaclav Hapla     } else {
6532e361344SVaclav Hapla       /* relative */
6542e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
6559566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
6562e361344SVaclav Hapla       gname = buf;
6572e361344SVaclav Hapla     }
6582e361344SVaclav Hapla   }
6599566063dSJacob Faibussowitsch   PetscCall(PetscNew(&groupNode));
6609566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
6615c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
6625c6c1daeSBarry Smith   hdf5->groups    = groupNode;
6635c6c1daeSBarry Smith   PetscFunctionReturn(0);
6645c6c1daeSBarry Smith }
6655c6c1daeSBarry Smith 
6663ef9c667SSatish Balay /*@
6675c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
6685c6c1daeSBarry Smith 
6695c6c1daeSBarry Smith   Not collective
6705c6c1daeSBarry Smith 
6715c6c1daeSBarry Smith   Input Parameter:
672811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
6735c6c1daeSBarry Smith 
6745c6c1daeSBarry Smith   Level: intermediate
6755c6c1daeSBarry Smith 
676811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
6775c6c1daeSBarry Smith @*/
6789371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) {
6795c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6807d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6815c6c1daeSBarry Smith 
6825c6c1daeSBarry Smith   PetscFunctionBegin;
6835c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
6845f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
6855c6c1daeSBarry Smith   groupNode    = hdf5->groups;
6865c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
6879566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
6889566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
6895c6c1daeSBarry Smith   PetscFunctionReturn(0);
6905c6c1daeSBarry Smith }
6915c6c1daeSBarry Smith 
6925c6c1daeSBarry Smith /*@C
693811af0c4SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
694874270d9SVaclav Hapla   If none has been assigned, returns NULL.
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith   Not collective
6975c6c1daeSBarry Smith 
6985c6c1daeSBarry Smith   Input Parameter:
699811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7005c6c1daeSBarry Smith 
7015c6c1daeSBarry Smith   Output Parameter:
7025c6c1daeSBarry Smith . name - The group name
7035c6c1daeSBarry Smith 
7045c6c1daeSBarry Smith   Level: intermediate
7055c6c1daeSBarry Smith 
706811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`
7075c6c1daeSBarry Smith @*/
7089371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) {
7095c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7105c6c1daeSBarry Smith 
7115c6c1daeSBarry Smith   PetscFunctionBegin;
7125c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
713c959eef4SJed Brown   PetscValidPointer(name, 2);
714a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
7150298fd71SBarry Smith   else *name = NULL;
7165c6c1daeSBarry Smith   PetscFunctionReturn(0);
7175c6c1daeSBarry Smith }
7185c6c1daeSBarry Smith 
7198c557b5aSVaclav Hapla /*@
720811af0c4SBarry Smith   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
721874270d9SVaclav Hapla   and return this group's ID and file ID.
722811af0c4SBarry Smith   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
723874270d9SVaclav Hapla 
724874270d9SVaclav Hapla   Not collective
725874270d9SVaclav Hapla 
726874270d9SVaclav Hapla   Input Parameter:
727811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
728874270d9SVaclav Hapla 
729d8d19677SJose E. Roman   Output Parameters:
730874270d9SVaclav Hapla + fileId - The HDF5 file ID
731874270d9SVaclav Hapla - groupId - The HDF5 group ID
732874270d9SVaclav Hapla 
733811af0c4SBarry Smith   Note:
734e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
735e74616beSVaclav Hapla 
736874270d9SVaclav Hapla   Level: intermediate
737874270d9SVaclav Hapla 
738811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
739874270d9SVaclav Hapla @*/
7409371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) {
74190e03537SVaclav Hapla   hid_t       file_id;
74290e03537SVaclav Hapla   H5O_type_t  type;
74376d59af2SVaclav Hapla   const char *groupName = NULL, *fileName = NULL;
74476d59af2SVaclav Hapla   PetscBool   writable, has;
74554dbf706SVaclav Hapla 
74654dbf706SVaclav Hapla   PetscFunctionBegin;
7479566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
7489566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
7499566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
7509566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName));
7519566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
75276d59af2SVaclav Hapla   if (!has) {
7535f80ce2aSJacob 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);
754f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
75576d59af2SVaclav Hapla   }
7565f80ce2aSJacob 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);
757792fecdfSBarry Smith   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName ? groupName : "/", H5P_DEFAULT));
75854dbf706SVaclav Hapla   *fileId = file_id;
75954dbf706SVaclav Hapla   PetscFunctionReturn(0);
76054dbf706SVaclav Hapla }
76154dbf706SVaclav Hapla 
7623ef9c667SSatish Balay /*@
763d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
7645c6c1daeSBarry Smith 
7655c6c1daeSBarry Smith   Not collective
7665c6c1daeSBarry Smith 
7675c6c1daeSBarry Smith   Input Parameter:
768811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7695c6c1daeSBarry Smith 
7705c6c1daeSBarry Smith   Level: intermediate
7715c6c1daeSBarry Smith 
772d7dd068bSVaclav Hapla   Notes:
773811af0c4SBarry Smith   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
774811af0c4SBarry Smith   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
775811af0c4SBarry Smith   Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. `VecView()`].
776811af0c4SBarry Smith   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
777811af0c4SBarry Smith   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
778d7dd068bSVaclav Hapla 
779d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
780d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
781d7dd068bSVaclav Hapla 
782811af0c4SBarry Smith   Developer note:
783d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
784d7dd068bSVaclav Hapla 
785811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
786d7dd068bSVaclav Hapla @*/
7879371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) {
788d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
789d7dd068bSVaclav Hapla 
790d7dd068bSVaclav Hapla   PetscFunctionBegin;
791d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7925f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
793d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
794d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
795d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
796d7dd068bSVaclav Hapla }
797d7dd068bSVaclav Hapla 
798d7dd068bSVaclav Hapla /*@
799d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
800d7dd068bSVaclav Hapla 
801d7dd068bSVaclav Hapla   Not collective
802d7dd068bSVaclav Hapla 
803d7dd068bSVaclav Hapla   Input Parameter:
804811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
805d7dd068bSVaclav Hapla 
806d7dd068bSVaclav Hapla   Level: intermediate
807d7dd068bSVaclav Hapla 
808811af0c4SBarry Smith   Note:
809811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
810d7dd068bSVaclav Hapla 
811811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
812d7dd068bSVaclav Hapla @*/
8139371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) {
814d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
815d7dd068bSVaclav Hapla 
816d7dd068bSVaclav Hapla   PetscFunctionBegin;
817d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8185f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
819d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
820d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
821d7dd068bSVaclav Hapla }
822d7dd068bSVaclav Hapla 
823d7dd068bSVaclav Hapla /*@
824d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
825d7dd068bSVaclav Hapla 
826d7dd068bSVaclav Hapla   Not collective
827d7dd068bSVaclav Hapla 
828d7dd068bSVaclav Hapla   Input Parameter:
829811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
830d7dd068bSVaclav Hapla 
831d7dd068bSVaclav Hapla   Output Parameter:
832d7dd068bSVaclav Hapla . flg - is timestepping active?
833d7dd068bSVaclav Hapla 
834d7dd068bSVaclav Hapla   Level: intermediate
835d7dd068bSVaclav Hapla 
836811af0c4SBarry Smith   Note:
837811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
838d7dd068bSVaclav Hapla 
839811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
840d7dd068bSVaclav Hapla @*/
8419371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) {
842d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
843d7dd068bSVaclav Hapla 
844d7dd068bSVaclav Hapla   PetscFunctionBegin;
845d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
846d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
847d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
848d7dd068bSVaclav Hapla }
849d7dd068bSVaclav Hapla 
850d7dd068bSVaclav Hapla /*@
851d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
852d7dd068bSVaclav Hapla 
853d7dd068bSVaclav Hapla   Not collective
854d7dd068bSVaclav Hapla 
855d7dd068bSVaclav Hapla   Input Parameter:
856811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
857d7dd068bSVaclav Hapla 
858d7dd068bSVaclav Hapla   Level: intermediate
859d7dd068bSVaclav Hapla 
860811af0c4SBarry Smith   Note:
861811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
862d7dd068bSVaclav Hapla 
863811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
8645c6c1daeSBarry Smith @*/
8659371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) {
8665c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
8675c6c1daeSBarry Smith 
8685c6c1daeSBarry Smith   PetscFunctionBegin;
8695c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8705f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
8715c6c1daeSBarry Smith   ++hdf5->timestep;
8725c6c1daeSBarry Smith   PetscFunctionReturn(0);
8735c6c1daeSBarry Smith }
8745c6c1daeSBarry Smith 
8753ef9c667SSatish Balay /*@
876d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
8775c6c1daeSBarry Smith 
878d7dd068bSVaclav Hapla   Logically collective
8795c6c1daeSBarry Smith 
8805c6c1daeSBarry Smith   Input Parameters:
881811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
882d7dd068bSVaclav Hapla - timestep - The timestep
8835c6c1daeSBarry Smith 
8845c6c1daeSBarry Smith   Level: intermediate
8855c6c1daeSBarry Smith 
886811af0c4SBarry Smith   Note:
887811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
888d7dd068bSVaclav Hapla 
889811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
8905c6c1daeSBarry Smith @*/
8919371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) {
8925c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
8935c6c1daeSBarry Smith 
8945c6c1daeSBarry Smith   PetscFunctionBegin;
8955c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
896d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
8975f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
8985f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
8995c6c1daeSBarry Smith   hdf5->timestep = timestep;
9005c6c1daeSBarry Smith   PetscFunctionReturn(0);
9015c6c1daeSBarry Smith }
9025c6c1daeSBarry Smith 
9033ef9c667SSatish Balay /*@
9045c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
9055c6c1daeSBarry Smith 
9065c6c1daeSBarry Smith   Not collective
9075c6c1daeSBarry Smith 
9085c6c1daeSBarry Smith   Input Parameter:
909811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
9105c6c1daeSBarry Smith 
9115c6c1daeSBarry Smith   Output Parameter:
912d7dd068bSVaclav Hapla . timestep - The timestep
9135c6c1daeSBarry Smith 
9145c6c1daeSBarry Smith   Level: intermediate
9155c6c1daeSBarry Smith 
916811af0c4SBarry Smith   Note:
917811af0c4SBarry Smith   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
918d7dd068bSVaclav Hapla 
919811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
9205c6c1daeSBarry Smith @*/
9219371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) {
9225c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9235c6c1daeSBarry Smith 
9245c6c1daeSBarry Smith   PetscFunctionBegin;
9255c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
926d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep, 2);
9275f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9285c6c1daeSBarry Smith   *timestep = hdf5->timestep;
9295c6c1daeSBarry Smith   PetscFunctionReturn(0);
9305c6c1daeSBarry Smith }
9315c6c1daeSBarry Smith 
93236bce6e8SMatthew G. Knepley /*@C
93336bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
93436bce6e8SMatthew G. Knepley 
93536bce6e8SMatthew G. Knepley   Not collective
93636bce6e8SMatthew G. Knepley 
93736bce6e8SMatthew G. Knepley   Input Parameter:
938811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
93936bce6e8SMatthew G. Knepley 
94036bce6e8SMatthew G. Knepley   Output Parameter:
941811af0c4SBarry Smith . mtype - the HDF5  datatype
94236bce6e8SMatthew G. Knepley 
94336bce6e8SMatthew G. Knepley   Level: advanced
94436bce6e8SMatthew G. Knepley 
945db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
94636bce6e8SMatthew G. Knepley @*/
9479371c9d4SSatish Balay PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) {
94836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
94936bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
95036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
95136bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_LLONG;
95236bce6e8SMatthew G. Knepley #else
95336bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_INT;
95436bce6e8SMatthew G. Knepley #endif
95536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
95636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
95736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
958c9450970SBarry Smith   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
959de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
96036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
96136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
96236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
9637e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
96436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
96536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
96636bce6e8SMatthew G. Knepley }
96736bce6e8SMatthew G. Knepley 
96836bce6e8SMatthew G. Knepley /*@C
96936bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
97036bce6e8SMatthew G. Knepley 
97136bce6e8SMatthew G. Knepley   Not collective
97236bce6e8SMatthew G. Knepley 
97336bce6e8SMatthew G. Knepley   Input Parameter:
974811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
97536bce6e8SMatthew G. Knepley 
97636bce6e8SMatthew G. Knepley   Output Parameter:
977811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
97836bce6e8SMatthew G. Knepley 
97936bce6e8SMatthew G. Knepley   Level: advanced
98036bce6e8SMatthew G. Knepley 
981db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
98236bce6e8SMatthew G. Knepley @*/
9839371c9d4SSatish Balay PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) {
98436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
98536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
98636bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
98736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
98836bce6e8SMatthew G. Knepley #else
98936bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
99036bce6e8SMatthew G. Knepley #endif
99136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
99236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
99336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
99436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
99536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
99636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
9977e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
99836bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
99936bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
100036bce6e8SMatthew G. Knepley }
100136bce6e8SMatthew G. Knepley 
1002df863907SAlex Fikl /*@C
1003b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
100436bce6e8SMatthew G. Knepley 
1005343a20b2SVaclav Hapla   Collective
1006343a20b2SVaclav Hapla 
100736bce6e8SMatthew G. Knepley   Input Parameters:
1008811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
10094302210dSVaclav Hapla . parent - The parent dataset/group name
101036bce6e8SMatthew G. Knepley . name   - The attribute name
101136bce6e8SMatthew G. Knepley . datatype - The attribute type
101236bce6e8SMatthew G. Knepley - value    - The attribute value
101336bce6e8SMatthew G. Knepley 
101436bce6e8SMatthew G. Knepley   Level: advanced
101536bce6e8SMatthew G. Knepley 
1016811af0c4SBarry Smith   Note:
1017343a20b2SVaclav 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.
10184302210dSVaclav Hapla 
1019811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1020811af0c4SBarry Smith           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
102136bce6e8SMatthew G. Knepley @*/
10229371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) {
10234302210dSVaclav Hapla   char     *parentAbsPath;
102460568592SMatthew G. Knepley   hid_t     h5, dataspace, obj, attribute, dtype;
1025080f144cSVaclav Hapla   PetscBool has;
102636bce6e8SMatthew G. Knepley 
102736bce6e8SMatthew G. Knepley   PetscFunctionBegin;
10285cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
10294302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1030c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
10314302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1032b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
10339566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
10349566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
10359566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
10369566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
10377e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
10387e97c476SMatthew G. Knepley     size_t len;
10399566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *)value, &len));
1040792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
10417e97c476SMatthew G. Knepley   }
10429566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1043792fecdfSBarry Smith   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1044792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1045080f144cSVaclav Hapla   if (has) {
1046792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1047080f144cSVaclav Hapla   } else {
1048792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1049080f144cSVaclav Hapla   }
1050792fecdfSBarry Smith   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1051792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1052792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1053792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
1054792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (dataspace));
10559566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
105636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
105736bce6e8SMatthew G. Knepley }
105836bce6e8SMatthew G. Knepley 
1059df863907SAlex Fikl /*@C
1060811af0c4SBarry Smith  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1061577e0e04SVaclav Hapla 
1062343a20b2SVaclav Hapla   Collective
1063343a20b2SVaclav Hapla 
1064577e0e04SVaclav Hapla   Input Parameters:
1065811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1066577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1067577e0e04SVaclav Hapla . name     - The attribute name
1068577e0e04SVaclav Hapla . datatype - The attribute type
1069577e0e04SVaclav Hapla - value    - The attribute value
1070577e0e04SVaclav Hapla 
1071811af0c4SBarry Smith   Note:
1072343a20b2SVaclav 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).
1073811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1074577e0e04SVaclav Hapla 
1075577e0e04SVaclav Hapla   Level: advanced
1076577e0e04SVaclav Hapla 
1077811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1078811af0c4SBarry Smith           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1079577e0e04SVaclav Hapla @*/
10809371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) {
1081577e0e04SVaclav Hapla   PetscFunctionBegin;
1082577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1083577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1084577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1085b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
10869566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
10879566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1088577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1089577e0e04SVaclav Hapla }
1090577e0e04SVaclav Hapla 
1091577e0e04SVaclav Hapla /*@C
1092b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1093ad0c61feSMatthew G. Knepley 
1094343a20b2SVaclav Hapla   Collective
1095343a20b2SVaclav Hapla 
1096ad0c61feSMatthew G. Knepley   Input Parameters:
1097811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
10984302210dSVaclav Hapla . parent - The parent dataset/group name
1099ad0c61feSMatthew G. Knepley . name   - The attribute name
1100a2d6be1bSVaclav Hapla . datatype - The attribute type
1101a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1102ad0c61feSMatthew G. Knepley 
1103ad0c61feSMatthew G. Knepley   Output Parameter:
1104a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1105ad0c61feSMatthew G. Knepley 
1106a2d6be1bSVaclav Hapla   Notes:
1107a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1108a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1109a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1110811af0c4SBarry Smith $  flg = `PETSC_FALSE`;
1111811af0c4SBarry Smith $  PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1112a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1113a2d6be1bSVaclav Hapla 
1114811af0c4SBarry Smith   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1115708d4cb3SBarry Smith 
1116343a20b2SVaclav 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.
1117ad0c61feSMatthew G. Knepley 
1118343a20b2SVaclav Hapla   Level: advanced
11194302210dSVaclav Hapla 
1120811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1121ad0c61feSMatthew G. Knepley @*/
11229371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) {
11234302210dSVaclav Hapla   char     *parentAbsPath;
1124a2d6be1bSVaclav Hapla   hid_t     h5, obj, attribute, dtype;
1125969748fdSVaclav Hapla   PetscBool has;
1126ad0c61feSMatthew G. Knepley 
1127ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
11285cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
11294302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1130c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1131a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1132a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
11339566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11349566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
11359566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
11369566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1137a2d6be1bSVaclav Hapla   if (!has) {
1138a2d6be1bSVaclav Hapla     if (defaultValue) {
1139a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1140a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
11419566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1142a2d6be1bSVaclav Hapla         } else {
1143a2d6be1bSVaclav Hapla           size_t len;
1144792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
11459566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1146a2d6be1bSVaclav Hapla         }
1147a2d6be1bSVaclav Hapla       }
11489566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
1149a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
115098921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1151a2d6be1bSVaclav Hapla   }
11529566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1153792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1154792fecdfSBarry Smith   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1155f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1156f0b58479SMatthew G. Knepley     size_t len;
1157a2d6be1bSVaclav Hapla     hid_t  atype;
1158792fecdfSBarry Smith     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1159792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
11609566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1161792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1162792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1163708d4cb3SBarry Smith   } else {
1164792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1165708d4cb3SBarry Smith   }
1166792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1167e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1168792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
11699566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
1170ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1171ad0c61feSMatthew G. Knepley }
1172ad0c61feSMatthew G. Knepley 
1173577e0e04SVaclav Hapla /*@C
1174811af0c4SBarry Smith  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1175577e0e04SVaclav Hapla 
1176343a20b2SVaclav Hapla   Collective
1177343a20b2SVaclav Hapla 
1178577e0e04SVaclav Hapla   Input Parameters:
1179811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1180577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1181577e0e04SVaclav Hapla . name     - The attribute name
1182577e0e04SVaclav Hapla - datatype - The attribute type
1183577e0e04SVaclav Hapla 
1184577e0e04SVaclav Hapla   Output Parameter:
1185577e0e04SVaclav Hapla . value    - The attribute value
1186577e0e04SVaclav Hapla 
1187811af0c4SBarry Smith   Note:
1188577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1189811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1190577e0e04SVaclav Hapla 
1191577e0e04SVaclav Hapla   Level: advanced
1192577e0e04SVaclav Hapla 
1193811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1194577e0e04SVaclav Hapla @*/
11959371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) {
1196577e0e04SVaclav Hapla   PetscFunctionBegin;
1197577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1198577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1199577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1200064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
12019566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12029566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1203577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1204577e0e04SVaclav Hapla }
1205577e0e04SVaclav Hapla 
12069371c9d4SSatish Balay static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) {
1207a75c3fd4SVaclav Hapla   htri_t exists;
1208a75c3fd4SVaclav Hapla   hid_t  group;
1209a75c3fd4SVaclav Hapla 
1210a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1211792fecdfSBarry Smith   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1212792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1213a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1214792fecdfSBarry Smith     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1215792fecdfSBarry Smith     PetscCallHDF5(H5Gclose, (group));
1216a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1217a75c3fd4SVaclav Hapla   }
1218a75c3fd4SVaclav Hapla   *exists_ = (PetscBool)exists;
1219a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1220a75c3fd4SVaclav Hapla }
1221a75c3fd4SVaclav Hapla 
12229371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) {
122390e03537SVaclav Hapla   const char rootGroupName[] = "/";
12245cdeb108SMatthew G. Knepley   hid_t      h5;
1225e5bf9ebcSVaclav Hapla   PetscBool  exists = PETSC_FALSE;
122615b861d2SVaclav Hapla   PetscInt   i;
122715b861d2SVaclav Hapla   int        n;
122885688ae2SVaclav Hapla   char     **hierarchy;
122985688ae2SVaclav Hapla   char       buf[PETSC_MAX_PATH_LEN] = "";
12305cdeb108SMatthew G. Knepley 
12315cdeb108SMatthew G. Knepley   PetscFunctionBegin;
12325cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
123390e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
123490e03537SVaclav Hapla   else name = rootGroupName;
1235ccd66a83SVaclav Hapla   if (has) {
1236064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
12375cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1238ccd66a83SVaclav Hapla   }
1239ccd66a83SVaclav Hapla   if (otype) {
1240064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
124156cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1242ccd66a83SVaclav Hapla   }
12439566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
124485688ae2SVaclav Hapla 
124585688ae2SVaclav Hapla   /*
124685688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
124785688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
124885688ae2SVaclav Hapla      1) whether it's a valid link
124985688ae2SVaclav Hapla      2) whether this link resolves to an object
125085688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
125185688ae2SVaclav Hapla   */
12529566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
125385688ae2SVaclav Hapla   if (!n) {
125485688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1255ccd66a83SVaclav Hapla     if (has) *has = PETSC_TRUE;
1256ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
12579566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
125885688ae2SVaclav Hapla     PetscFunctionReturn(0);
125985688ae2SVaclav Hapla   }
126085688ae2SVaclav Hapla   for (i = 0; i < n; i++) {
12619566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
12629566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, hierarchy[i]));
12639566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1264a75c3fd4SVaclav Hapla     if (!exists) break;
126585688ae2SVaclav Hapla   }
12669566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
126785688ae2SVaclav Hapla 
126885688ae2SVaclav Hapla   /* If the object exists, get its type */
1269ccd66a83SVaclav Hapla   if (exists && otype) {
12705cdeb108SMatthew G. Knepley     H5O_info_t info;
12715cdeb108SMatthew G. Knepley 
127276276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1273792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
127456cc0592SVaclav Hapla     *otype = info.type;
12755cdeb108SMatthew G. Knepley   }
1276ccd66a83SVaclav Hapla   if (has) *has = exists;
12775cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
12785cdeb108SMatthew G. Knepley }
12795cdeb108SMatthew G. Knepley 
12804302210dSVaclav Hapla /*@C
12810a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
12820a9f38efSVaclav Hapla 
1283343a20b2SVaclav Hapla   Collective
1284343a20b2SVaclav Hapla 
12850a9f38efSVaclav Hapla   Input Parameters:
1286811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
12874302210dSVaclav Hapla - path - The group path
12880a9f38efSVaclav Hapla 
12890a9f38efSVaclav Hapla   Output Parameter:
12900a9f38efSVaclav Hapla . has - Flag for group existence
12910a9f38efSVaclav Hapla 
12920a9f38efSVaclav Hapla   Level: advanced
12930a9f38efSVaclav Hapla 
12944302210dSVaclav Hapla   Notes:
1295785c443eSVaclav 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.
1296785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
12974302210dSVaclav Hapla 
1298811af0c4SBarry Smith   If path exists but is not a group, `PETSC_FALSE` is returned.
1299811af0c4SBarry Smith 
1300811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
13010a9f38efSVaclav Hapla @*/
13029371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) {
13030a9f38efSVaclav Hapla   H5O_type_t type;
13044302210dSVaclav Hapla   char      *abspath;
13050a9f38efSVaclav Hapla 
13060a9f38efSVaclav Hapla   PetscFunctionBegin;
13070a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
13084302210dSVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
13094302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
13109566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
13119566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
13124302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
13139566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
13140a9f38efSVaclav Hapla   PetscFunctionReturn(0);
13150a9f38efSVaclav Hapla }
13160a9f38efSVaclav Hapla 
131789e0ef10SVaclav Hapla /*@C
131889e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
131989e0ef10SVaclav Hapla 
1320343a20b2SVaclav Hapla   Collective
1321343a20b2SVaclav Hapla 
132289e0ef10SVaclav Hapla   Input Parameters:
1323811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
132489e0ef10SVaclav Hapla - path - The dataset path
132589e0ef10SVaclav Hapla 
132689e0ef10SVaclav Hapla   Output Parameter:
132789e0ef10SVaclav Hapla . has - Flag whether dataset exists
132889e0ef10SVaclav Hapla 
132989e0ef10SVaclav Hapla   Level: advanced
133089e0ef10SVaclav Hapla 
133189e0ef10SVaclav Hapla   Notes:
1332343a20b2SVaclav 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.
133389e0ef10SVaclav Hapla 
1334811af0c4SBarry Smith   If path is NULL or empty, has is set to `PETSC_FALSE`.
1335811af0c4SBarry Smith 
1336811af0c4SBarry Smith   If path exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1337811af0c4SBarry Smith 
1338811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
133989e0ef10SVaclav Hapla @*/
13409371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) {
134189e0ef10SVaclav Hapla   H5O_type_t type;
134289e0ef10SVaclav Hapla   char      *abspath;
134389e0ef10SVaclav Hapla 
134489e0ef10SVaclav Hapla   PetscFunctionBegin;
134589e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
134689e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
134789e0ef10SVaclav Hapla   PetscValidBoolPointer(has, 3);
13489566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
13499566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
135089e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
13519566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
135289e0ef10SVaclav Hapla   PetscFunctionReturn(0);
135389e0ef10SVaclav Hapla }
135489e0ef10SVaclav Hapla 
13550a9f38efSVaclav Hapla /*@
1356e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1357ecce1506SVaclav Hapla 
1358343a20b2SVaclav Hapla   Collective
1359343a20b2SVaclav Hapla 
1360ecce1506SVaclav Hapla   Input Parameters:
1361811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1362ecce1506SVaclav Hapla - obj    - The named object
1363ecce1506SVaclav Hapla 
1364ecce1506SVaclav Hapla   Output Parameter:
136589e0ef10SVaclav Hapla . has    - Flag for dataset existence
1366ecce1506SVaclav Hapla 
1367e3f143f7SVaclav Hapla   Notes:
136889e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1369811af0c4SBarry Smith 
1370811af0c4SBarry Smith   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1371e3f143f7SVaclav Hapla 
1372ecce1506SVaclav Hapla   Level: advanced
1373ecce1506SVaclav Hapla 
1374811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1375ecce1506SVaclav Hapla @*/
13769371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) {
137789e0ef10SVaclav Hapla   size_t len;
1378ecce1506SVaclav Hapla 
1379ecce1506SVaclav Hapla   PetscFunctionBegin;
1380c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1381c57a1a9aSVaclav Hapla   PetscValidHeader(obj, 2);
13824302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
13839566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
13845f80ce2aSJacob Faibussowitsch   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
13859566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1386ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1387ecce1506SVaclav Hapla }
1388ecce1506SVaclav Hapla 
1389df863907SAlex Fikl /*@C
1390ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1391e7bdbf8eSMatthew G. Knepley 
1392343a20b2SVaclav Hapla   Collective
1393343a20b2SVaclav Hapla 
1394e7bdbf8eSMatthew G. Knepley   Input Parameters:
1395811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
13964302210dSVaclav Hapla . parent - The parent dataset/group name
1397e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1398e7bdbf8eSMatthew G. Knepley 
1399e7bdbf8eSMatthew G. Knepley   Output Parameter:
1400e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1401e7bdbf8eSMatthew G. Knepley 
1402e7bdbf8eSMatthew G. Knepley   Level: advanced
1403e7bdbf8eSMatthew G. Knepley 
1404811af0c4SBarry Smith   Note:
1405343a20b2SVaclav 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.
14064302210dSVaclav Hapla 
1407811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1408e7bdbf8eSMatthew G. Knepley @*/
14099371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
14104302210dSVaclav Hapla   char *parentAbsPath;
1411e7bdbf8eSMatthew G. Knepley 
1412e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
14135cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14144302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1415c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
14164302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
14179566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
14189566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
14199566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
14209566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
142106db490cSVaclav Hapla   PetscFunctionReturn(0);
142206db490cSVaclav Hapla }
142306db490cSVaclav Hapla 
1424577e0e04SVaclav Hapla /*@C
1425811af0c4SBarry Smith  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1426577e0e04SVaclav Hapla 
1427343a20b2SVaclav Hapla   Collective
1428343a20b2SVaclav Hapla 
1429577e0e04SVaclav Hapla   Input Parameters:
1430811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1431577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1432577e0e04SVaclav Hapla - name   - The attribute name
1433577e0e04SVaclav Hapla 
1434577e0e04SVaclav Hapla   Output Parameter:
1435577e0e04SVaclav Hapla . has    - Flag for attribute existence
1436577e0e04SVaclav Hapla 
1437811af0c4SBarry Smith   Note:
1438577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1439811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1440577e0e04SVaclav Hapla 
1441577e0e04SVaclav Hapla   Level: advanced
1442577e0e04SVaclav Hapla 
1443811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1444577e0e04SVaclav Hapla @*/
14459371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) {
1446577e0e04SVaclav Hapla   PetscFunctionBegin;
1447577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1448577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1449577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
14504302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
14519566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
14529566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1453577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1454577e0e04SVaclav Hapla }
1455577e0e04SVaclav Hapla 
14569371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
145706db490cSVaclav Hapla   hid_t  h5;
145806db490cSVaclav Hapla   htri_t hhas;
145906db490cSVaclav Hapla 
146006db490cSVaclav Hapla   PetscFunctionBegin;
14619566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1462792fecdfSBarry Smith   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
14635cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1464e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1465e7bdbf8eSMatthew G. Knepley }
1466e7bdbf8eSMatthew G. Knepley 
1467a75e6a4aSMatthew G. Knepley /*
1468a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1469a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1470a75e6a4aSMatthew G. Knepley */
1471d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1472a75e6a4aSMatthew G. Knepley 
1473a75e6a4aSMatthew G. Knepley /*@C
1474811af0c4SBarry Smith   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1475a75e6a4aSMatthew G. Knepley 
1476d083f849SBarry Smith   Collective
1477a75e6a4aSMatthew G. Knepley 
1478a75e6a4aSMatthew G. Knepley   Input Parameter:
1479811af0c4SBarry Smith . comm - the MPI communicator to share the HDF5 `PetscViewer`
1480a75e6a4aSMatthew G. Knepley 
1481a75e6a4aSMatthew G. Knepley   Level: intermediate
1482a75e6a4aSMatthew G. Knepley 
1483811af0c4SBarry Smith   Options Database Key:
148410699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1485a75e6a4aSMatthew G. Knepley 
1486811af0c4SBarry Smith   Environmental variable:
1487811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1488a75e6a4aSMatthew G. Knepley 
1489811af0c4SBarry Smith   Note:
1490811af0c4SBarry Smith   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1491811af0c4SBarry Smith   an error code.  The HDF5 `PetscViewer` is usually used in the form
1492a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1493a75e6a4aSMatthew G. Knepley 
1494811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1495a75e6a4aSMatthew G. Knepley @*/
14969371c9d4SSatish Balay PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) {
1497a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1498a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1499a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1500a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1501a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1502a75e6a4aSMatthew G. Knepley 
1503a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
15049371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
15059371c9d4SSatish Balay   if (ierr) {
15069371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15079371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15089371c9d4SSatish Balay   }
1509a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1510908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
15119371c9d4SSatish Balay     if (ierr) {
15129371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15139371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15149371c9d4SSatish Balay     }
1515a75e6a4aSMatthew G. Knepley   }
151647435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
15179371c9d4SSatish Balay   if (ierr) {
15189371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15199371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15209371c9d4SSatish Balay   }
1521a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1522a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
15239371c9d4SSatish Balay     if (ierr) {
15249371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15259371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15269371c9d4SSatish Balay     }
1527a75e6a4aSMatthew G. Knepley     if (!flg) {
1528a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname, "output.h5");
15299371c9d4SSatish Balay       if (ierr) {
15309371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15319371c9d4SSatish Balay         PetscFunctionReturn(NULL);
15329371c9d4SSatish Balay       }
1533a75e6a4aSMatthew G. Knepley     }
1534a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
15359371c9d4SSatish Balay     if (ierr) {
15369371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15379371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15389371c9d4SSatish Balay     }
1539a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15409371c9d4SSatish Balay     if (ierr) {
15419371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15429371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15439371c9d4SSatish Balay     }
154447435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
15459371c9d4SSatish Balay     if (ierr) {
15469371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15479371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15489371c9d4SSatish Balay     }
1549a75e6a4aSMatthew G. Knepley   }
1550a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
15519371c9d4SSatish Balay   if (ierr) {
15529371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15539371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15549371c9d4SSatish Balay   }
1555a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1556a75e6a4aSMatthew G. Knepley }
1557