xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision 811af0c4b09a35de4306c442f88bd09fdc09897d)
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 
141*811af0c4SBarry Smith     Logically Collective on viewer
14282ea9b62SBarry Smith 
14382ea9b62SBarry Smith   Input Parameters:
144*811af0c4SBarry Smith +  viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
145*811af0c4SBarry 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 
147*811af0c4SBarry 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 
150*811af0c4SBarry 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 
156*811af0c4SBarry 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 
169*811af0c4SBarry Smith     Logically Collective on viewer
17082ea9b62SBarry Smith 
17182ea9b62SBarry Smith   Input Parameter:
172*811af0c4SBarry Smith .  viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`
17382ea9b62SBarry Smith 
17482ea9b62SBarry Smith   Output Parameter:
175*811af0c4SBarry 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 
177*811af0c4SBarry 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 
183*811af0c4SBarry 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
204*811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2059a0502c6SHåkon Strandenes 
206*811af0c4SBarry Smith     Logically Collective on viewer
2079a0502c6SHåkon Strandenes 
2089a0502c6SHåkon Strandenes   Input Parameters:
209*811af0c4SBarry Smith +  viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
210*811af0c4SBarry Smith -  flg - if `PETSC_TRUE` the data will be written to disk with single precision
2119a0502c6SHåkon Strandenes 
212*811af0c4SBarry Smith   Options Database Key:
2139a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2149a0502c6SHåkon Strandenes 
215*811af0c4SBarry 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 
221*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
222*811af0c4SBarry 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
233*811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2349a0502c6SHåkon Strandenes 
235*811af0c4SBarry Smith     Logically Collective on viewer
2369a0502c6SHåkon Strandenes 
2379a0502c6SHåkon Strandenes   Input Parameter:
238*811af0c4SBarry Smith .  viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`
2399a0502c6SHåkon Strandenes 
2409a0502c6SHåkon Strandenes   Output Parameter:
241*811af0c4SBarry 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()`,
246*811af0c4SBarry 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:
278*811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
279*811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)
280ee8b9732SVaclav Hapla 
281*811af0c4SBarry Smith   Options Database Key:
282ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
283ee8b9732SVaclav Hapla 
284*811af0c4SBarry 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 
288*811af0c4SBarry Smith   Developer Note:
289*811af0c4SBarry 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 
295*811af0c4SBarry 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:
327*811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer`
328ee8b9732SVaclav Hapla 
329ee8b9732SVaclav Hapla   Output Parameters:
330ee8b9732SVaclav Hapla . flg - the flag
331ee8b9732SVaclav Hapla 
332ee8b9732SVaclav Hapla   Level: intermediate
333ee8b9732SVaclav Hapla 
334*811af0c4SBarry Smith   Note:
335*811af0c4SBarry 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.
336*811af0c4SBarry Smith   For more details, see `PetscViewerHDF5SetCollective()`.
337ee8b9732SVaclav Hapla 
338*811af0c4SBarry 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 
431*811af0c4SBarry Smith   Logically Collective on viewer
43219a20e4cSMatthew G. Knepley 
43319a20e4cSMatthew G. Knepley   Input Parameters:
434*811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
435*811af0c4SBarry Smith - flg    - if `PETSC_TRUE` we will assume that timestepping is on
43619a20e4cSMatthew G. Knepley 
437*811af0c4SBarry Smith   Options Database Key:
43819a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
43919a20e4cSMatthew G. Knepley 
440*811af0c4SBarry 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 
445*811af0c4SBarry 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:
460*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
46119a20e4cSMatthew G. Knepley 
46219a20e4cSMatthew G. Knepley   Output Parameter:
463*811af0c4SBarry Smith . flg    - if `PETSC_TRUE` we will assume that timestepping is on
46419a20e4cSMatthew G. Knepley 
46519a20e4cSMatthew G. Knepley   Level: intermediate
46619a20e4cSMatthew G. Knepley 
467*811af0c4SBarry 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 
479*811af0c4SBarry Smith   Level: beginner
480*811af0c4SBarry 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 
4999566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(v, &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
528*811af0c4SBarry 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:
538*811af0c4SBarry Smith .  hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
5395c6c1daeSBarry Smith 
540*811af0c4SBarry 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
548*811af0c4SBarry Smith .vb
549*811af0c4SBarry Smith   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
550*811af0c4SBarry Smith   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
551*811af0c4SBarry 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]
552*811af0c4SBarry Smith   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
553*811af0c4SBarry Smith .ve
5547e4fd573SVaclav Hapla 
555*811af0c4SBarry 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 
559*811af0c4SBarry 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:
579*811af0c4SBarry 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 
586*811af0c4SBarry 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:
603*811af0c4SBarry 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.
610*811af0c4SBarry Smith .vb
611*811af0c4SBarry 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.
612*811af0c4SBarry Smith   NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
613*811af0c4SBarry Smith   "." means the current group is pushed again.
614*811af0c4SBarry 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 
623*811af0c4SBarry Smith   Developer Note:
6242e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
625820fc2d1SVaclav Hapla 
626*811af0c4SBarry 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 
6661ad4d0e4SMatthew G. Knepley /*@C
6671ad4d0e4SMatthew G. Knepley   PetscViewerHDF5PushGroupRelative - Set the current HDF5 group for output, appending the path in the argument
6681ad4d0e4SMatthew G. Knepley 
6691ad4d0e4SMatthew G. Knepley   Not collective
6701ad4d0e4SMatthew G. Knepley 
6711ad4d0e4SMatthew G. Knepley   Input Parameters:
672*811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
6731ad4d0e4SMatthew G. Knepley - name - The group name to append
6741ad4d0e4SMatthew G. Knepley 
6751ad4d0e4SMatthew G. Knepley   Level: intermediate
6761ad4d0e4SMatthew G. Knepley 
677*811af0c4SBarry Smith   Note:
6781ad4d0e4SMatthew G. Knepley   This is designed to mnemonically resemble the Unix cd command.
679*811af0c4SBarry Smith .vb
680*811af0c4SBarry 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.
681*811af0c4SBarry Smith   NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
682*811af0c4SBarry Smith   "." means the current group is pushed again.
683*811af0c4SBarry Smith .ve
6841ad4d0e4SMatthew G. Knepley 
6851ad4d0e4SMatthew G. Knepley   Example:
6861ad4d0e4SMatthew G. Knepley   Suppose the current group is "/a".
6871ad4d0e4SMatthew G. Knepley   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
6881ad4d0e4SMatthew G. Knepley   . If name is ".", then the new group will be "/a".
6891ad4d0e4SMatthew G. Knepley   . If name is "b", then the new group will be "/a/b".
6901ad4d0e4SMatthew G. Knepley   - If name is "/b", then the new group will be "/b".
6911ad4d0e4SMatthew G. Knepley 
692*811af0c4SBarry Smith   Developer Note:
6931ad4d0e4SMatthew G. Knepley   The root group "/" is internally stored as NULL.
6941ad4d0e4SMatthew G. Knepley 
695*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5Open()`
6961ad4d0e4SMatthew G. Knepley @*/
6979371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushGroupRelative(PetscViewer viewer, const char name[]) {
6981ad4d0e4SMatthew G. Knepley   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6991ad4d0e4SMatthew G. Knepley   PetscViewerHDF5GroupList *groupNode;
7001ad4d0e4SMatthew G. Knepley   char                      groupname[PETSC_MAX_PATH_LEN];
7011ad4d0e4SMatthew G. Knepley 
7021ad4d0e4SMatthew G. Knepley   PetscFunctionBegin;
7031ad4d0e4SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7041ad4d0e4SMatthew G. Knepley   if (name) PetscValidCharPointer(name, 2);
7051ad4d0e4SMatthew G. Knepley   if (name && name[0]) {
7061ad4d0e4SMatthew G. Knepley     size_t i, len;
7071ad4d0e4SMatthew G. Knepley     PetscCall(PetscStrlen(name, &len));
7089371c9d4SSatish Balay     for (i = 0; i < len; i++)
7099371c9d4SSatish Balay       if (name[i] != '/') break;
7101ad4d0e4SMatthew G. Knepley     if (i == len) name = NULL;
7111ad4d0e4SMatthew G. Knepley   } else name = NULL;
7121ad4d0e4SMatthew G. Knepley   PetscCall(PetscNew(&groupNode));
7131ad4d0e4SMatthew G. Knepley   PetscCall(PetscStrncpy(groupname, hdf5->groups->name, sizeof(groupname)));
7141ad4d0e4SMatthew G. Knepley   if (name[0] != '/') PetscCall(PetscStrlcat(groupname, "/", sizeof(groupname)));
7151ad4d0e4SMatthew G. Knepley   PetscCall(PetscStrlcat(groupname, name, sizeof(groupname)));
7161ad4d0e4SMatthew G. Knepley   PetscCall(PetscStrallocpy(groupname, (char **)&groupNode->name));
7171ad4d0e4SMatthew G. Knepley   groupNode->next = hdf5->groups;
7181ad4d0e4SMatthew G. Knepley   hdf5->groups    = groupNode;
7191ad4d0e4SMatthew G. Knepley   PetscFunctionReturn(0);
7201ad4d0e4SMatthew G. Knepley }
7211ad4d0e4SMatthew G. Knepley 
7223ef9c667SSatish Balay /*@
7235c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
7245c6c1daeSBarry Smith 
7255c6c1daeSBarry Smith   Not collective
7265c6c1daeSBarry Smith 
7275c6c1daeSBarry Smith   Input Parameter:
728*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith   Level: intermediate
7315c6c1daeSBarry Smith 
732*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
7335c6c1daeSBarry Smith @*/
7349371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) {
7355c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7367d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
7375c6c1daeSBarry Smith 
7385c6c1daeSBarry Smith   PetscFunctionBegin;
7395c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7405f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
7415c6c1daeSBarry Smith   groupNode    = hdf5->groups;
7425c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
7439566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
7449566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
7455c6c1daeSBarry Smith   PetscFunctionReturn(0);
7465c6c1daeSBarry Smith }
7475c6c1daeSBarry Smith 
7485c6c1daeSBarry Smith /*@C
749*811af0c4SBarry Smith   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
750874270d9SVaclav Hapla   If none has been assigned, returns NULL.
7515c6c1daeSBarry Smith 
7525c6c1daeSBarry Smith   Not collective
7535c6c1daeSBarry Smith 
7545c6c1daeSBarry Smith   Input Parameter:
755*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7565c6c1daeSBarry Smith 
7575c6c1daeSBarry Smith   Output Parameter:
7585c6c1daeSBarry Smith . name - The group name
7595c6c1daeSBarry Smith 
7605c6c1daeSBarry Smith   Level: intermediate
7615c6c1daeSBarry Smith 
762*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`
7635c6c1daeSBarry Smith @*/
7649371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) {
7655c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7665c6c1daeSBarry Smith 
7675c6c1daeSBarry Smith   PetscFunctionBegin;
7685c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
769c959eef4SJed Brown   PetscValidPointer(name, 2);
770a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
7710298fd71SBarry Smith   else *name = NULL;
7725c6c1daeSBarry Smith   PetscFunctionReturn(0);
7735c6c1daeSBarry Smith }
7745c6c1daeSBarry Smith 
7758c557b5aSVaclav Hapla /*@
776*811af0c4SBarry Smith   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
777874270d9SVaclav Hapla   and return this group's ID and file ID.
778*811af0c4SBarry Smith   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
779874270d9SVaclav Hapla 
780874270d9SVaclav Hapla   Not collective
781874270d9SVaclav Hapla 
782874270d9SVaclav Hapla   Input Parameter:
783*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
784874270d9SVaclav Hapla 
785d8d19677SJose E. Roman   Output Parameters:
786874270d9SVaclav Hapla + fileId - The HDF5 file ID
787874270d9SVaclav Hapla - groupId - The HDF5 group ID
788874270d9SVaclav Hapla 
789*811af0c4SBarry Smith   Note:
790e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
791e74616beSVaclav Hapla 
792874270d9SVaclav Hapla   Level: intermediate
793874270d9SVaclav Hapla 
794*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
795874270d9SVaclav Hapla @*/
7969371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) {
79790e03537SVaclav Hapla   hid_t       file_id;
79890e03537SVaclav Hapla   H5O_type_t  type;
79976d59af2SVaclav Hapla   const char *groupName = NULL, *fileName = NULL;
80076d59af2SVaclav Hapla   PetscBool   writable, has;
80154dbf706SVaclav Hapla 
80254dbf706SVaclav Hapla   PetscFunctionBegin;
8039566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
8049566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
8059566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
8069566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName));
8079566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
80876d59af2SVaclav Hapla   if (!has) {
8095f80ce2aSJacob 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);
810f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
81176d59af2SVaclav Hapla   }
8125f80ce2aSJacob 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);
813792fecdfSBarry Smith   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName ? groupName : "/", H5P_DEFAULT));
81454dbf706SVaclav Hapla   *fileId = file_id;
81554dbf706SVaclav Hapla   PetscFunctionReturn(0);
81654dbf706SVaclav Hapla }
81754dbf706SVaclav Hapla 
8183ef9c667SSatish Balay /*@
819d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
8205c6c1daeSBarry Smith 
8215c6c1daeSBarry Smith   Not collective
8225c6c1daeSBarry Smith 
8235c6c1daeSBarry Smith   Input Parameter:
824*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
8255c6c1daeSBarry Smith 
8265c6c1daeSBarry Smith   Level: intermediate
8275c6c1daeSBarry Smith 
828d7dd068bSVaclav Hapla   Notes:
829*811af0c4SBarry Smith   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
830*811af0c4SBarry Smith   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
831*811af0c4SBarry 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()`].
832*811af0c4SBarry Smith   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
833*811af0c4SBarry Smith   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
834d7dd068bSVaclav Hapla 
835d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
836d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
837d7dd068bSVaclav Hapla 
838*811af0c4SBarry Smith   Developer note:
839d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
840d7dd068bSVaclav Hapla 
841*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
842d7dd068bSVaclav Hapla @*/
8439371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) {
844d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
845d7dd068bSVaclav Hapla 
846d7dd068bSVaclav Hapla   PetscFunctionBegin;
847d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8485f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
849d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
850d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
851d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
852d7dd068bSVaclav Hapla }
853d7dd068bSVaclav Hapla 
854d7dd068bSVaclav Hapla /*@
855d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
856d7dd068bSVaclav Hapla 
857d7dd068bSVaclav Hapla   Not collective
858d7dd068bSVaclav Hapla 
859d7dd068bSVaclav Hapla   Input Parameter:
860*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
861d7dd068bSVaclav Hapla 
862d7dd068bSVaclav Hapla   Level: intermediate
863d7dd068bSVaclav Hapla 
864*811af0c4SBarry Smith   Note:
865*811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
866d7dd068bSVaclav Hapla 
867*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
868d7dd068bSVaclav Hapla @*/
8699371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) {
870d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
871d7dd068bSVaclav Hapla 
872d7dd068bSVaclav Hapla   PetscFunctionBegin;
873d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8745f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
875d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
876d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
877d7dd068bSVaclav Hapla }
878d7dd068bSVaclav Hapla 
879d7dd068bSVaclav Hapla /*@
880d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
881d7dd068bSVaclav Hapla 
882d7dd068bSVaclav Hapla   Not collective
883d7dd068bSVaclav Hapla 
884d7dd068bSVaclav Hapla   Input Parameter:
885*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
886d7dd068bSVaclav Hapla 
887d7dd068bSVaclav Hapla   Output Parameter:
888d7dd068bSVaclav Hapla . flg - is timestepping active?
889d7dd068bSVaclav Hapla 
890d7dd068bSVaclav Hapla   Level: intermediate
891d7dd068bSVaclav Hapla 
892*811af0c4SBarry Smith   Note:
893*811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
894d7dd068bSVaclav Hapla 
895*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
896d7dd068bSVaclav Hapla @*/
8979371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) {
898d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
899d7dd068bSVaclav Hapla 
900d7dd068bSVaclav Hapla   PetscFunctionBegin;
901d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
902d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
903d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
904d7dd068bSVaclav Hapla }
905d7dd068bSVaclav Hapla 
906d7dd068bSVaclav Hapla /*@
907d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
908d7dd068bSVaclav Hapla 
909d7dd068bSVaclav Hapla   Not collective
910d7dd068bSVaclav Hapla 
911d7dd068bSVaclav Hapla   Input Parameter:
912*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
913d7dd068bSVaclav Hapla 
914d7dd068bSVaclav Hapla   Level: intermediate
915d7dd068bSVaclav Hapla 
916*811af0c4SBarry Smith   Note:
917*811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
918d7dd068bSVaclav Hapla 
919*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
9205c6c1daeSBarry Smith @*/
9219371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) {
9225c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9235c6c1daeSBarry Smith 
9245c6c1daeSBarry Smith   PetscFunctionBegin;
9255c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9265f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9275c6c1daeSBarry Smith   ++hdf5->timestep;
9285c6c1daeSBarry Smith   PetscFunctionReturn(0);
9295c6c1daeSBarry Smith }
9305c6c1daeSBarry Smith 
9313ef9c667SSatish Balay /*@
932d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
9335c6c1daeSBarry Smith 
934d7dd068bSVaclav Hapla   Logically collective
9355c6c1daeSBarry Smith 
9365c6c1daeSBarry Smith   Input Parameters:
937*811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
938d7dd068bSVaclav Hapla - timestep - The timestep
9395c6c1daeSBarry Smith 
9405c6c1daeSBarry Smith   Level: intermediate
9415c6c1daeSBarry Smith 
942*811af0c4SBarry Smith   Note:
943*811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
944d7dd068bSVaclav Hapla 
945*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
9465c6c1daeSBarry Smith @*/
9479371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) {
9485c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9495c6c1daeSBarry Smith 
9505c6c1daeSBarry Smith   PetscFunctionBegin;
9515c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
952d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
9535f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
9545f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9555c6c1daeSBarry Smith   hdf5->timestep = timestep;
9565c6c1daeSBarry Smith   PetscFunctionReturn(0);
9575c6c1daeSBarry Smith }
9585c6c1daeSBarry Smith 
9593ef9c667SSatish Balay /*@
9605c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
9615c6c1daeSBarry Smith 
9625c6c1daeSBarry Smith   Not collective
9635c6c1daeSBarry Smith 
9645c6c1daeSBarry Smith   Input Parameter:
965*811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
9665c6c1daeSBarry Smith 
9675c6c1daeSBarry Smith   Output Parameter:
968d7dd068bSVaclav Hapla . timestep - The timestep
9695c6c1daeSBarry Smith 
9705c6c1daeSBarry Smith   Level: intermediate
9715c6c1daeSBarry Smith 
972*811af0c4SBarry Smith   Note:
973*811af0c4SBarry Smith   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
974d7dd068bSVaclav Hapla 
975*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
9765c6c1daeSBarry Smith @*/
9779371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) {
9785c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9795c6c1daeSBarry Smith 
9805c6c1daeSBarry Smith   PetscFunctionBegin;
9815c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
982d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep, 2);
9835f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9845c6c1daeSBarry Smith   *timestep = hdf5->timestep;
9855c6c1daeSBarry Smith   PetscFunctionReturn(0);
9865c6c1daeSBarry Smith }
9875c6c1daeSBarry Smith 
98836bce6e8SMatthew G. Knepley /*@C
98936bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
99036bce6e8SMatthew G. Knepley 
99136bce6e8SMatthew G. Knepley   Not collective
99236bce6e8SMatthew G. Knepley 
99336bce6e8SMatthew G. Knepley   Input Parameter:
994*811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
99536bce6e8SMatthew G. Knepley 
99636bce6e8SMatthew G. Knepley   Output Parameter:
997*811af0c4SBarry Smith . mtype - the HDF5  datatype
99836bce6e8SMatthew G. Knepley 
99936bce6e8SMatthew G. Knepley   Level: advanced
100036bce6e8SMatthew G. Knepley 
1001db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
100236bce6e8SMatthew G. Knepley @*/
10039371c9d4SSatish Balay PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) {
100436bce6e8SMatthew G. Knepley   PetscFunctionBegin;
100536bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
100636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
100736bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_LLONG;
100836bce6e8SMatthew G. Knepley #else
100936bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_INT;
101036bce6e8SMatthew G. Knepley #endif
101136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
101236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
101336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1014c9450970SBarry Smith   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
1015de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
101636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
101736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
101836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
10197e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
102036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
102136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
102236bce6e8SMatthew G. Knepley }
102336bce6e8SMatthew G. Knepley 
102436bce6e8SMatthew G. Knepley /*@C
102536bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
102636bce6e8SMatthew G. Knepley 
102736bce6e8SMatthew G. Knepley   Not collective
102836bce6e8SMatthew G. Knepley 
102936bce6e8SMatthew G. Knepley   Input Parameter:
1030*811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
103136bce6e8SMatthew G. Knepley 
103236bce6e8SMatthew G. Knepley   Output Parameter:
1033*811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
103436bce6e8SMatthew G. Knepley 
103536bce6e8SMatthew G. Knepley   Level: advanced
103636bce6e8SMatthew G. Knepley 
1037db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
103836bce6e8SMatthew G. Knepley @*/
10399371c9d4SSatish Balay PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) {
104036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
104136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
104236bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
104336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
104436bce6e8SMatthew G. Knepley #else
104536bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
104636bce6e8SMatthew G. Knepley #endif
104736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
104836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
104936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
105036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
105136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
105236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
10537e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
105436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
105536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
105636bce6e8SMatthew G. Knepley }
105736bce6e8SMatthew G. Knepley 
1058df863907SAlex Fikl /*@C
1059b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
106036bce6e8SMatthew G. Knepley 
1061343a20b2SVaclav Hapla   Collective
1062343a20b2SVaclav Hapla 
106336bce6e8SMatthew G. Knepley   Input Parameters:
1064*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
10654302210dSVaclav Hapla . parent - The parent dataset/group name
106636bce6e8SMatthew G. Knepley . name   - The attribute name
106736bce6e8SMatthew G. Knepley . datatype - The attribute type
106836bce6e8SMatthew G. Knepley - value    - The attribute value
106936bce6e8SMatthew G. Knepley 
107036bce6e8SMatthew G. Knepley   Level: advanced
107136bce6e8SMatthew G. Knepley 
1072*811af0c4SBarry Smith   Note:
1073343a20b2SVaclav 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.
10744302210dSVaclav Hapla 
1075*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1076*811af0c4SBarry Smith           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
107736bce6e8SMatthew G. Knepley @*/
10789371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) {
10794302210dSVaclav Hapla   char     *parentAbsPath;
108060568592SMatthew G. Knepley   hid_t     h5, dataspace, obj, attribute, dtype;
1081080f144cSVaclav Hapla   PetscBool has;
108236bce6e8SMatthew G. Knepley 
108336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
10845cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
10854302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1086c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
10874302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1088b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
10899566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
10909566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
10919566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
10929566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
10937e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
10947e97c476SMatthew G. Knepley     size_t len;
10959566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *)value, &len));
1096792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
10977e97c476SMatthew G. Knepley   }
10989566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1099792fecdfSBarry Smith   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1100792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1101080f144cSVaclav Hapla   if (has) {
1102792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1103080f144cSVaclav Hapla   } else {
1104792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1105080f144cSVaclav Hapla   }
1106792fecdfSBarry Smith   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1107792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1108792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1109792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
1110792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (dataspace));
11119566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
111236bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
111336bce6e8SMatthew G. Knepley }
111436bce6e8SMatthew G. Knepley 
1115df863907SAlex Fikl /*@C
1116*811af0c4SBarry Smith  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1117577e0e04SVaclav Hapla 
1118343a20b2SVaclav Hapla   Collective
1119343a20b2SVaclav Hapla 
1120577e0e04SVaclav Hapla   Input Parameters:
1121*811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1122577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1123577e0e04SVaclav Hapla . name     - The attribute name
1124577e0e04SVaclav Hapla . datatype - The attribute type
1125577e0e04SVaclav Hapla - value    - The attribute value
1126577e0e04SVaclav Hapla 
1127*811af0c4SBarry Smith   Note:
1128343a20b2SVaclav 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).
1129*811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1130577e0e04SVaclav Hapla 
1131577e0e04SVaclav Hapla   Level: advanced
1132577e0e04SVaclav Hapla 
1133*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1134*811af0c4SBarry Smith           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1135577e0e04SVaclav Hapla @*/
11369371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) {
1137577e0e04SVaclav Hapla   PetscFunctionBegin;
1138577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1139577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1140577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1141b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
11429566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
11439566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1144577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1145577e0e04SVaclav Hapla }
1146577e0e04SVaclav Hapla 
1147577e0e04SVaclav Hapla /*@C
1148b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1149ad0c61feSMatthew G. Knepley 
1150343a20b2SVaclav Hapla   Collective
1151343a20b2SVaclav Hapla 
1152ad0c61feSMatthew G. Knepley   Input Parameters:
1153*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
11544302210dSVaclav Hapla . parent - The parent dataset/group name
1155ad0c61feSMatthew G. Knepley . name   - The attribute name
1156a2d6be1bSVaclav Hapla . datatype - The attribute type
1157a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1158ad0c61feSMatthew G. Knepley 
1159ad0c61feSMatthew G. Knepley   Output Parameter:
1160a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1161ad0c61feSMatthew G. Knepley 
1162a2d6be1bSVaclav Hapla   Notes:
1163a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1164a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1165a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1166*811af0c4SBarry Smith $  flg = `PETSC_FALSE`;
1167*811af0c4SBarry Smith $  PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1168a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1169a2d6be1bSVaclav Hapla 
1170*811af0c4SBarry Smith   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1171708d4cb3SBarry Smith 
1172343a20b2SVaclav 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.
1173ad0c61feSMatthew G. Knepley 
1174343a20b2SVaclav Hapla   Level: advanced
11754302210dSVaclav Hapla 
1176*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1177ad0c61feSMatthew G. Knepley @*/
11789371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) {
11794302210dSVaclav Hapla   char     *parentAbsPath;
1180a2d6be1bSVaclav Hapla   hid_t     h5, obj, attribute, dtype;
1181969748fdSVaclav Hapla   PetscBool has;
1182ad0c61feSMatthew G. Knepley 
1183ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
11845cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
11854302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1186c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1187a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1188a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
11899566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11909566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
11919566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
11929566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1193a2d6be1bSVaclav Hapla   if (!has) {
1194a2d6be1bSVaclav Hapla     if (defaultValue) {
1195a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1196a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
11979566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1198a2d6be1bSVaclav Hapla         } else {
1199a2d6be1bSVaclav Hapla           size_t len;
1200792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
12019566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1202a2d6be1bSVaclav Hapla         }
1203a2d6be1bSVaclav Hapla       }
12049566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
1205a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
120698921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1207a2d6be1bSVaclav Hapla   }
12089566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1209792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1210792fecdfSBarry Smith   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1211f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1212f0b58479SMatthew G. Knepley     size_t len;
1213a2d6be1bSVaclav Hapla     hid_t  atype;
1214792fecdfSBarry Smith     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1215792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
12169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1217792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1218792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1219708d4cb3SBarry Smith   } else {
1220792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1221708d4cb3SBarry Smith   }
1222792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1223e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1224792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
12259566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
1226ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1227ad0c61feSMatthew G. Knepley }
1228ad0c61feSMatthew G. Knepley 
1229577e0e04SVaclav Hapla /*@C
1230*811af0c4SBarry Smith  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1231577e0e04SVaclav Hapla 
1232343a20b2SVaclav Hapla   Collective
1233343a20b2SVaclav Hapla 
1234577e0e04SVaclav Hapla   Input Parameters:
1235*811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1236577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1237577e0e04SVaclav Hapla . name     - The attribute name
1238577e0e04SVaclav Hapla - datatype - The attribute type
1239577e0e04SVaclav Hapla 
1240577e0e04SVaclav Hapla   Output Parameter:
1241577e0e04SVaclav Hapla . value    - The attribute value
1242577e0e04SVaclav Hapla 
1243*811af0c4SBarry Smith   Note:
1244577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1245*811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1246577e0e04SVaclav Hapla 
1247577e0e04SVaclav Hapla   Level: advanced
1248577e0e04SVaclav Hapla 
1249*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1250577e0e04SVaclav Hapla @*/
12519371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) {
1252577e0e04SVaclav Hapla   PetscFunctionBegin;
1253577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1254577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1255577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1256064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
12579566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12589566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1259577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1260577e0e04SVaclav Hapla }
1261577e0e04SVaclav Hapla 
12629371c9d4SSatish Balay static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) {
1263a75c3fd4SVaclav Hapla   htri_t exists;
1264a75c3fd4SVaclav Hapla   hid_t  group;
1265a75c3fd4SVaclav Hapla 
1266a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1267792fecdfSBarry Smith   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1268792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1269a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1270792fecdfSBarry Smith     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1271792fecdfSBarry Smith     PetscCallHDF5(H5Gclose, (group));
1272a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1273a75c3fd4SVaclav Hapla   }
1274a75c3fd4SVaclav Hapla   *exists_ = (PetscBool)exists;
1275a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1276a75c3fd4SVaclav Hapla }
1277a75c3fd4SVaclav Hapla 
12789371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) {
127990e03537SVaclav Hapla   const char rootGroupName[] = "/";
12805cdeb108SMatthew G. Knepley   hid_t      h5;
1281e5bf9ebcSVaclav Hapla   PetscBool  exists = PETSC_FALSE;
128215b861d2SVaclav Hapla   PetscInt   i;
128315b861d2SVaclav Hapla   int        n;
128485688ae2SVaclav Hapla   char     **hierarchy;
128585688ae2SVaclav Hapla   char       buf[PETSC_MAX_PATH_LEN] = "";
12865cdeb108SMatthew G. Knepley 
12875cdeb108SMatthew G. Knepley   PetscFunctionBegin;
12885cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
128990e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
129090e03537SVaclav Hapla   else name = rootGroupName;
1291ccd66a83SVaclav Hapla   if (has) {
1292064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
12935cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1294ccd66a83SVaclav Hapla   }
1295ccd66a83SVaclav Hapla   if (otype) {
1296064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
129756cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1298ccd66a83SVaclav Hapla   }
12999566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
130085688ae2SVaclav Hapla 
130185688ae2SVaclav Hapla   /*
130285688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
130385688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
130485688ae2SVaclav Hapla      1) whether it's a valid link
130585688ae2SVaclav Hapla      2) whether this link resolves to an object
130685688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
130785688ae2SVaclav Hapla   */
13089566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
130985688ae2SVaclav Hapla   if (!n) {
131085688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1311ccd66a83SVaclav Hapla     if (has) *has = PETSC_TRUE;
1312ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
13139566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
131485688ae2SVaclav Hapla     PetscFunctionReturn(0);
131585688ae2SVaclav Hapla   }
131685688ae2SVaclav Hapla   for (i = 0; i < n; i++) {
13179566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
13189566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, hierarchy[i]));
13199566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1320a75c3fd4SVaclav Hapla     if (!exists) break;
132185688ae2SVaclav Hapla   }
13229566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
132385688ae2SVaclav Hapla 
132485688ae2SVaclav Hapla   /* If the object exists, get its type */
1325ccd66a83SVaclav Hapla   if (exists && otype) {
13265cdeb108SMatthew G. Knepley     H5O_info_t info;
13275cdeb108SMatthew G. Knepley 
132876276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1329792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
133056cc0592SVaclav Hapla     *otype = info.type;
13315cdeb108SMatthew G. Knepley   }
1332ccd66a83SVaclav Hapla   if (has) *has = exists;
13335cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
13345cdeb108SMatthew G. Knepley }
13355cdeb108SMatthew G. Knepley 
13364302210dSVaclav Hapla /*@C
13370a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
13380a9f38efSVaclav Hapla 
1339343a20b2SVaclav Hapla   Collective
1340343a20b2SVaclav Hapla 
13410a9f38efSVaclav Hapla   Input Parameters:
1342*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
13434302210dSVaclav Hapla - path - The group path
13440a9f38efSVaclav Hapla 
13450a9f38efSVaclav Hapla   Output Parameter:
13460a9f38efSVaclav Hapla . has - Flag for group existence
13470a9f38efSVaclav Hapla 
13480a9f38efSVaclav Hapla   Level: advanced
13490a9f38efSVaclav Hapla 
13504302210dSVaclav Hapla   Notes:
1351785c443eSVaclav 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.
1352785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
13534302210dSVaclav Hapla 
1354*811af0c4SBarry Smith   If path exists but is not a group, `PETSC_FALSE` is returned.
1355*811af0c4SBarry Smith 
1356*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
13570a9f38efSVaclav Hapla @*/
13589371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) {
13590a9f38efSVaclav Hapla   H5O_type_t type;
13604302210dSVaclav Hapla   char      *abspath;
13610a9f38efSVaclav Hapla 
13620a9f38efSVaclav Hapla   PetscFunctionBegin;
13630a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
13644302210dSVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
13654302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
13669566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
13679566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
13684302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
13699566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
13700a9f38efSVaclav Hapla   PetscFunctionReturn(0);
13710a9f38efSVaclav Hapla }
13720a9f38efSVaclav Hapla 
137389e0ef10SVaclav Hapla /*@C
137489e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
137589e0ef10SVaclav Hapla 
1376343a20b2SVaclav Hapla   Collective
1377343a20b2SVaclav Hapla 
137889e0ef10SVaclav Hapla   Input Parameters:
1379*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
138089e0ef10SVaclav Hapla - path - The dataset path
138189e0ef10SVaclav Hapla 
138289e0ef10SVaclav Hapla   Output Parameter:
138389e0ef10SVaclav Hapla . has - Flag whether dataset exists
138489e0ef10SVaclav Hapla 
138589e0ef10SVaclav Hapla   Level: advanced
138689e0ef10SVaclav Hapla 
138789e0ef10SVaclav Hapla   Notes:
1388343a20b2SVaclav 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.
138989e0ef10SVaclav Hapla 
1390*811af0c4SBarry Smith   If path is NULL or empty, has is set to `PETSC_FALSE`.
1391*811af0c4SBarry Smith 
1392*811af0c4SBarry Smith   If path exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1393*811af0c4SBarry Smith 
1394*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
139589e0ef10SVaclav Hapla @*/
13969371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) {
139789e0ef10SVaclav Hapla   H5O_type_t type;
139889e0ef10SVaclav Hapla   char      *abspath;
139989e0ef10SVaclav Hapla 
140089e0ef10SVaclav Hapla   PetscFunctionBegin;
140189e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
140289e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
140389e0ef10SVaclav Hapla   PetscValidBoolPointer(has, 3);
14049566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath));
14059566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
140689e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
14079566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
140889e0ef10SVaclav Hapla   PetscFunctionReturn(0);
140989e0ef10SVaclav Hapla }
141089e0ef10SVaclav Hapla 
14110a9f38efSVaclav Hapla /*@
1412e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1413ecce1506SVaclav Hapla 
1414343a20b2SVaclav Hapla   Collective
1415343a20b2SVaclav Hapla 
1416ecce1506SVaclav Hapla   Input Parameters:
1417*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1418ecce1506SVaclav Hapla - obj    - The named object
1419ecce1506SVaclav Hapla 
1420ecce1506SVaclav Hapla   Output Parameter:
142189e0ef10SVaclav Hapla . has    - Flag for dataset existence
1422ecce1506SVaclav Hapla 
1423e3f143f7SVaclav Hapla   Notes:
142489e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1425*811af0c4SBarry Smith 
1426*811af0c4SBarry Smith   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1427e3f143f7SVaclav Hapla 
1428ecce1506SVaclav Hapla   Level: advanced
1429ecce1506SVaclav Hapla 
1430*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1431ecce1506SVaclav Hapla @*/
14329371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) {
143389e0ef10SVaclav Hapla   size_t len;
1434ecce1506SVaclav Hapla 
1435ecce1506SVaclav Hapla   PetscFunctionBegin;
1436c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1437c57a1a9aSVaclav Hapla   PetscValidHeader(obj, 2);
14384302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
14399566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
14405f80ce2aSJacob Faibussowitsch   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
14419566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1442ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1443ecce1506SVaclav Hapla }
1444ecce1506SVaclav Hapla 
1445df863907SAlex Fikl /*@C
1446ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1447e7bdbf8eSMatthew G. Knepley 
1448343a20b2SVaclav Hapla   Collective
1449343a20b2SVaclav Hapla 
1450e7bdbf8eSMatthew G. Knepley   Input Parameters:
1451*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
14524302210dSVaclav Hapla . parent - The parent dataset/group name
1453e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1454e7bdbf8eSMatthew G. Knepley 
1455e7bdbf8eSMatthew G. Knepley   Output Parameter:
1456e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1457e7bdbf8eSMatthew G. Knepley 
1458e7bdbf8eSMatthew G. Knepley   Level: advanced
1459e7bdbf8eSMatthew G. Knepley 
1460*811af0c4SBarry Smith   Note:
1461343a20b2SVaclav 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.
14624302210dSVaclav Hapla 
1463*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1464e7bdbf8eSMatthew G. Knepley @*/
14659371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
14664302210dSVaclav Hapla   char *parentAbsPath;
1467e7bdbf8eSMatthew G. Knepley 
1468e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
14695cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14704302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1471c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
14724302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
14739566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath));
14749566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
14759566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
14769566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
147706db490cSVaclav Hapla   PetscFunctionReturn(0);
147806db490cSVaclav Hapla }
147906db490cSVaclav Hapla 
1480577e0e04SVaclav Hapla /*@C
1481*811af0c4SBarry Smith  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1482577e0e04SVaclav Hapla 
1483343a20b2SVaclav Hapla   Collective
1484343a20b2SVaclav Hapla 
1485577e0e04SVaclav Hapla   Input Parameters:
1486*811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1487577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1488577e0e04SVaclav Hapla - name   - The attribute name
1489577e0e04SVaclav Hapla 
1490577e0e04SVaclav Hapla   Output Parameter:
1491577e0e04SVaclav Hapla . has    - Flag for attribute existence
1492577e0e04SVaclav Hapla 
1493*811af0c4SBarry Smith   Note:
1494577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1495*811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1496577e0e04SVaclav Hapla 
1497577e0e04SVaclav Hapla   Level: advanced
1498577e0e04SVaclav Hapla 
1499*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1500577e0e04SVaclav Hapla @*/
15019371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) {
1502577e0e04SVaclav Hapla   PetscFunctionBegin;
1503577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1504577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1505577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
15064302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
15079566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
15089566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1509577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1510577e0e04SVaclav Hapla }
1511577e0e04SVaclav Hapla 
15129371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) {
151306db490cSVaclav Hapla   hid_t  h5;
151406db490cSVaclav Hapla   htri_t hhas;
151506db490cSVaclav Hapla 
151606db490cSVaclav Hapla   PetscFunctionBegin;
15179566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1518792fecdfSBarry Smith   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
15195cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1520e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1521e7bdbf8eSMatthew G. Knepley }
1522e7bdbf8eSMatthew G. Knepley 
1523a75e6a4aSMatthew G. Knepley /*
1524a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1525a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1526a75e6a4aSMatthew G. Knepley */
1527d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1528a75e6a4aSMatthew G. Knepley 
1529a75e6a4aSMatthew G. Knepley /*@C
1530*811af0c4SBarry Smith   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1531a75e6a4aSMatthew G. Knepley 
1532d083f849SBarry Smith   Collective
1533a75e6a4aSMatthew G. Knepley 
1534a75e6a4aSMatthew G. Knepley   Input Parameter:
1535*811af0c4SBarry Smith . comm - the MPI communicator to share the HDF5 `PetscViewer`
1536a75e6a4aSMatthew G. Knepley 
1537a75e6a4aSMatthew G. Knepley   Level: intermediate
1538a75e6a4aSMatthew G. Knepley 
1539*811af0c4SBarry Smith   Options Database Key:
154010699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1541a75e6a4aSMatthew G. Knepley 
1542*811af0c4SBarry Smith   Environmental variable:
1543*811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1544a75e6a4aSMatthew G. Knepley 
1545*811af0c4SBarry Smith   Note:
1546*811af0c4SBarry Smith   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1547*811af0c4SBarry Smith   an error code.  The HDF5 `PetscViewer` is usually used in the form
1548a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1549a75e6a4aSMatthew G. Knepley 
1550*811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1551a75e6a4aSMatthew G. Knepley @*/
15529371c9d4SSatish Balay PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) {
1553a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1554a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1555a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1556a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1557a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1558a75e6a4aSMatthew G. Knepley 
1559a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
15609371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
15619371c9d4SSatish Balay   if (ierr) {
15629371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15639371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15649371c9d4SSatish Balay   }
1565a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1566908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
15679371c9d4SSatish Balay     if (ierr) {
15689371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15699371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15709371c9d4SSatish Balay     }
1571a75e6a4aSMatthew G. Knepley   }
157247435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
15739371c9d4SSatish Balay   if (ierr) {
15749371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
15759371c9d4SSatish Balay     PetscFunctionReturn(NULL);
15769371c9d4SSatish Balay   }
1577a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1578a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
15799371c9d4SSatish Balay     if (ierr) {
15809371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15819371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15829371c9d4SSatish Balay     }
1583a75e6a4aSMatthew G. Knepley     if (!flg) {
1584a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname, "output.h5");
15859371c9d4SSatish Balay       if (ierr) {
15869371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15879371c9d4SSatish Balay         PetscFunctionReturn(NULL);
15889371c9d4SSatish Balay       }
1589a75e6a4aSMatthew G. Knepley     }
1590a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
15919371c9d4SSatish Balay     if (ierr) {
15929371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15939371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15949371c9d4SSatish Balay     }
1595a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
15969371c9d4SSatish Balay     if (ierr) {
15979371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
15989371c9d4SSatish Balay       PetscFunctionReturn(NULL);
15999371c9d4SSatish Balay     }
160047435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
16019371c9d4SSatish Balay     if (ierr) {
16029371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16039371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16049371c9d4SSatish Balay     }
1605a75e6a4aSMatthew G. Knepley   }
1606a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
16079371c9d4SSatish Balay   if (ierr) {
16089371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16099371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16109371c9d4SSatish Balay   }
1611a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1612a75e6a4aSMatthew G. Knepley }
1613