xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision c3339decea92175325d9368fa13196bcd0e0e58b)
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 
677717648SVaclav Hapla /*@C
777717648SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
877717648SVaclav Hapla 
977717648SVaclav Hapla   Not collective
1077717648SVaclav Hapla 
1177717648SVaclav Hapla   Input Parameters:
1277717648SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
1377717648SVaclav Hapla - path - (Optional) The path relative to the pushed group
1477717648SVaclav Hapla 
1577717648SVaclav Hapla   Output Parameter:
1677717648SVaclav Hapla . abspath - The absolute HDF5 path (group)
1777717648SVaclav Hapla 
1877717648SVaclav Hapla   Level: intermediate
1977717648SVaclav Hapla 
2077717648SVaclav Hapla   Notes:
2177717648SVaclav 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.
2277717648SVaclav Hapla   So NULL or empty path means the current pushed group.
2377717648SVaclav Hapla 
2477717648SVaclav Hapla   The output abspath is newly allocated so needs to be freed.
2577717648SVaclav Hapla 
2663cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
2777717648SVaclav Hapla @*/
28d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[])
29d71ae5a4SJacob Faibussowitsch {
3077717648SVaclav Hapla   size_t      len;
314302210dSVaclav Hapla   PetscBool   relative = PETSC_FALSE;
326c132bc1SVaclav Hapla   const char *group;
336c132bc1SVaclav Hapla   char        buf[PETSC_MAX_PATH_LEN] = "";
346c132bc1SVaclav Hapla 
356c132bc1SVaclav Hapla   PetscFunctionBegin;
3677717648SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
3777717648SVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
3877717648SVaclav Hapla   PetscValidPointer(abspath, 3);
3977717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
4077717648SVaclav Hapla   PetscCall(PetscStrlen(path, &len));
4177717648SVaclav Hapla   relative = (PetscBool)(!len || path[0] != '/');
424302210dSVaclav Hapla   if (relative) {
439566063dSJacob Faibussowitsch     PetscCall(PetscStrcpy(buf, group));
4477717648SVaclav Hapla     if (!group || len) PetscCall(PetscStrcat(buf, "/"));
459566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, path));
469566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(buf, abspath));
474302210dSVaclav Hapla   } else {
489566063dSJacob Faibussowitsch     PetscCall(PetscStrallocpy(path, abspath));
494302210dSVaclav Hapla   }
506c132bc1SVaclav Hapla   PetscFunctionReturn(0);
516c132bc1SVaclav Hapla }
526c132bc1SVaclav Hapla 
53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
54d71ae5a4SJacob Faibussowitsch {
55577e0e04SVaclav Hapla   PetscBool has;
56577e0e04SVaclav Hapla 
57577e0e04SVaclav Hapla   PetscFunctionBegin;
589566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
59577e0e04SVaclav Hapla   if (!has) {
6077717648SVaclav Hapla     char *group;
6177717648SVaclav Hapla     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
6277717648SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group);
63577e0e04SVaclav Hapla   }
64577e0e04SVaclav Hapla   PetscFunctionReturn(0);
65577e0e04SVaclav Hapla }
66577e0e04SVaclav Hapla 
67d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject)
68d71ae5a4SJacob Faibussowitsch {
69ee8b9732SVaclav Hapla   PetscBool         flg  = PETSC_FALSE, set;
7082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
7182ea9b62SBarry Smith 
7282ea9b62SBarry Smith   PetscFunctionBegin;
73d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
779566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
7819a20e4cSMatthew G. Knepley   flg = PETSC_FALSE;
799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
809566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
81d0609cedSBarry Smith   PetscOptionsHeadEnd();
8282ea9b62SBarry Smith   PetscFunctionReturn(0);
8382ea9b62SBarry Smith }
8482ea9b62SBarry Smith 
85d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
86d71ae5a4SJacob Faibussowitsch {
871b793a25SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
8803fa8834SVaclav Hapla   PetscBool         flg;
891b793a25SVaclav Hapla 
901b793a25SVaclav Hapla   PetscFunctionBegin;
9148a46eb9SPierre Jolivet   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
929566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
939566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
971b793a25SVaclav Hapla   PetscFunctionReturn(0);
981b793a25SVaclav Hapla }
991b793a25SVaclav Hapla 
100d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
101d71ae5a4SJacob Faibussowitsch {
1025c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1035c6c1daeSBarry Smith 
1045c6c1daeSBarry Smith   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5->filename));
106792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
1075c6c1daeSBarry Smith   PetscFunctionReturn(0);
1085c6c1daeSBarry Smith }
1095c6c1daeSBarry Smith 
110d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
111d71ae5a4SJacob Faibussowitsch {
1126226335aSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1136226335aSVaclav Hapla 
1146226335aSVaclav Hapla   PetscFunctionBegin;
115792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
1166226335aSVaclav Hapla   PetscFunctionReturn(0);
1176226335aSVaclav Hapla }
1186226335aSVaclav Hapla 
119d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
120d71ae5a4SJacob Faibussowitsch {
1215c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1225c6c1daeSBarry Smith 
1235c6c1daeSBarry Smith   PetscFunctionBegin;
124792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
1259566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_HDF5(viewer));
1265c6c1daeSBarry Smith   while (hdf5->groups) {
1277d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
1285c6c1daeSBarry Smith 
1299566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups->name));
1309566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups));
1315c6c1daeSBarry Smith     hdf5->groups = tmp;
1325c6c1daeSBarry Smith   }
1339566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5));
1349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
1359566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
1369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
1372e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
1402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
1412e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
1439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
1445c6c1daeSBarry Smith   PetscFunctionReturn(0);
1455c6c1daeSBarry Smith }
1465c6c1daeSBarry Smith 
147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
148d71ae5a4SJacob Faibussowitsch {
1495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1505c6c1daeSBarry Smith 
1515c6c1daeSBarry Smith   PetscFunctionBegin;
1525c6c1daeSBarry Smith   hdf5->btype = type;
1535c6c1daeSBarry Smith   PetscFunctionReturn(0);
1545c6c1daeSBarry Smith }
1555c6c1daeSBarry Smith 
156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
157d71ae5a4SJacob Faibussowitsch {
1580b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1590b62783dSVaclav Hapla 
1600b62783dSVaclav Hapla   PetscFunctionBegin;
1610b62783dSVaclav Hapla   *type = hdf5->btype;
1620b62783dSVaclav Hapla   PetscFunctionReturn(0);
1630b62783dSVaclav Hapla }
1640b62783dSVaclav Hapla 
165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
166d71ae5a4SJacob Faibussowitsch {
16782ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
16882ea9b62SBarry Smith 
16982ea9b62SBarry Smith   PetscFunctionBegin;
17082ea9b62SBarry Smith   hdf5->basedimension2 = flg;
17182ea9b62SBarry Smith   PetscFunctionReturn(0);
17282ea9b62SBarry Smith }
17382ea9b62SBarry Smith 
174df863907SAlex Fikl /*@
17582ea9b62SBarry Smith      PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
17682ea9b62SBarry Smith        dimension of 2.
17782ea9b62SBarry Smith 
178*c3339decSBarry Smith     Logically Collective
17982ea9b62SBarry Smith 
18082ea9b62SBarry Smith   Input Parameters:
181811af0c4SBarry Smith +  viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
182811af0c4SBarry Smith -  flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1
18382ea9b62SBarry Smith 
184811af0c4SBarry Smith   Options Database Key:
18582ea9b62SBarry 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
18682ea9b62SBarry Smith 
187811af0c4SBarry Smith   Note:
18895452b02SPatrick 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
18982ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
19082ea9b62SBarry Smith 
19182ea9b62SBarry Smith   Level: intermediate
19282ea9b62SBarry Smith 
193811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
19482ea9b62SBarry Smith @*/
195d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
196d71ae5a4SJacob Faibussowitsch {
19782ea9b62SBarry Smith   PetscFunctionBegin;
19882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
199cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
20082ea9b62SBarry Smith   PetscFunctionReturn(0);
20182ea9b62SBarry Smith }
20282ea9b62SBarry Smith 
203df863907SAlex Fikl /*@
20482ea9b62SBarry Smith      PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
20582ea9b62SBarry Smith        dimension of 2.
20682ea9b62SBarry Smith 
207*c3339decSBarry Smith     Logically Collective
20882ea9b62SBarry Smith 
20982ea9b62SBarry Smith   Input Parameter:
210811af0c4SBarry Smith .  viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`
21182ea9b62SBarry Smith 
21282ea9b62SBarry Smith   Output Parameter:
213811af0c4SBarry 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
21482ea9b62SBarry Smith 
215811af0c4SBarry Smith   Note:
21695452b02SPatrick 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
21782ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
21882ea9b62SBarry Smith 
21982ea9b62SBarry Smith   Level: intermediate
22082ea9b62SBarry Smith 
221811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
22282ea9b62SBarry Smith @*/
223d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
224d71ae5a4SJacob Faibussowitsch {
22582ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
22682ea9b62SBarry Smith 
22782ea9b62SBarry Smith   PetscFunctionBegin;
22882ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
22982ea9b62SBarry Smith   *flg = hdf5->basedimension2;
23082ea9b62SBarry Smith   PetscFunctionReturn(0);
23182ea9b62SBarry Smith }
23282ea9b62SBarry Smith 
233d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
234d71ae5a4SJacob Faibussowitsch {
2359a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
2369a0502c6SHåkon Strandenes 
2379a0502c6SHåkon Strandenes   PetscFunctionBegin;
2389a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2399a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2409a0502c6SHåkon Strandenes }
2419a0502c6SHåkon Strandenes 
242df863907SAlex Fikl /*@
2439a0502c6SHåkon Strandenes      PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
244811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2459a0502c6SHåkon Strandenes 
246*c3339decSBarry Smith     Logically Collective
2479a0502c6SHåkon Strandenes 
2489a0502c6SHåkon Strandenes   Input Parameters:
249811af0c4SBarry Smith +  viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
250811af0c4SBarry Smith -  flg - if `PETSC_TRUE` the data will be written to disk with single precision
2519a0502c6SHåkon Strandenes 
252811af0c4SBarry Smith   Options Database Key:
2539a0502c6SHåkon Strandenes .  -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2549a0502c6SHåkon Strandenes 
255811af0c4SBarry Smith   Note:
25695452b02SPatrick Sanan     Setting this option does not make any difference if PETSc is compiled with single precision
2579a0502c6SHåkon Strandenes          in the first place. It does not affect reading datasets (HDF5 handle this internally).
2589a0502c6SHåkon Strandenes 
2599a0502c6SHåkon Strandenes   Level: intermediate
2609a0502c6SHåkon Strandenes 
261811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
262811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
2639a0502c6SHåkon Strandenes @*/
264d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
265d71ae5a4SJacob Faibussowitsch {
2669a0502c6SHåkon Strandenes   PetscFunctionBegin;
2679a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
268cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
2699a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2709a0502c6SHåkon Strandenes }
2719a0502c6SHåkon Strandenes 
272df863907SAlex Fikl /*@
2739a0502c6SHåkon Strandenes      PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
274811af0c4SBarry Smith        compiled with double precision `PetscReal`.
2759a0502c6SHåkon Strandenes 
276*c3339decSBarry Smith     Logically Collective
2779a0502c6SHåkon Strandenes 
2789a0502c6SHåkon Strandenes   Input Parameter:
279811af0c4SBarry Smith .  viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`
2809a0502c6SHåkon Strandenes 
2819a0502c6SHåkon Strandenes   Output Parameter:
282811af0c4SBarry Smith .  flg - if `PETSC_TRUE` the data will be written to disk with single precision
2839a0502c6SHåkon Strandenes 
2849a0502c6SHåkon Strandenes   Level: intermediate
2859a0502c6SHåkon Strandenes 
286db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
287811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
2889a0502c6SHåkon Strandenes @*/
289d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
290d71ae5a4SJacob Faibussowitsch {
2919a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
2929a0502c6SHåkon Strandenes 
2939a0502c6SHåkon Strandenes   PetscFunctionBegin;
2949a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
2959a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
2969a0502c6SHåkon Strandenes   PetscFunctionReturn(0);
2979a0502c6SHåkon Strandenes }
2989a0502c6SHåkon Strandenes 
299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
300d71ae5a4SJacob Faibussowitsch {
301ee8b9732SVaclav Hapla   PetscFunctionBegin;
302ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
303ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
304c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
305ee8b9732SVaclav Hapla   {
306ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
307792fecdfSBarry Smith     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
308ee8b9732SVaclav Hapla   }
309ee8b9732SVaclav Hapla #else
3109566063dSJacob 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"));
311ee8b9732SVaclav Hapla #endif
312ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
313ee8b9732SVaclav Hapla }
314ee8b9732SVaclav Hapla 
315ee8b9732SVaclav Hapla /*@
316ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
317ee8b9732SVaclav Hapla 
318ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
319ee8b9732SVaclav Hapla 
320ee8b9732SVaclav Hapla   Input Parameters:
321811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
322811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)
323ee8b9732SVaclav Hapla 
324811af0c4SBarry Smith   Options Database Key:
325ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
326ee8b9732SVaclav Hapla 
327811af0c4SBarry Smith   Note:
328ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
32953effed7SBarry 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.
330ee8b9732SVaclav Hapla 
331811af0c4SBarry Smith   Developer Note:
332811af0c4SBarry Smith   In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively.
333ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
334ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
335ee8b9732SVaclav Hapla 
336ee8b9732SVaclav Hapla   Level: intermediate
337ee8b9732SVaclav Hapla 
338811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
339ee8b9732SVaclav Hapla @*/
340d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
341d71ae5a4SJacob Faibussowitsch {
342ee8b9732SVaclav Hapla   PetscFunctionBegin;
343ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
344ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer, flg, 2);
345cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
346ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
347ee8b9732SVaclav Hapla }
348ee8b9732SVaclav Hapla 
349d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
350d71ae5a4SJacob Faibussowitsch {
351c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
352ee8b9732SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
353ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
354c796909eSBarry Smith #endif
355ee8b9732SVaclav Hapla 
356ee8b9732SVaclav Hapla   PetscFunctionBegin;
357c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
358c796909eSBarry Smith   *flg = PETSC_FALSE;
359c796909eSBarry Smith #else
360792fecdfSBarry Smith   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
361ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
362c796909eSBarry Smith #endif
363ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
364ee8b9732SVaclav Hapla }
365ee8b9732SVaclav Hapla 
366ee8b9732SVaclav Hapla /*@
367ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
368ee8b9732SVaclav Hapla 
369ee8b9732SVaclav Hapla   Not Collective
370ee8b9732SVaclav Hapla 
371ee8b9732SVaclav Hapla   Input Parameters:
372811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer`
373ee8b9732SVaclav Hapla 
374ee8b9732SVaclav Hapla   Output Parameters:
375ee8b9732SVaclav Hapla . flg - the flag
376ee8b9732SVaclav Hapla 
377ee8b9732SVaclav Hapla   Level: intermediate
378ee8b9732SVaclav Hapla 
379811af0c4SBarry Smith   Note:
380811af0c4SBarry 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.
381811af0c4SBarry Smith   For more details, see `PetscViewerHDF5SetCollective()`.
382ee8b9732SVaclav Hapla 
383811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
384ee8b9732SVaclav Hapla @*/
385d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
386d71ae5a4SJacob Faibussowitsch {
387ee8b9732SVaclav Hapla   PetscFunctionBegin;
388ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
389534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
390ee8b9732SVaclav Hapla 
391cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
392ee8b9732SVaclav Hapla   PetscFunctionReturn(0);
393ee8b9732SVaclav Hapla }
394ee8b9732SVaclav Hapla 
395d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
396d71ae5a4SJacob Faibussowitsch {
3975c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
3985c6c1daeSBarry Smith   hid_t             plist_id;
3995c6c1daeSBarry Smith 
4005c6c1daeSBarry Smith   PetscFunctionBegin;
401792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
4029566063dSJacob Faibussowitsch   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
4039566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &hdf5->filename));
4045c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
405792fecdfSBarry Smith   PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS));
406c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
407792fecdfSBarry Smith   PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
408c796909eSBarry Smith #endif
4095c6c1daeSBarry Smith   /* Create or open the file collectively */
4105c6c1daeSBarry Smith   switch (hdf5->btype) {
4115c6c1daeSBarry Smith   case FILE_MODE_READ:
4128a2871f6SBarry Smith     if (PetscDefined(USE_DEBUG)) {
4138a2871f6SBarry Smith       PetscMPIInt rank;
4148a2871f6SBarry Smith       PetscBool   flg;
4158a2871f6SBarry Smith 
4169566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
4178a2871f6SBarry Smith       if (rank == 0) {
4189566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
4195f80ce2aSJacob Faibussowitsch         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
4208a2871f6SBarry Smith       }
4219566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
4228a2871f6SBarry Smith     }
423792fecdfSBarry Smith     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id));
4245c6c1daeSBarry Smith     break;
4255c6c1daeSBarry Smith   case FILE_MODE_APPEND:
4269371c9d4SSatish Balay   case FILE_MODE_UPDATE: {
4277e4fd573SVaclav Hapla     PetscBool flg;
4289566063dSJacob Faibussowitsch     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
429792fecdfSBarry Smith     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id));
430792fecdfSBarry Smith     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4315c6c1daeSBarry Smith     break;
4327e4fd573SVaclav Hapla   }
433d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
434d71ae5a4SJacob Faibussowitsch     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
435d71ae5a4SJacob Faibussowitsch     break;
436d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
437d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
438d71ae5a4SJacob Faibussowitsch   default:
439d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
4405c6c1daeSBarry Smith   }
4415f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
442792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (plist_id));
443671abe24SVaclav Hapla   PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
4445c6c1daeSBarry Smith   PetscFunctionReturn(0);
4455c6c1daeSBarry Smith }
4465c6c1daeSBarry Smith 
447d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
448d71ae5a4SJacob Faibussowitsch {
449d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
450d1232d7fSToby Isaac 
451d1232d7fSToby Isaac   PetscFunctionBegin;
452d1232d7fSToby Isaac   *name = vhdf5->filename;
453d1232d7fSToby Isaac   PetscFunctionReturn(0);
454d1232d7fSToby Isaac }
455d1232d7fSToby Isaac 
456d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
457d71ae5a4SJacob Faibussowitsch {
4585dc64a97SVaclav Hapla   /*
459b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4605dc64a97SVaclav Hapla   */
461b723ab35SVaclav Hapla 
462b723ab35SVaclav Hapla   PetscFunctionBegin;
463b723ab35SVaclav Hapla   PetscFunctionReturn(0);
464b723ab35SVaclav Hapla }
465b723ab35SVaclav Hapla 
466d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
467d71ae5a4SJacob Faibussowitsch {
46819a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
46919a20e4cSMatthew G. Knepley 
47019a20e4cSMatthew G. Knepley   PetscFunctionBegin;
47119a20e4cSMatthew G. Knepley   hdf5->defTimestepping = flg;
47219a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
47319a20e4cSMatthew G. Knepley }
47419a20e4cSMatthew G. Knepley 
475d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
476d71ae5a4SJacob Faibussowitsch {
47719a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
47819a20e4cSMatthew G. Knepley 
47919a20e4cSMatthew G. Knepley   PetscFunctionBegin;
48019a20e4cSMatthew G. Knepley   *flg = hdf5->defTimestepping;
48119a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
48219a20e4cSMatthew G. Knepley }
48319a20e4cSMatthew G. Knepley 
48419a20e4cSMatthew G. Knepley /*@
48519a20e4cSMatthew G. Knepley   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
48619a20e4cSMatthew G. Knepley 
487*c3339decSBarry Smith   Logically Collective
48819a20e4cSMatthew G. Knepley 
48919a20e4cSMatthew G. Knepley   Input Parameters:
490811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
491811af0c4SBarry Smith - flg    - if `PETSC_TRUE` we will assume that timestepping is on
49219a20e4cSMatthew G. Knepley 
493811af0c4SBarry Smith   Options Database Key:
49419a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
49519a20e4cSMatthew G. Knepley 
496811af0c4SBarry Smith   Note:
49719a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
49819a20e4cSMatthew G. Knepley 
49919a20e4cSMatthew G. Knepley   Level: intermediate
50019a20e4cSMatthew G. Knepley 
501811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
50219a20e4cSMatthew G. Knepley @*/
503d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
504d71ae5a4SJacob Faibussowitsch {
50519a20e4cSMatthew G. Knepley   PetscFunctionBegin;
50619a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
507cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
50819a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
50919a20e4cSMatthew G. Knepley }
51019a20e4cSMatthew G. Knepley 
51119a20e4cSMatthew G. Knepley /*@
51219a20e4cSMatthew G. Knepley   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
51319a20e4cSMatthew G. Knepley 
51419a20e4cSMatthew G. Knepley   Not collective
51519a20e4cSMatthew G. Knepley 
51619a20e4cSMatthew G. Knepley   Input Parameter:
517811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
51819a20e4cSMatthew G. Knepley 
51919a20e4cSMatthew G. Knepley   Output Parameter:
520811af0c4SBarry Smith . flg    - if `PETSC_TRUE` we will assume that timestepping is on
52119a20e4cSMatthew G. Knepley 
52219a20e4cSMatthew G. Knepley   Level: intermediate
52319a20e4cSMatthew G. Knepley 
524811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
52519a20e4cSMatthew G. Knepley @*/
526d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
527d71ae5a4SJacob Faibussowitsch {
52819a20e4cSMatthew G. Knepley   PetscFunctionBegin;
52919a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
530cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
53119a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
53219a20e4cSMatthew G. Knepley }
53319a20e4cSMatthew G. Knepley 
5348556b5ebSBarry Smith /*MC
5358556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
5368556b5ebSBarry Smith 
537811af0c4SBarry Smith   Level: beginner
538811af0c4SBarry Smith 
539db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
540db781477SPatrick Sanan           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
541db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
542db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
5438556b5ebSBarry Smith M*/
544d1232d7fSToby Isaac 
545d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
546d71ae5a4SJacob Faibussowitsch {
5475c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
5485c6c1daeSBarry Smith 
5495c6c1daeSBarry Smith   PetscFunctionBegin;
55099335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
55199335c34SVaclav Hapla   {
55299335c34SVaclav Hapla     PetscMPIInt size;
5539566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
5545f80ce2aSJacob 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)");
55599335c34SVaclav Hapla   }
55699335c34SVaclav Hapla #endif
55799335c34SVaclav Hapla 
5584dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hdf5));
5595c6c1daeSBarry Smith 
5605c6c1daeSBarry Smith   v->data                = (void *)hdf5;
5615c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
56282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
563b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
5641b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
5656226335aSVaclav Hapla   v->ops->flush          = PetscViewerFlush_HDF5;
5667e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
567908793a3SLisandro Dalcin   hdf5->filename         = NULL;
5685c6c1daeSBarry Smith   hdf5->timestep         = -1;
5690298fd71SBarry Smith   hdf5->groups           = NULL;
5705c6c1daeSBarry Smith 
571792fecdfSBarry Smith   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
5729c5072fbSVaclav Hapla 
5739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
5749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
5759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
5769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
5789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
5799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
5809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
5819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
5829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
5835c6c1daeSBarry Smith   PetscFunctionReturn(0);
5845c6c1daeSBarry Smith }
5855c6c1daeSBarry Smith 
5865c6c1daeSBarry Smith /*@C
587811af0c4SBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`
5885c6c1daeSBarry Smith 
589d083f849SBarry Smith    Collective
5905c6c1daeSBarry Smith 
5915c6c1daeSBarry Smith    Input Parameters:
5925c6c1daeSBarry Smith +  comm - MPI communicator
5935c6c1daeSBarry Smith .  name - name of file
5945c6c1daeSBarry Smith -  type - type of file
5955c6c1daeSBarry Smith 
5965c6c1daeSBarry Smith    Output Parameter:
597811af0c4SBarry Smith .  hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
5985c6c1daeSBarry Smith 
599811af0c4SBarry Smith   Options Database Keys:
600a2b725a8SWilliam 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
601a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
60282ea9b62SBarry Smith 
6035c6c1daeSBarry Smith    Level: beginner
6045c6c1daeSBarry Smith 
6057e4fd573SVaclav Hapla    Notes:
6067e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
607811af0c4SBarry Smith .vb
608811af0c4SBarry Smith   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
609811af0c4SBarry Smith   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
610811af0c4SBarry 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]
611811af0c4SBarry Smith   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
612811af0c4SBarry Smith .ve
6137e4fd573SVaclav Hapla 
614811af0c4SBarry 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.
6157e4fd573SVaclav Hapla 
6165c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
6175c6c1daeSBarry Smith 
618811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
619db781477SPatrick Sanan           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
620db781477SPatrick Sanan           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
6215c6c1daeSBarry Smith @*/
622d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
623d71ae5a4SJacob Faibussowitsch {
6245c6c1daeSBarry Smith   PetscFunctionBegin;
6259566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, hdf5v));
6269566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
6279566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
6289566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*hdf5v, name));
6299566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*hdf5v));
6305c6c1daeSBarry Smith   PetscFunctionReturn(0);
6315c6c1daeSBarry Smith }
6325c6c1daeSBarry Smith 
6335c6c1daeSBarry Smith /*@C
6345c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
6355c6c1daeSBarry Smith 
6365c6c1daeSBarry Smith   Not collective
6375c6c1daeSBarry Smith 
6385c6c1daeSBarry Smith   Input Parameter:
639811af0c4SBarry Smith . viewer - the `PetscViewer`
6405c6c1daeSBarry Smith 
6415c6c1daeSBarry Smith   Output Parameter:
6425c6c1daeSBarry Smith . file_id - The file id
6435c6c1daeSBarry Smith 
6445c6c1daeSBarry Smith   Level: intermediate
6455c6c1daeSBarry Smith 
646811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
6475c6c1daeSBarry Smith @*/
648d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
649d71ae5a4SJacob Faibussowitsch {
6505c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6515c6c1daeSBarry Smith 
6525c6c1daeSBarry Smith   PetscFunctionBegin;
6535c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
6545c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
6555c6c1daeSBarry Smith   PetscFunctionReturn(0);
6565c6c1daeSBarry Smith }
6575c6c1daeSBarry Smith 
6585c6c1daeSBarry Smith /*@C
6595c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
6605c6c1daeSBarry Smith 
6615c6c1daeSBarry Smith   Not collective
6625c6c1daeSBarry Smith 
6635c6c1daeSBarry Smith   Input Parameters:
664811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
6655c6c1daeSBarry Smith - name - The group name
6665c6c1daeSBarry Smith 
6675c6c1daeSBarry Smith   Level: intermediate
6685c6c1daeSBarry Smith 
6692e361344SVaclav Hapla   Notes:
6702e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
671811af0c4SBarry Smith .vb
672811af0c4SBarry 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.
673811af0c4SBarry Smith   NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
674811af0c4SBarry Smith   "." means the current group is pushed again.
675811af0c4SBarry Smith .ve
6762e361344SVaclav Hapla 
6772e361344SVaclav Hapla   Example:
6782e361344SVaclav Hapla   Suppose the current group is "/a".
6792e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
6802e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
6812e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
6822e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
6832e361344SVaclav Hapla 
684811af0c4SBarry Smith   Developer Note:
6852e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
686820fc2d1SVaclav Hapla 
68763cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
6885c6c1daeSBarry Smith @*/
689d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
690d71ae5a4SJacob Faibussowitsch {
6915c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6927d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6932e361344SVaclav Hapla   size_t                    i, len;
6942e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
6952e361344SVaclav Hapla   const char               *gname;
6965c6c1daeSBarry Smith 
6975c6c1daeSBarry Smith   PetscFunctionBegin;
6985c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
699820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
7009566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
7012e361344SVaclav Hapla   gname = NULL;
7022e361344SVaclav Hapla   if (len) {
7032e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
7042e361344SVaclav Hapla       /* use current name */
7052e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
7062e361344SVaclav Hapla     } else if (name[0] == '/') {
7072e361344SVaclav Hapla       /* absolute */
7082e361344SVaclav Hapla       for (i = 1; i < len; i++) {
7092e361344SVaclav Hapla         if (name[i] != '/') {
7102e361344SVaclav Hapla           gname = name;
7112e361344SVaclav Hapla           break;
7122e361344SVaclav Hapla         }
7132e361344SVaclav Hapla       }
7142e361344SVaclav Hapla     } else {
7152e361344SVaclav Hapla       /* relative */
7162e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
7179566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
7182e361344SVaclav Hapla       gname = buf;
7192e361344SVaclav Hapla     }
7202e361344SVaclav Hapla   }
7219566063dSJacob Faibussowitsch   PetscCall(PetscNew(&groupNode));
7229566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
7235c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
7245c6c1daeSBarry Smith   hdf5->groups    = groupNode;
7255c6c1daeSBarry Smith   PetscFunctionReturn(0);
7265c6c1daeSBarry Smith }
7275c6c1daeSBarry Smith 
7283ef9c667SSatish Balay /*@
7295c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
7305c6c1daeSBarry Smith 
7315c6c1daeSBarry Smith   Not collective
7325c6c1daeSBarry Smith 
7335c6c1daeSBarry Smith   Input Parameter:
734811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7355c6c1daeSBarry Smith 
7365c6c1daeSBarry Smith   Level: intermediate
7375c6c1daeSBarry Smith 
73863cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
7395c6c1daeSBarry Smith @*/
740d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
741d71ae5a4SJacob Faibussowitsch {
7425c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7437d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
7445c6c1daeSBarry Smith 
7455c6c1daeSBarry Smith   PetscFunctionBegin;
7465c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7475f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
7485c6c1daeSBarry Smith   groupNode    = hdf5->groups;
7495c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
7509566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
7519566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
7525c6c1daeSBarry Smith   PetscFunctionReturn(0);
7535c6c1daeSBarry Smith }
7545c6c1daeSBarry Smith 
755d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
756d71ae5a4SJacob Faibussowitsch {
7575c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7585c6c1daeSBarry Smith 
7595c6c1daeSBarry Smith   PetscFunctionBegin;
7605c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
761c959eef4SJed Brown   PetscValidPointer(name, 2);
762a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
7630298fd71SBarry Smith   else *name = NULL;
7645c6c1daeSBarry Smith   PetscFunctionReturn(0);
7655c6c1daeSBarry Smith }
7665c6c1daeSBarry Smith 
7673014b61aSVaclav Hapla /*@C
768811af0c4SBarry Smith   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
769874270d9SVaclav Hapla   and return this group's ID and file ID.
770811af0c4SBarry Smith   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
771874270d9SVaclav Hapla 
772874270d9SVaclav Hapla   Not collective
773874270d9SVaclav Hapla 
7743014b61aSVaclav Hapla   Input Parameters:
7753014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7763014b61aSVaclav Hapla - path - (Optional) The path relative to the pushed group
777874270d9SVaclav Hapla 
778d8d19677SJose E. Roman   Output Parameters:
779874270d9SVaclav Hapla + fileId - The HDF5 file ID
780874270d9SVaclav Hapla - groupId - The HDF5 group ID
781874270d9SVaclav Hapla 
782811af0c4SBarry Smith   Note:
7833014b61aSVaclav 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.
7843014b61aSVaclav Hapla   So NULL or empty path means the current pushed group.
7853014b61aSVaclav Hapla 
786e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
787e74616beSVaclav Hapla 
788874270d9SVaclav Hapla   Level: intermediate
789874270d9SVaclav Hapla 
79063cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
791874270d9SVaclav Hapla @*/
792d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
793d71ae5a4SJacob Faibussowitsch {
79490e03537SVaclav Hapla   hid_t       file_id;
79590e03537SVaclav Hapla   H5O_type_t  type;
7963014b61aSVaclav Hapla   const char *fileName  = NULL;
7973014b61aSVaclav Hapla   char       *groupName = NULL;
79876d59af2SVaclav Hapla   PetscBool   writable, has;
79954dbf706SVaclav Hapla 
80054dbf706SVaclav Hapla   PetscFunctionBegin;
8019566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
8029566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
8039566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
8043014b61aSVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
8059566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
80676d59af2SVaclav Hapla   if (!has) {
8075f80ce2aSJacob 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);
808f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
80976d59af2SVaclav Hapla   }
8105f80ce2aSJacob 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);
8113014b61aSVaclav Hapla   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
8123014b61aSVaclav Hapla   PetscCall(PetscFree(groupName));
81354dbf706SVaclav Hapla   *fileId = file_id;
81454dbf706SVaclav Hapla   PetscFunctionReturn(0);
81554dbf706SVaclav Hapla }
81654dbf706SVaclav Hapla 
81763cb69f5SVaclav Hapla /*@C
81863cb69f5SVaclav Hapla   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file
81963cb69f5SVaclav Hapla 
82063cb69f5SVaclav Hapla   Not collective
82163cb69f5SVaclav Hapla 
82263cb69f5SVaclav Hapla   Input Parameters:
82363cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
82463cb69f5SVaclav Hapla - path - (Optional) The path relative to the pushed group
82563cb69f5SVaclav Hapla 
82663cb69f5SVaclav Hapla   Note:
82763cb69f5SVaclav 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.
82863cb69f5SVaclav Hapla   So NULL or empty path means the current pushed group.
82963cb69f5SVaclav Hapla 
83063cb69f5SVaclav Hapla   This will fail if the viewer is not writable.
83163cb69f5SVaclav Hapla 
83263cb69f5SVaclav Hapla   Level: intermediate
83363cb69f5SVaclav Hapla 
83463cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
83563cb69f5SVaclav Hapla @*/
836d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
837d71ae5a4SJacob Faibussowitsch {
83863cb69f5SVaclav Hapla   hid_t fileId, groupId;
83963cb69f5SVaclav Hapla 
84063cb69f5SVaclav Hapla   PetscFunctionBegin;
84163cb69f5SVaclav Hapla   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
84263cb69f5SVaclav Hapla   PetscCallHDF5(H5Gclose, (groupId));
84363cb69f5SVaclav Hapla   PetscFunctionReturn(0);
84463cb69f5SVaclav Hapla }
84563cb69f5SVaclav Hapla 
8463ef9c667SSatish Balay /*@
847d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
8485c6c1daeSBarry Smith 
8495c6c1daeSBarry Smith   Not collective
8505c6c1daeSBarry Smith 
8515c6c1daeSBarry Smith   Input Parameter:
852811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
8535c6c1daeSBarry Smith 
8545c6c1daeSBarry Smith   Level: intermediate
8555c6c1daeSBarry Smith 
856d7dd068bSVaclav Hapla   Notes:
857811af0c4SBarry Smith   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
858811af0c4SBarry Smith   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
859811af0c4SBarry 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()`].
860811af0c4SBarry Smith   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
861811af0c4SBarry Smith   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
862d7dd068bSVaclav Hapla 
863d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
864d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
865d7dd068bSVaclav Hapla 
866811af0c4SBarry Smith   Developer note:
867d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
868d7dd068bSVaclav Hapla 
869811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
870d7dd068bSVaclav Hapla @*/
871d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
872d71ae5a4SJacob Faibussowitsch {
873d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
874d7dd068bSVaclav Hapla 
875d7dd068bSVaclav Hapla   PetscFunctionBegin;
876d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8775f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
878d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
879d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
880d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
881d7dd068bSVaclav Hapla }
882d7dd068bSVaclav Hapla 
883d7dd068bSVaclav Hapla /*@
884d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
885d7dd068bSVaclav Hapla 
886d7dd068bSVaclav Hapla   Not collective
887d7dd068bSVaclav Hapla 
888d7dd068bSVaclav Hapla   Input Parameter:
889811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
890d7dd068bSVaclav Hapla 
891d7dd068bSVaclav Hapla   Level: intermediate
892d7dd068bSVaclav Hapla 
893811af0c4SBarry Smith   Note:
894811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
895d7dd068bSVaclav Hapla 
896811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
897d7dd068bSVaclav Hapla @*/
898d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
899d71ae5a4SJacob Faibussowitsch {
900d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
901d7dd068bSVaclav Hapla 
902d7dd068bSVaclav Hapla   PetscFunctionBegin;
903d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9045f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
905d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
906d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
907d7dd068bSVaclav Hapla }
908d7dd068bSVaclav Hapla 
909d7dd068bSVaclav Hapla /*@
910d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
911d7dd068bSVaclav Hapla 
912d7dd068bSVaclav Hapla   Not collective
913d7dd068bSVaclav Hapla 
914d7dd068bSVaclav Hapla   Input Parameter:
915811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
916d7dd068bSVaclav Hapla 
917d7dd068bSVaclav Hapla   Output Parameter:
918d7dd068bSVaclav Hapla . flg - is timestepping active?
919d7dd068bSVaclav Hapla 
920d7dd068bSVaclav Hapla   Level: intermediate
921d7dd068bSVaclav Hapla 
922811af0c4SBarry Smith   Note:
923811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
924d7dd068bSVaclav Hapla 
925811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
926d7dd068bSVaclav Hapla @*/
927d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
928d71ae5a4SJacob Faibussowitsch {
929d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
930d7dd068bSVaclav Hapla 
931d7dd068bSVaclav Hapla   PetscFunctionBegin;
932d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
933d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
934d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
935d7dd068bSVaclav Hapla }
936d7dd068bSVaclav Hapla 
937d7dd068bSVaclav Hapla /*@
938d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
939d7dd068bSVaclav Hapla 
940d7dd068bSVaclav Hapla   Not collective
941d7dd068bSVaclav Hapla 
942d7dd068bSVaclav Hapla   Input Parameter:
943811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
944d7dd068bSVaclav Hapla 
945d7dd068bSVaclav Hapla   Level: intermediate
946d7dd068bSVaclav Hapla 
947811af0c4SBarry Smith   Note:
948811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
949d7dd068bSVaclav Hapla 
950811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
9515c6c1daeSBarry Smith @*/
952d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
953d71ae5a4SJacob Faibussowitsch {
9545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9555c6c1daeSBarry Smith 
9565c6c1daeSBarry Smith   PetscFunctionBegin;
9575c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9585f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9595c6c1daeSBarry Smith   ++hdf5->timestep;
9605c6c1daeSBarry Smith   PetscFunctionReturn(0);
9615c6c1daeSBarry Smith }
9625c6c1daeSBarry Smith 
9633ef9c667SSatish Balay /*@
964d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
9655c6c1daeSBarry Smith 
966d7dd068bSVaclav Hapla   Logically collective
9675c6c1daeSBarry Smith 
9685c6c1daeSBarry Smith   Input Parameters:
969811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
970d7dd068bSVaclav Hapla - timestep - The timestep
9715c6c1daeSBarry Smith 
9725c6c1daeSBarry Smith   Level: intermediate
9735c6c1daeSBarry Smith 
974811af0c4SBarry Smith   Note:
975811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
976d7dd068bSVaclav Hapla 
977811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
9785c6c1daeSBarry Smith @*/
979d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
980d71ae5a4SJacob Faibussowitsch {
9815c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9825c6c1daeSBarry Smith 
9835c6c1daeSBarry Smith   PetscFunctionBegin;
9845c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
985d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
9865f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
9875f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9885c6c1daeSBarry Smith   hdf5->timestep = timestep;
9895c6c1daeSBarry Smith   PetscFunctionReturn(0);
9905c6c1daeSBarry Smith }
9915c6c1daeSBarry Smith 
9923ef9c667SSatish Balay /*@
9935c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
9945c6c1daeSBarry Smith 
9955c6c1daeSBarry Smith   Not collective
9965c6c1daeSBarry Smith 
9975c6c1daeSBarry Smith   Input Parameter:
998811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
9995c6c1daeSBarry Smith 
10005c6c1daeSBarry Smith   Output Parameter:
1001d7dd068bSVaclav Hapla . timestep - The timestep
10025c6c1daeSBarry Smith 
10035c6c1daeSBarry Smith   Level: intermediate
10045c6c1daeSBarry Smith 
1005811af0c4SBarry Smith   Note:
1006811af0c4SBarry Smith   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1007d7dd068bSVaclav Hapla 
1008811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
10095c6c1daeSBarry Smith @*/
1010d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
1011d71ae5a4SJacob Faibussowitsch {
10125c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
10135c6c1daeSBarry Smith 
10145c6c1daeSBarry Smith   PetscFunctionBegin;
10155c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1016d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep, 2);
10175f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
10185c6c1daeSBarry Smith   *timestep = hdf5->timestep;
10195c6c1daeSBarry Smith   PetscFunctionReturn(0);
10205c6c1daeSBarry Smith }
10215c6c1daeSBarry Smith 
102236bce6e8SMatthew G. Knepley /*@C
102336bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
102436bce6e8SMatthew G. Knepley 
102536bce6e8SMatthew G. Knepley   Not collective
102636bce6e8SMatthew G. Knepley 
102736bce6e8SMatthew G. Knepley   Input Parameter:
1028811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
102936bce6e8SMatthew G. Knepley 
103036bce6e8SMatthew G. Knepley   Output Parameter:
1031811af0c4SBarry Smith . mtype - the HDF5  datatype
103236bce6e8SMatthew G. Knepley 
103336bce6e8SMatthew G. Knepley   Level: advanced
103436bce6e8SMatthew G. Knepley 
1035db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
103636bce6e8SMatthew G. Knepley @*/
1037d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
1038d71ae5a4SJacob Faibussowitsch {
103936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
104036bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
104136bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
104236bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_LLONG;
104336bce6e8SMatthew G. Knepley #else
104436bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_INT;
104536bce6e8SMatthew G. Knepley #endif
104636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
104736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
104836bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1049c9450970SBarry Smith   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
1050de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
105136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
105236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
105336bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
10547e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
105536bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
105636bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
105736bce6e8SMatthew G. Knepley }
105836bce6e8SMatthew G. Knepley 
105936bce6e8SMatthew G. Knepley /*@C
106036bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
106136bce6e8SMatthew G. Knepley 
106236bce6e8SMatthew G. Knepley   Not collective
106336bce6e8SMatthew G. Knepley 
106436bce6e8SMatthew G. Knepley   Input Parameter:
1065811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
106636bce6e8SMatthew G. Knepley 
106736bce6e8SMatthew G. Knepley   Output Parameter:
1068811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
106936bce6e8SMatthew G. Knepley 
107036bce6e8SMatthew G. Knepley   Level: advanced
107136bce6e8SMatthew G. Knepley 
1072db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
107336bce6e8SMatthew G. Knepley @*/
1074d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
1075d71ae5a4SJacob Faibussowitsch {
107636bce6e8SMatthew G. Knepley   PetscFunctionBegin;
107736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
107836bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
107936bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
108036bce6e8SMatthew G. Knepley #else
108136bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
108236bce6e8SMatthew G. Knepley #endif
108336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
108436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
108536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
108636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
108736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
108836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
10897e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
109036bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
109136bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
109236bce6e8SMatthew G. Knepley }
109336bce6e8SMatthew G. Knepley 
1094df863907SAlex Fikl /*@C
1095b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
109636bce6e8SMatthew G. Knepley 
1097343a20b2SVaclav Hapla   Collective
1098343a20b2SVaclav Hapla 
109936bce6e8SMatthew G. Knepley   Input Parameters:
1100811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
11014302210dSVaclav Hapla . parent - The parent dataset/group name
110236bce6e8SMatthew G. Knepley . name   - The attribute name
110336bce6e8SMatthew G. Knepley . datatype - The attribute type
110436bce6e8SMatthew G. Knepley - value    - The attribute value
110536bce6e8SMatthew G. Knepley 
110636bce6e8SMatthew G. Knepley   Level: advanced
110736bce6e8SMatthew G. Knepley 
1108811af0c4SBarry Smith   Note:
1109343a20b2SVaclav 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.
11104302210dSVaclav Hapla 
1111811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1112811af0c4SBarry Smith           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
111336bce6e8SMatthew G. Knepley @*/
1114d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1115d71ae5a4SJacob Faibussowitsch {
11164302210dSVaclav Hapla   char     *parentAbsPath;
111760568592SMatthew G. Knepley   hid_t     h5, dataspace, obj, attribute, dtype;
1118080f144cSVaclav Hapla   PetscBool has;
111936bce6e8SMatthew G. Knepley 
112036bce6e8SMatthew G. Knepley   PetscFunctionBegin;
11215cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
11224302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1123c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
11244302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1125b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
112677717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
11279566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
11289566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
11299566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11307e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
11317e97c476SMatthew G. Knepley     size_t len;
11329566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *)value, &len));
1133792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
11347e97c476SMatthew G. Knepley   }
11359566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1136792fecdfSBarry Smith   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1137792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1138080f144cSVaclav Hapla   if (has) {
1139792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1140080f144cSVaclav Hapla   } else {
1141792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1142080f144cSVaclav Hapla   }
1143792fecdfSBarry Smith   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1144792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1145792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1146792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
1147792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (dataspace));
11489566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
114936bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
115036bce6e8SMatthew G. Knepley }
115136bce6e8SMatthew G. Knepley 
1152df863907SAlex Fikl /*@C
1153811af0c4SBarry Smith  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1154577e0e04SVaclav Hapla 
1155343a20b2SVaclav Hapla   Collective
1156343a20b2SVaclav Hapla 
1157577e0e04SVaclav Hapla   Input Parameters:
1158811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1159577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1160577e0e04SVaclav Hapla . name     - The attribute name
1161577e0e04SVaclav Hapla . datatype - The attribute type
1162577e0e04SVaclav Hapla - value    - The attribute value
1163577e0e04SVaclav Hapla 
1164811af0c4SBarry Smith   Note:
1165343a20b2SVaclav 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).
1166811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1167577e0e04SVaclav Hapla 
1168577e0e04SVaclav Hapla   Level: advanced
1169577e0e04SVaclav Hapla 
1170811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1171811af0c4SBarry Smith           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1172577e0e04SVaclav Hapla @*/
1173d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1174d71ae5a4SJacob Faibussowitsch {
1175577e0e04SVaclav Hapla   PetscFunctionBegin;
1176577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1177577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1178577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1179b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
11809566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
11819566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1182577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1183577e0e04SVaclav Hapla }
1184577e0e04SVaclav Hapla 
1185577e0e04SVaclav Hapla /*@C
1186b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1187ad0c61feSMatthew G. Knepley 
1188343a20b2SVaclav Hapla   Collective
1189343a20b2SVaclav Hapla 
1190ad0c61feSMatthew G. Knepley   Input Parameters:
1191811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
11924302210dSVaclav Hapla . parent - The parent dataset/group name
1193ad0c61feSMatthew G. Knepley . name   - The attribute name
1194a2d6be1bSVaclav Hapla . datatype - The attribute type
1195a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1196ad0c61feSMatthew G. Knepley 
1197ad0c61feSMatthew G. Knepley   Output Parameter:
1198a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1199ad0c61feSMatthew G. Knepley 
1200a2d6be1bSVaclav Hapla   Notes:
1201a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1202a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1203a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1204811af0c4SBarry Smith $  flg = `PETSC_FALSE`;
1205811af0c4SBarry Smith $  PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1206a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1207a2d6be1bSVaclav Hapla 
1208811af0c4SBarry Smith   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1209708d4cb3SBarry Smith 
1210343a20b2SVaclav 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.
1211ad0c61feSMatthew G. Knepley 
1212343a20b2SVaclav Hapla   Level: advanced
12134302210dSVaclav Hapla 
1214811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1215ad0c61feSMatthew G. Knepley @*/
1216d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1217d71ae5a4SJacob Faibussowitsch {
12184302210dSVaclav Hapla   char     *parentAbsPath;
1219a2d6be1bSVaclav Hapla   hid_t     h5, obj, attribute, dtype;
1220969748fdSVaclav Hapla   PetscBool has;
1221ad0c61feSMatthew G. Knepley 
1222ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
12235cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12244302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1225c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1226a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1227a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
12289566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
122977717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
12309566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
12319566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1232a2d6be1bSVaclav Hapla   if (!has) {
1233a2d6be1bSVaclav Hapla     if (defaultValue) {
1234a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1235a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
12369566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1237a2d6be1bSVaclav Hapla         } else {
1238a2d6be1bSVaclav Hapla           size_t len;
1239792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
12409566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1241a2d6be1bSVaclav Hapla         }
1242a2d6be1bSVaclav Hapla       }
12439566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
1244a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
124598921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1246a2d6be1bSVaclav Hapla   }
12479566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1248792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1249792fecdfSBarry Smith   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1250f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1251f0b58479SMatthew G. Knepley     size_t len;
1252a2d6be1bSVaclav Hapla     hid_t  atype;
1253792fecdfSBarry Smith     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1254792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
12559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1256792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1257792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1258708d4cb3SBarry Smith   } else {
1259792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1260708d4cb3SBarry Smith   }
1261792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1262e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1263792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
12649566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
1265ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1266ad0c61feSMatthew G. Knepley }
1267ad0c61feSMatthew G. Knepley 
1268577e0e04SVaclav Hapla /*@C
1269811af0c4SBarry Smith  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1270577e0e04SVaclav Hapla 
1271343a20b2SVaclav Hapla   Collective
1272343a20b2SVaclav Hapla 
1273577e0e04SVaclav Hapla   Input Parameters:
1274811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1275577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1276577e0e04SVaclav Hapla . name     - The attribute name
1277577e0e04SVaclav Hapla - datatype - The attribute type
1278577e0e04SVaclav Hapla 
1279577e0e04SVaclav Hapla   Output Parameter:
1280577e0e04SVaclav Hapla . value    - The attribute value
1281577e0e04SVaclav Hapla 
1282811af0c4SBarry Smith   Note:
1283577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1284811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1285577e0e04SVaclav Hapla 
1286577e0e04SVaclav Hapla   Level: advanced
1287577e0e04SVaclav Hapla 
1288811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1289577e0e04SVaclav Hapla @*/
1290d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1291d71ae5a4SJacob Faibussowitsch {
1292577e0e04SVaclav Hapla   PetscFunctionBegin;
1293577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1294577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1295577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1296064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
12979566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12989566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1299577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1300577e0e04SVaclav Hapla }
1301577e0e04SVaclav Hapla 
1302d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1303d71ae5a4SJacob Faibussowitsch {
1304a75c3fd4SVaclav Hapla   htri_t exists;
1305a75c3fd4SVaclav Hapla   hid_t  group;
1306a75c3fd4SVaclav Hapla 
1307a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1308792fecdfSBarry Smith   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1309792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1310a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1311792fecdfSBarry Smith     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1312792fecdfSBarry Smith     PetscCallHDF5(H5Gclose, (group));
1313a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1314a75c3fd4SVaclav Hapla   }
1315a75c3fd4SVaclav Hapla   *exists_ = (PetscBool)exists;
1316a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1317a75c3fd4SVaclav Hapla }
1318a75c3fd4SVaclav Hapla 
1319d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1320d71ae5a4SJacob Faibussowitsch {
132190e03537SVaclav Hapla   const char rootGroupName[] = "/";
13225cdeb108SMatthew G. Knepley   hid_t      h5;
1323e5bf9ebcSVaclav Hapla   PetscBool  exists = PETSC_FALSE;
132415b861d2SVaclav Hapla   PetscInt   i;
132515b861d2SVaclav Hapla   int        n;
132685688ae2SVaclav Hapla   char     **hierarchy;
132785688ae2SVaclav Hapla   char       buf[PETSC_MAX_PATH_LEN] = "";
13285cdeb108SMatthew G. Knepley 
13295cdeb108SMatthew G. Knepley   PetscFunctionBegin;
13305cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
133190e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
133290e03537SVaclav Hapla   else name = rootGroupName;
1333ccd66a83SVaclav Hapla   if (has) {
1334064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
13355cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1336ccd66a83SVaclav Hapla   }
1337ccd66a83SVaclav Hapla   if (otype) {
1338064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
133956cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1340ccd66a83SVaclav Hapla   }
13419566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
134285688ae2SVaclav Hapla 
134385688ae2SVaclav Hapla   /*
134485688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
134585688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
134685688ae2SVaclav Hapla      1) whether it's a valid link
134785688ae2SVaclav Hapla      2) whether this link resolves to an object
134885688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
134985688ae2SVaclav Hapla   */
13509566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
135185688ae2SVaclav Hapla   if (!n) {
135285688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1353ccd66a83SVaclav Hapla     if (has) *has = PETSC_TRUE;
1354ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
13559566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
135685688ae2SVaclav Hapla     PetscFunctionReturn(0);
135785688ae2SVaclav Hapla   }
135885688ae2SVaclav Hapla   for (i = 0; i < n; i++) {
13599566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
13609566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, hierarchy[i]));
13619566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1362a75c3fd4SVaclav Hapla     if (!exists) break;
136385688ae2SVaclav Hapla   }
13649566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
136585688ae2SVaclav Hapla 
136685688ae2SVaclav Hapla   /* If the object exists, get its type */
1367ccd66a83SVaclav Hapla   if (exists && otype) {
13685cdeb108SMatthew G. Knepley     H5O_info_t info;
13695cdeb108SMatthew G. Knepley 
137076276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1371792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
137256cc0592SVaclav Hapla     *otype = info.type;
13735cdeb108SMatthew G. Knepley   }
1374ccd66a83SVaclav Hapla   if (has) *has = exists;
13755cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
13765cdeb108SMatthew G. Knepley }
13775cdeb108SMatthew G. Knepley 
13784302210dSVaclav Hapla /*@C
13790a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
13800a9f38efSVaclav Hapla 
1381343a20b2SVaclav Hapla   Collective
1382343a20b2SVaclav Hapla 
13830a9f38efSVaclav Hapla   Input Parameters:
1384811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1385a0558868SVaclav Hapla - path - (Optional) The path relative to the pushed group
13860a9f38efSVaclav Hapla 
13870a9f38efSVaclav Hapla   Output Parameter:
13880a9f38efSVaclav Hapla . has - Flag for group existence
13890a9f38efSVaclav Hapla 
13900a9f38efSVaclav Hapla   Level: advanced
13910a9f38efSVaclav Hapla 
13924302210dSVaclav Hapla   Notes:
1393785c443eSVaclav 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.
1394785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
13954302210dSVaclav Hapla 
1396811af0c4SBarry Smith   If path exists but is not a group, `PETSC_FALSE` is returned.
1397811af0c4SBarry Smith 
1398811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
13990a9f38efSVaclav Hapla @*/
1400d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1401d71ae5a4SJacob Faibussowitsch {
14020a9f38efSVaclav Hapla   H5O_type_t type;
14034302210dSVaclav Hapla   char      *abspath;
14040a9f38efSVaclav Hapla 
14050a9f38efSVaclav Hapla   PetscFunctionBegin;
14060a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14074302210dSVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
14084302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
140977717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
14109566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
14114302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
14129566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
14130a9f38efSVaclav Hapla   PetscFunctionReturn(0);
14140a9f38efSVaclav Hapla }
14150a9f38efSVaclav Hapla 
141689e0ef10SVaclav Hapla /*@C
141789e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
141889e0ef10SVaclav Hapla 
1419343a20b2SVaclav Hapla   Collective
1420343a20b2SVaclav Hapla 
142189e0ef10SVaclav Hapla   Input Parameters:
1422811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
142389e0ef10SVaclav Hapla - path - The dataset path
142489e0ef10SVaclav Hapla 
142589e0ef10SVaclav Hapla   Output Parameter:
142689e0ef10SVaclav Hapla . has - Flag whether dataset exists
142789e0ef10SVaclav Hapla 
142889e0ef10SVaclav Hapla   Level: advanced
142989e0ef10SVaclav Hapla 
143089e0ef10SVaclav Hapla   Notes:
1431343a20b2SVaclav 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.
143289e0ef10SVaclav Hapla 
1433811af0c4SBarry Smith   If path is NULL or empty, has is set to `PETSC_FALSE`.
1434811af0c4SBarry Smith 
1435811af0c4SBarry Smith   If path exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1436811af0c4SBarry Smith 
1437811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
143889e0ef10SVaclav Hapla @*/
1439d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1440d71ae5a4SJacob Faibussowitsch {
144189e0ef10SVaclav Hapla   H5O_type_t type;
144289e0ef10SVaclav Hapla   char      *abspath;
144389e0ef10SVaclav Hapla 
144489e0ef10SVaclav Hapla   PetscFunctionBegin;
144589e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
144689e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
144789e0ef10SVaclav Hapla   PetscValidBoolPointer(has, 3);
144877717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
14499566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
145089e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
14519566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
145289e0ef10SVaclav Hapla   PetscFunctionReturn(0);
145389e0ef10SVaclav Hapla }
145489e0ef10SVaclav Hapla 
14550a9f38efSVaclav Hapla /*@
1456e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1457ecce1506SVaclav Hapla 
1458343a20b2SVaclav Hapla   Collective
1459343a20b2SVaclav Hapla 
1460ecce1506SVaclav Hapla   Input Parameters:
1461811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1462ecce1506SVaclav Hapla - obj    - The named object
1463ecce1506SVaclav Hapla 
1464ecce1506SVaclav Hapla   Output Parameter:
146589e0ef10SVaclav Hapla . has    - Flag for dataset existence
1466ecce1506SVaclav Hapla 
1467e3f143f7SVaclav Hapla   Notes:
146889e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1469811af0c4SBarry Smith 
1470811af0c4SBarry Smith   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1471e3f143f7SVaclav Hapla 
1472ecce1506SVaclav Hapla   Level: advanced
1473ecce1506SVaclav Hapla 
1474811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1475ecce1506SVaclav Hapla @*/
1476d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1477d71ae5a4SJacob Faibussowitsch {
147889e0ef10SVaclav Hapla   size_t len;
1479ecce1506SVaclav Hapla 
1480ecce1506SVaclav Hapla   PetscFunctionBegin;
1481c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1482c57a1a9aSVaclav Hapla   PetscValidHeader(obj, 2);
14834302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
14849566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
14855f80ce2aSJacob Faibussowitsch   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
14869566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1487ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1488ecce1506SVaclav Hapla }
1489ecce1506SVaclav Hapla 
1490df863907SAlex Fikl /*@C
1491ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1492e7bdbf8eSMatthew G. Knepley 
1493343a20b2SVaclav Hapla   Collective
1494343a20b2SVaclav Hapla 
1495e7bdbf8eSMatthew G. Knepley   Input Parameters:
1496811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
14974302210dSVaclav Hapla . parent - The parent dataset/group name
1498e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1499e7bdbf8eSMatthew G. Knepley 
1500e7bdbf8eSMatthew G. Knepley   Output Parameter:
1501e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1502e7bdbf8eSMatthew G. Knepley 
1503e7bdbf8eSMatthew G. Knepley   Level: advanced
1504e7bdbf8eSMatthew G. Knepley 
1505811af0c4SBarry Smith   Note:
1506343a20b2SVaclav 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.
15074302210dSVaclav Hapla 
1508811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1509e7bdbf8eSMatthew G. Knepley @*/
1510d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1511d71ae5a4SJacob Faibussowitsch {
15124302210dSVaclav Hapla   char *parentAbsPath;
1513e7bdbf8eSMatthew G. Knepley 
1514e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
15155cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
15164302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1517c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
15184302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
151977717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
15209566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
15219566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
15229566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
152306db490cSVaclav Hapla   PetscFunctionReturn(0);
152406db490cSVaclav Hapla }
152506db490cSVaclav Hapla 
1526577e0e04SVaclav Hapla /*@C
1527811af0c4SBarry Smith  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1528577e0e04SVaclav Hapla 
1529343a20b2SVaclav Hapla   Collective
1530343a20b2SVaclav Hapla 
1531577e0e04SVaclav Hapla   Input Parameters:
1532811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1533577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1534577e0e04SVaclav Hapla - name   - The attribute name
1535577e0e04SVaclav Hapla 
1536577e0e04SVaclav Hapla   Output Parameter:
1537577e0e04SVaclav Hapla . has    - Flag for attribute existence
1538577e0e04SVaclav Hapla 
1539811af0c4SBarry Smith   Note:
1540577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1541811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1542577e0e04SVaclav Hapla 
1543577e0e04SVaclav Hapla   Level: advanced
1544577e0e04SVaclav Hapla 
1545811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1546577e0e04SVaclav Hapla @*/
1547d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1548d71ae5a4SJacob Faibussowitsch {
1549577e0e04SVaclav Hapla   PetscFunctionBegin;
1550577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1551577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1552577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
15534302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
15549566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
15559566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1556577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1557577e0e04SVaclav Hapla }
1558577e0e04SVaclav Hapla 
1559d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1560d71ae5a4SJacob Faibussowitsch {
156106db490cSVaclav Hapla   hid_t  h5;
156206db490cSVaclav Hapla   htri_t hhas;
156306db490cSVaclav Hapla 
156406db490cSVaclav Hapla   PetscFunctionBegin;
15659566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1566792fecdfSBarry Smith   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
15675cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1568e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1569e7bdbf8eSMatthew G. Knepley }
1570e7bdbf8eSMatthew G. Knepley 
1571a75e6a4aSMatthew G. Knepley /*
1572a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1573a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1574a75e6a4aSMatthew G. Knepley */
1575d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1576a75e6a4aSMatthew G. Knepley 
1577a75e6a4aSMatthew G. Knepley /*@C
1578811af0c4SBarry Smith   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1579a75e6a4aSMatthew G. Knepley 
1580d083f849SBarry Smith   Collective
1581a75e6a4aSMatthew G. Knepley 
1582a75e6a4aSMatthew G. Knepley   Input Parameter:
1583811af0c4SBarry Smith . comm - the MPI communicator to share the HDF5 `PetscViewer`
1584a75e6a4aSMatthew G. Knepley 
1585a75e6a4aSMatthew G. Knepley   Level: intermediate
1586a75e6a4aSMatthew G. Knepley 
1587811af0c4SBarry Smith   Options Database Key:
158810699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1589a75e6a4aSMatthew G. Knepley 
1590811af0c4SBarry Smith   Environmental variable:
1591811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1592a75e6a4aSMatthew G. Knepley 
1593811af0c4SBarry Smith   Note:
1594811af0c4SBarry Smith   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1595811af0c4SBarry Smith   an error code.  The HDF5 `PetscViewer` is usually used in the form
1596a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1597a75e6a4aSMatthew G. Knepley 
1598811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1599a75e6a4aSMatthew G. Knepley @*/
1600d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1601d71ae5a4SJacob Faibussowitsch {
1602a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1603a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1604a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1605a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1606a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1607a75e6a4aSMatthew G. Knepley 
1608a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
16099371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16109371c9d4SSatish Balay   if (ierr) {
16119371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16129371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16139371c9d4SSatish Balay   }
1614a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1615908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
16169371c9d4SSatish Balay     if (ierr) {
16179371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16189371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16199371c9d4SSatish Balay     }
1620a75e6a4aSMatthew G. Knepley   }
162147435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
16229371c9d4SSatish Balay   if (ierr) {
16239371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16249371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16259371c9d4SSatish Balay   }
1626a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1627a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16289371c9d4SSatish Balay     if (ierr) {
16299371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16309371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16319371c9d4SSatish Balay     }
1632a75e6a4aSMatthew G. Knepley     if (!flg) {
1633a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname, "output.h5");
16349371c9d4SSatish Balay       if (ierr) {
16359371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16369371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16379371c9d4SSatish Balay       }
1638a75e6a4aSMatthew G. Knepley     }
1639a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
16409371c9d4SSatish Balay     if (ierr) {
16419371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16429371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16439371c9d4SSatish Balay     }
1644a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16459371c9d4SSatish Balay     if (ierr) {
16469371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16479371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16489371c9d4SSatish Balay     }
164947435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
16509371c9d4SSatish Balay     if (ierr) {
16519371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16529371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16539371c9d4SSatish Balay     }
1654a75e6a4aSMatthew G. Knepley   }
1655a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
16569371c9d4SSatish Balay   if (ierr) {
16579371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16589371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16599371c9d4SSatish Balay   }
1660a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1661a75e6a4aSMatthew G. Knepley }
1662