xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision d71ae5a4db6382e7f06317b8d368875286fe9008)
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 @*/
28*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[])
29*d71ae5a4SJacob 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 
53*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
54*d71ae5a4SJacob 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 
67*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject)
68*d71ae5a4SJacob 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 
85*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
86*d71ae5a4SJacob 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 
100*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
101*d71ae5a4SJacob 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 
110*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
111*d71ae5a4SJacob 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 
119*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
120*d71ae5a4SJacob 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 
147*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
148*d71ae5a4SJacob 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 
156*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
157*d71ae5a4SJacob 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 
165*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
166*d71ae5a4SJacob 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 
178811af0c4SBarry Smith     Logically Collective on viewer
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 @*/
195*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
196*d71ae5a4SJacob 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 
207811af0c4SBarry Smith     Logically Collective on viewer
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 @*/
223*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
224*d71ae5a4SJacob 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 
233*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
234*d71ae5a4SJacob 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 
246811af0c4SBarry Smith     Logically Collective on viewer
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 @*/
264*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
265*d71ae5a4SJacob 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 
276811af0c4SBarry Smith     Logically Collective on viewer
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 @*/
289*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
290*d71ae5a4SJacob 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 
299*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
300*d71ae5a4SJacob 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 @*/
340*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
341*d71ae5a4SJacob 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 
349*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
350*d71ae5a4SJacob 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 @*/
385*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
386*d71ae5a4SJacob 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 
395*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
396*d71ae5a4SJacob 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   }
433*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
434*d71ae5a4SJacob Faibussowitsch     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
435*d71ae5a4SJacob Faibussowitsch     break;
436*d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
437*d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
438*d71ae5a4SJacob Faibussowitsch   default:
439*d71ae5a4SJacob 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));
4435c6c1daeSBarry Smith   PetscFunctionReturn(0);
4445c6c1daeSBarry Smith }
4455c6c1daeSBarry Smith 
446*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
447*d71ae5a4SJacob Faibussowitsch {
448d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
449d1232d7fSToby Isaac 
450d1232d7fSToby Isaac   PetscFunctionBegin;
451d1232d7fSToby Isaac   *name = vhdf5->filename;
452d1232d7fSToby Isaac   PetscFunctionReturn(0);
453d1232d7fSToby Isaac }
454d1232d7fSToby Isaac 
455*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
456*d71ae5a4SJacob Faibussowitsch {
4575dc64a97SVaclav Hapla   /*
458b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4595dc64a97SVaclav Hapla   */
460b723ab35SVaclav Hapla 
461b723ab35SVaclav Hapla   PetscFunctionBegin;
462b723ab35SVaclav Hapla   PetscFunctionReturn(0);
463b723ab35SVaclav Hapla }
464b723ab35SVaclav Hapla 
465*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
466*d71ae5a4SJacob Faibussowitsch {
46719a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
46819a20e4cSMatthew G. Knepley 
46919a20e4cSMatthew G. Knepley   PetscFunctionBegin;
47019a20e4cSMatthew G. Knepley   hdf5->defTimestepping = flg;
47119a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
47219a20e4cSMatthew G. Knepley }
47319a20e4cSMatthew G. Knepley 
474*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
475*d71ae5a4SJacob Faibussowitsch {
47619a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
47719a20e4cSMatthew G. Knepley 
47819a20e4cSMatthew G. Knepley   PetscFunctionBegin;
47919a20e4cSMatthew G. Knepley   *flg = hdf5->defTimestepping;
48019a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
48119a20e4cSMatthew G. Knepley }
48219a20e4cSMatthew G. Knepley 
48319a20e4cSMatthew G. Knepley /*@
48419a20e4cSMatthew G. Knepley   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
48519a20e4cSMatthew G. Knepley 
486811af0c4SBarry Smith   Logically Collective on viewer
48719a20e4cSMatthew G. Knepley 
48819a20e4cSMatthew G. Knepley   Input Parameters:
489811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
490811af0c4SBarry Smith - flg    - if `PETSC_TRUE` we will assume that timestepping is on
49119a20e4cSMatthew G. Knepley 
492811af0c4SBarry Smith   Options Database Key:
49319a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
49419a20e4cSMatthew G. Knepley 
495811af0c4SBarry Smith   Note:
49619a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
49719a20e4cSMatthew G. Knepley 
49819a20e4cSMatthew G. Knepley   Level: intermediate
49919a20e4cSMatthew G. Knepley 
500811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
50119a20e4cSMatthew G. Knepley @*/
502*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
503*d71ae5a4SJacob Faibussowitsch {
50419a20e4cSMatthew G. Knepley   PetscFunctionBegin;
50519a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
506cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
50719a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
50819a20e4cSMatthew G. Knepley }
50919a20e4cSMatthew G. Knepley 
51019a20e4cSMatthew G. Knepley /*@
51119a20e4cSMatthew G. Knepley   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
51219a20e4cSMatthew G. Knepley 
51319a20e4cSMatthew G. Knepley   Not collective
51419a20e4cSMatthew G. Knepley 
51519a20e4cSMatthew G. Knepley   Input Parameter:
516811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
51719a20e4cSMatthew G. Knepley 
51819a20e4cSMatthew G. Knepley   Output Parameter:
519811af0c4SBarry Smith . flg    - if `PETSC_TRUE` we will assume that timestepping is on
52019a20e4cSMatthew G. Knepley 
52119a20e4cSMatthew G. Knepley   Level: intermediate
52219a20e4cSMatthew G. Knepley 
523811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
52419a20e4cSMatthew G. Knepley @*/
525*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
526*d71ae5a4SJacob Faibussowitsch {
52719a20e4cSMatthew G. Knepley   PetscFunctionBegin;
52819a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
529cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
53019a20e4cSMatthew G. Knepley   PetscFunctionReturn(0);
53119a20e4cSMatthew G. Knepley }
53219a20e4cSMatthew G. Knepley 
5338556b5ebSBarry Smith /*MC
5348556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
5358556b5ebSBarry Smith 
536811af0c4SBarry Smith   Level: beginner
537811af0c4SBarry Smith 
538db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
539db781477SPatrick Sanan           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
540db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
541db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
5428556b5ebSBarry Smith M*/
543d1232d7fSToby Isaac 
544*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
545*d71ae5a4SJacob Faibussowitsch {
5465c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
5475c6c1daeSBarry Smith 
5485c6c1daeSBarry Smith   PetscFunctionBegin;
54999335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
55099335c34SVaclav Hapla   {
55199335c34SVaclav Hapla     PetscMPIInt size;
5529566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
5535f80ce2aSJacob 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)");
55499335c34SVaclav Hapla   }
55599335c34SVaclav Hapla #endif
55699335c34SVaclav Hapla 
5574dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hdf5));
5585c6c1daeSBarry Smith 
5595c6c1daeSBarry Smith   v->data                = (void *)hdf5;
5605c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
56182ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
562b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
5631b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
5646226335aSVaclav Hapla   v->ops->flush          = PetscViewerFlush_HDF5;
5657e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
566908793a3SLisandro Dalcin   hdf5->filename         = NULL;
5675c6c1daeSBarry Smith   hdf5->timestep         = -1;
5680298fd71SBarry Smith   hdf5->groups           = NULL;
5695c6c1daeSBarry Smith 
570792fecdfSBarry Smith   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
5719c5072fbSVaclav Hapla 
5729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
5739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
5749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
5759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
5769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
5779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
5789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
5799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
5809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
5819566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
5825c6c1daeSBarry Smith   PetscFunctionReturn(0);
5835c6c1daeSBarry Smith }
5845c6c1daeSBarry Smith 
5855c6c1daeSBarry Smith /*@C
586811af0c4SBarry Smith    PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`
5875c6c1daeSBarry Smith 
588d083f849SBarry Smith    Collective
5895c6c1daeSBarry Smith 
5905c6c1daeSBarry Smith    Input Parameters:
5915c6c1daeSBarry Smith +  comm - MPI communicator
5925c6c1daeSBarry Smith .  name - name of file
5935c6c1daeSBarry Smith -  type - type of file
5945c6c1daeSBarry Smith 
5955c6c1daeSBarry Smith    Output Parameter:
596811af0c4SBarry Smith .  hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
5975c6c1daeSBarry Smith 
598811af0c4SBarry Smith   Options Database Keys:
599a2b725a8SWilliam 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
600a2b725a8SWilliam Gropp -  -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
60182ea9b62SBarry Smith 
6025c6c1daeSBarry Smith    Level: beginner
6035c6c1daeSBarry Smith 
6047e4fd573SVaclav Hapla    Notes:
6057e4fd573SVaclav Hapla    Reading is always available, regardless of the mode. Available modes are
606811af0c4SBarry Smith .vb
607811af0c4SBarry Smith   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
608811af0c4SBarry Smith   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
609811af0c4SBarry 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]
610811af0c4SBarry Smith   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
611811af0c4SBarry Smith .ve
6127e4fd573SVaclav Hapla 
613811af0c4SBarry 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.
6147e4fd573SVaclav Hapla 
6155c6c1daeSBarry Smith    This PetscViewer should be destroyed with PetscViewerDestroy().
6165c6c1daeSBarry Smith 
617811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
618db781477SPatrick Sanan           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
619db781477SPatrick Sanan           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
6205c6c1daeSBarry Smith @*/
621*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
622*d71ae5a4SJacob Faibussowitsch {
6235c6c1daeSBarry Smith   PetscFunctionBegin;
6249566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, hdf5v));
6259566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
6269566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
6279566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*hdf5v, name));
6289566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*hdf5v));
6295c6c1daeSBarry Smith   PetscFunctionReturn(0);
6305c6c1daeSBarry Smith }
6315c6c1daeSBarry Smith 
6325c6c1daeSBarry Smith /*@C
6335c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
6345c6c1daeSBarry Smith 
6355c6c1daeSBarry Smith   Not collective
6365c6c1daeSBarry Smith 
6375c6c1daeSBarry Smith   Input Parameter:
638811af0c4SBarry Smith . viewer - the `PetscViewer`
6395c6c1daeSBarry Smith 
6405c6c1daeSBarry Smith   Output Parameter:
6415c6c1daeSBarry Smith . file_id - The file id
6425c6c1daeSBarry Smith 
6435c6c1daeSBarry Smith   Level: intermediate
6445c6c1daeSBarry Smith 
645811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
6465c6c1daeSBarry Smith @*/
647*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
648*d71ae5a4SJacob Faibussowitsch {
6495c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6505c6c1daeSBarry Smith 
6515c6c1daeSBarry Smith   PetscFunctionBegin;
6525c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
6535c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
6545c6c1daeSBarry Smith   PetscFunctionReturn(0);
6555c6c1daeSBarry Smith }
6565c6c1daeSBarry Smith 
6575c6c1daeSBarry Smith /*@C
6585c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
6595c6c1daeSBarry Smith 
6605c6c1daeSBarry Smith   Not collective
6615c6c1daeSBarry Smith 
6625c6c1daeSBarry Smith   Input Parameters:
663811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
6645c6c1daeSBarry Smith - name - The group name
6655c6c1daeSBarry Smith 
6665c6c1daeSBarry Smith   Level: intermediate
6675c6c1daeSBarry Smith 
6682e361344SVaclav Hapla   Notes:
6692e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
670811af0c4SBarry Smith .vb
671811af0c4SBarry 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.
672811af0c4SBarry Smith   NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
673811af0c4SBarry Smith   "." means the current group is pushed again.
674811af0c4SBarry Smith .ve
6752e361344SVaclav Hapla 
6762e361344SVaclav Hapla   Example:
6772e361344SVaclav Hapla   Suppose the current group is "/a".
6782e361344SVaclav Hapla   + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
6792e361344SVaclav Hapla   . If name is ".", then the new group will be "/a".
6802e361344SVaclav Hapla   . If name is "b", then the new group will be "/a/b".
6812e361344SVaclav Hapla   - If name is "/b", then the new group will be "/b".
6822e361344SVaclav Hapla 
683811af0c4SBarry Smith   Developer Note:
6842e361344SVaclav Hapla   The root group "/" is internally stored as NULL.
685820fc2d1SVaclav Hapla 
68663cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
6875c6c1daeSBarry Smith @*/
688*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
689*d71ae5a4SJacob Faibussowitsch {
6905c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
6917d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
6922e361344SVaclav Hapla   size_t                    i, len;
6932e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
6942e361344SVaclav Hapla   const char               *gname;
6955c6c1daeSBarry Smith 
6965c6c1daeSBarry Smith   PetscFunctionBegin;
6975c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
698820fc2d1SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
6999566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
7002e361344SVaclav Hapla   gname = NULL;
7012e361344SVaclav Hapla   if (len) {
7022e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
7032e361344SVaclav Hapla       /* use current name */
7042e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
7052e361344SVaclav Hapla     } else if (name[0] == '/') {
7062e361344SVaclav Hapla       /* absolute */
7072e361344SVaclav Hapla       for (i = 1; i < len; i++) {
7082e361344SVaclav Hapla         if (name[i] != '/') {
7092e361344SVaclav Hapla           gname = name;
7102e361344SVaclav Hapla           break;
7112e361344SVaclav Hapla         }
7122e361344SVaclav Hapla       }
7132e361344SVaclav Hapla     } else {
7142e361344SVaclav Hapla       /* relative */
7152e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
7169566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
7172e361344SVaclav Hapla       gname = buf;
7182e361344SVaclav Hapla     }
7192e361344SVaclav Hapla   }
7209566063dSJacob Faibussowitsch   PetscCall(PetscNew(&groupNode));
7219566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
7225c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
7235c6c1daeSBarry Smith   hdf5->groups    = groupNode;
7245c6c1daeSBarry Smith   PetscFunctionReturn(0);
7255c6c1daeSBarry Smith }
7265c6c1daeSBarry Smith 
7273ef9c667SSatish Balay /*@
7285c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith   Not collective
7315c6c1daeSBarry Smith 
7325c6c1daeSBarry Smith   Input Parameter:
733811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   Level: intermediate
7365c6c1daeSBarry Smith 
73763cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
7385c6c1daeSBarry Smith @*/
739*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
740*d71ae5a4SJacob Faibussowitsch {
7415c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7427d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
7435c6c1daeSBarry Smith 
7445c6c1daeSBarry Smith   PetscFunctionBegin;
7455c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7465f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
7475c6c1daeSBarry Smith   groupNode    = hdf5->groups;
7485c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
7499566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
7509566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
7515c6c1daeSBarry Smith   PetscFunctionReturn(0);
7525c6c1daeSBarry Smith }
7535c6c1daeSBarry Smith 
754*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
755*d71ae5a4SJacob Faibussowitsch {
7565c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7575c6c1daeSBarry Smith 
7585c6c1daeSBarry Smith   PetscFunctionBegin;
7595c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
760c959eef4SJed Brown   PetscValidPointer(name, 2);
761a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
7620298fd71SBarry Smith   else *name = NULL;
7635c6c1daeSBarry Smith   PetscFunctionReturn(0);
7645c6c1daeSBarry Smith }
7655c6c1daeSBarry Smith 
7663014b61aSVaclav Hapla /*@C
767811af0c4SBarry Smith   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
768874270d9SVaclav Hapla   and return this group's ID and file ID.
769811af0c4SBarry Smith   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
770874270d9SVaclav Hapla 
771874270d9SVaclav Hapla   Not collective
772874270d9SVaclav Hapla 
7733014b61aSVaclav Hapla   Input Parameters:
7743014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7753014b61aSVaclav Hapla - path - (Optional) The path relative to the pushed group
776874270d9SVaclav Hapla 
777d8d19677SJose E. Roman   Output Parameters:
778874270d9SVaclav Hapla + fileId - The HDF5 file ID
779874270d9SVaclav Hapla - groupId - The HDF5 group ID
780874270d9SVaclav Hapla 
781811af0c4SBarry Smith   Note:
7823014b61aSVaclav 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.
7833014b61aSVaclav Hapla   So NULL or empty path means the current pushed group.
7843014b61aSVaclav Hapla 
785e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
786e74616beSVaclav Hapla 
787874270d9SVaclav Hapla   Level: intermediate
788874270d9SVaclav Hapla 
78963cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
790874270d9SVaclav Hapla @*/
791*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
792*d71ae5a4SJacob Faibussowitsch {
79390e03537SVaclav Hapla   hid_t       file_id;
79490e03537SVaclav Hapla   H5O_type_t  type;
7953014b61aSVaclav Hapla   const char *fileName  = NULL;
7963014b61aSVaclav Hapla   char       *groupName = NULL;
79776d59af2SVaclav Hapla   PetscBool   writable, has;
79854dbf706SVaclav Hapla 
79954dbf706SVaclav Hapla   PetscFunctionBegin;
8009566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
8019566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
8029566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
8033014b61aSVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
8049566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
80576d59af2SVaclav Hapla   if (!has) {
8065f80ce2aSJacob 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);
807f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
80876d59af2SVaclav Hapla   }
8095f80ce2aSJacob 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);
8103014b61aSVaclav Hapla   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
8113014b61aSVaclav Hapla   PetscCall(PetscFree(groupName));
81254dbf706SVaclav Hapla   *fileId = file_id;
81354dbf706SVaclav Hapla   PetscFunctionReturn(0);
81454dbf706SVaclav Hapla }
81554dbf706SVaclav Hapla 
81663cb69f5SVaclav Hapla /*@C
81763cb69f5SVaclav Hapla   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file
81863cb69f5SVaclav Hapla 
81963cb69f5SVaclav Hapla   Not collective
82063cb69f5SVaclav Hapla 
82163cb69f5SVaclav Hapla   Input Parameters:
82263cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
82363cb69f5SVaclav Hapla - path - (Optional) The path relative to the pushed group
82463cb69f5SVaclav Hapla 
82563cb69f5SVaclav Hapla   Note:
82663cb69f5SVaclav 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.
82763cb69f5SVaclav Hapla   So NULL or empty path means the current pushed group.
82863cb69f5SVaclav Hapla 
82963cb69f5SVaclav Hapla   This will fail if the viewer is not writable.
83063cb69f5SVaclav Hapla 
83163cb69f5SVaclav Hapla   Level: intermediate
83263cb69f5SVaclav Hapla 
83363cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
83463cb69f5SVaclav Hapla @*/
835*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
836*d71ae5a4SJacob Faibussowitsch {
83763cb69f5SVaclav Hapla   hid_t fileId, groupId;
83863cb69f5SVaclav Hapla 
83963cb69f5SVaclav Hapla   PetscFunctionBegin;
84063cb69f5SVaclav Hapla   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
84163cb69f5SVaclav Hapla   PetscCallHDF5(H5Gclose, (groupId));
84263cb69f5SVaclav Hapla   PetscFunctionReturn(0);
84363cb69f5SVaclav Hapla }
84463cb69f5SVaclav Hapla 
8453ef9c667SSatish Balay /*@
846d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
8475c6c1daeSBarry Smith 
8485c6c1daeSBarry Smith   Not collective
8495c6c1daeSBarry Smith 
8505c6c1daeSBarry Smith   Input Parameter:
851811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
8525c6c1daeSBarry Smith 
8535c6c1daeSBarry Smith   Level: intermediate
8545c6c1daeSBarry Smith 
855d7dd068bSVaclav Hapla   Notes:
856811af0c4SBarry Smith   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
857811af0c4SBarry Smith   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
858811af0c4SBarry 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()`].
859811af0c4SBarry Smith   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
860811af0c4SBarry Smith   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
861d7dd068bSVaclav Hapla 
862d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
863d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
864d7dd068bSVaclav Hapla 
865811af0c4SBarry Smith   Developer note:
866d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
867d7dd068bSVaclav Hapla 
868811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
869d7dd068bSVaclav Hapla @*/
870*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
871*d71ae5a4SJacob Faibussowitsch {
872d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
873d7dd068bSVaclav Hapla 
874d7dd068bSVaclav Hapla   PetscFunctionBegin;
875d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8765f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
877d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
878d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
879d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
880d7dd068bSVaclav Hapla }
881d7dd068bSVaclav Hapla 
882d7dd068bSVaclav Hapla /*@
883d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
884d7dd068bSVaclav Hapla 
885d7dd068bSVaclav Hapla   Not collective
886d7dd068bSVaclav Hapla 
887d7dd068bSVaclav Hapla   Input Parameter:
888811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
889d7dd068bSVaclav Hapla 
890d7dd068bSVaclav Hapla   Level: intermediate
891d7dd068bSVaclav Hapla 
892811af0c4SBarry Smith   Note:
893811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
894d7dd068bSVaclav Hapla 
895811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
896d7dd068bSVaclav Hapla @*/
897*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
898*d71ae5a4SJacob Faibussowitsch {
899d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
900d7dd068bSVaclav Hapla 
901d7dd068bSVaclav Hapla   PetscFunctionBegin;
902d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9035f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
904d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
905d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
906d7dd068bSVaclav Hapla }
907d7dd068bSVaclav Hapla 
908d7dd068bSVaclav Hapla /*@
909d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
910d7dd068bSVaclav Hapla 
911d7dd068bSVaclav Hapla   Not collective
912d7dd068bSVaclav Hapla 
913d7dd068bSVaclav Hapla   Input Parameter:
914811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
915d7dd068bSVaclav Hapla 
916d7dd068bSVaclav Hapla   Output Parameter:
917d7dd068bSVaclav Hapla . flg - is timestepping active?
918d7dd068bSVaclav Hapla 
919d7dd068bSVaclav Hapla   Level: intermediate
920d7dd068bSVaclav Hapla 
921811af0c4SBarry Smith   Note:
922811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
923d7dd068bSVaclav Hapla 
924811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
925d7dd068bSVaclav Hapla @*/
926*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
927*d71ae5a4SJacob Faibussowitsch {
928d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
929d7dd068bSVaclav Hapla 
930d7dd068bSVaclav Hapla   PetscFunctionBegin;
931d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
932d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
933d7dd068bSVaclav Hapla   PetscFunctionReturn(0);
934d7dd068bSVaclav Hapla }
935d7dd068bSVaclav Hapla 
936d7dd068bSVaclav Hapla /*@
937d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
938d7dd068bSVaclav Hapla 
939d7dd068bSVaclav Hapla   Not collective
940d7dd068bSVaclav Hapla 
941d7dd068bSVaclav Hapla   Input Parameter:
942811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
943d7dd068bSVaclav Hapla 
944d7dd068bSVaclav Hapla   Level: intermediate
945d7dd068bSVaclav Hapla 
946811af0c4SBarry Smith   Note:
947811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
948d7dd068bSVaclav Hapla 
949811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
9505c6c1daeSBarry Smith @*/
951*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
952*d71ae5a4SJacob Faibussowitsch {
9535c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9545c6c1daeSBarry Smith 
9555c6c1daeSBarry Smith   PetscFunctionBegin;
9565c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9575f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9585c6c1daeSBarry Smith   ++hdf5->timestep;
9595c6c1daeSBarry Smith   PetscFunctionReturn(0);
9605c6c1daeSBarry Smith }
9615c6c1daeSBarry Smith 
9623ef9c667SSatish Balay /*@
963d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
9645c6c1daeSBarry Smith 
965d7dd068bSVaclav Hapla   Logically collective
9665c6c1daeSBarry Smith 
9675c6c1daeSBarry Smith   Input Parameters:
968811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
969d7dd068bSVaclav Hapla - timestep - The timestep
9705c6c1daeSBarry Smith 
9715c6c1daeSBarry Smith   Level: intermediate
9725c6c1daeSBarry Smith 
973811af0c4SBarry Smith   Note:
974811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
975d7dd068bSVaclav Hapla 
976811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
9775c6c1daeSBarry Smith @*/
978*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
979*d71ae5a4SJacob Faibussowitsch {
9805c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
9815c6c1daeSBarry Smith 
9825c6c1daeSBarry Smith   PetscFunctionBegin;
9835c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
984d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
9855f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
9865f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
9875c6c1daeSBarry Smith   hdf5->timestep = timestep;
9885c6c1daeSBarry Smith   PetscFunctionReturn(0);
9895c6c1daeSBarry Smith }
9905c6c1daeSBarry Smith 
9913ef9c667SSatish Balay /*@
9925c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
9935c6c1daeSBarry Smith 
9945c6c1daeSBarry Smith   Not collective
9955c6c1daeSBarry Smith 
9965c6c1daeSBarry Smith   Input Parameter:
997811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
9985c6c1daeSBarry Smith 
9995c6c1daeSBarry Smith   Output Parameter:
1000d7dd068bSVaclav Hapla . timestep - The timestep
10015c6c1daeSBarry Smith 
10025c6c1daeSBarry Smith   Level: intermediate
10035c6c1daeSBarry Smith 
1004811af0c4SBarry Smith   Note:
1005811af0c4SBarry Smith   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1006d7dd068bSVaclav Hapla 
1007811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
10085c6c1daeSBarry Smith @*/
1009*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
1010*d71ae5a4SJacob Faibussowitsch {
10115c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
10125c6c1daeSBarry Smith 
10135c6c1daeSBarry Smith   PetscFunctionBegin;
10145c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1015d7dd068bSVaclav Hapla   PetscValidIntPointer(timestep, 2);
10165f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
10175c6c1daeSBarry Smith   *timestep = hdf5->timestep;
10185c6c1daeSBarry Smith   PetscFunctionReturn(0);
10195c6c1daeSBarry Smith }
10205c6c1daeSBarry Smith 
102136bce6e8SMatthew G. Knepley /*@C
102236bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
102336bce6e8SMatthew G. Knepley 
102436bce6e8SMatthew G. Knepley   Not collective
102536bce6e8SMatthew G. Knepley 
102636bce6e8SMatthew G. Knepley   Input Parameter:
1027811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
102836bce6e8SMatthew G. Knepley 
102936bce6e8SMatthew G. Knepley   Output Parameter:
1030811af0c4SBarry Smith . mtype - the HDF5  datatype
103136bce6e8SMatthew G. Knepley 
103236bce6e8SMatthew G. Knepley   Level: advanced
103336bce6e8SMatthew G. Knepley 
1034db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
103536bce6e8SMatthew G. Knepley @*/
1036*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
1037*d71ae5a4SJacob Faibussowitsch {
103836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
103936bce6e8SMatthew G. Knepley   if (ptype == PETSC_INT)
104036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
104136bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_LLONG;
104236bce6e8SMatthew G. Knepley #else
104336bce6e8SMatthew G. Knepley     *htype = H5T_NATIVE_INT;
104436bce6e8SMatthew G. Knepley #endif
104536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
104636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
104736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1048c9450970SBarry Smith   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
1049de3ba262SVaclav Hapla   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT;
105036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
105136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
105236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
10537e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
105436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
105536bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
105636bce6e8SMatthew G. Knepley }
105736bce6e8SMatthew G. Knepley 
105836bce6e8SMatthew G. Knepley /*@C
105936bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
106036bce6e8SMatthew G. Knepley 
106136bce6e8SMatthew G. Knepley   Not collective
106236bce6e8SMatthew G. Knepley 
106336bce6e8SMatthew G. Knepley   Input Parameter:
1064811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
106536bce6e8SMatthew G. Knepley 
106636bce6e8SMatthew G. Knepley   Output Parameter:
1067811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
106836bce6e8SMatthew G. Knepley 
106936bce6e8SMatthew G. Knepley   Level: advanced
107036bce6e8SMatthew G. Knepley 
1071db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
107236bce6e8SMatthew G. Knepley @*/
1073*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
1074*d71ae5a4SJacob Faibussowitsch {
107536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
107636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
107736bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG;
107836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
107936bce6e8SMatthew G. Knepley #else
108036bce6e8SMatthew G. Knepley   if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT;
108136bce6e8SMatthew G. Knepley #endif
108236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
108336bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
108436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
108536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
108636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
108736bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
10887e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
108936bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
109036bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
109136bce6e8SMatthew G. Knepley }
109236bce6e8SMatthew G. Knepley 
1093df863907SAlex Fikl /*@C
1094b8ee85a1SVaclav Hapla  PetscViewerHDF5WriteAttribute - Write an attribute
109536bce6e8SMatthew G. Knepley 
1096343a20b2SVaclav Hapla   Collective
1097343a20b2SVaclav Hapla 
109836bce6e8SMatthew G. Knepley   Input Parameters:
1099811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
11004302210dSVaclav Hapla . parent - The parent dataset/group name
110136bce6e8SMatthew G. Knepley . name   - The attribute name
110236bce6e8SMatthew G. Knepley . datatype - The attribute type
110336bce6e8SMatthew G. Knepley - value    - The attribute value
110436bce6e8SMatthew G. Knepley 
110536bce6e8SMatthew G. Knepley   Level: advanced
110636bce6e8SMatthew G. Knepley 
1107811af0c4SBarry Smith   Note:
1108343a20b2SVaclav 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.
11094302210dSVaclav Hapla 
1110811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1111811af0c4SBarry Smith           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
111236bce6e8SMatthew G. Knepley @*/
1113*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1114*d71ae5a4SJacob Faibussowitsch {
11154302210dSVaclav Hapla   char     *parentAbsPath;
111660568592SMatthew G. Knepley   hid_t     h5, dataspace, obj, attribute, dtype;
1117080f144cSVaclav Hapla   PetscBool has;
111836bce6e8SMatthew G. Knepley 
111936bce6e8SMatthew G. Knepley   PetscFunctionBegin;
11205cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
11214302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1122c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
11234302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
1124b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
112577717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
11269566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
11279566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
11289566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11297e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
11307e97c476SMatthew G. Knepley     size_t len;
11319566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *)value, &len));
1132792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
11337e97c476SMatthew G. Knepley   }
11349566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1135792fecdfSBarry Smith   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1136792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1137080f144cSVaclav Hapla   if (has) {
1138792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1139080f144cSVaclav Hapla   } else {
1140792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1141080f144cSVaclav Hapla   }
1142792fecdfSBarry Smith   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1143792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1144792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1145792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
1146792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (dataspace));
11479566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
114836bce6e8SMatthew G. Knepley   PetscFunctionReturn(0);
114936bce6e8SMatthew G. Knepley }
115036bce6e8SMatthew G. Knepley 
1151df863907SAlex Fikl /*@C
1152811af0c4SBarry Smith  PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1153577e0e04SVaclav Hapla 
1154343a20b2SVaclav Hapla   Collective
1155343a20b2SVaclav Hapla 
1156577e0e04SVaclav Hapla   Input Parameters:
1157811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1158577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1159577e0e04SVaclav Hapla . name     - The attribute name
1160577e0e04SVaclav Hapla . datatype - The attribute type
1161577e0e04SVaclav Hapla - value    - The attribute value
1162577e0e04SVaclav Hapla 
1163811af0c4SBarry Smith   Note:
1164343a20b2SVaclav 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).
1165811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1166577e0e04SVaclav Hapla 
1167577e0e04SVaclav Hapla   Level: advanced
1168577e0e04SVaclav Hapla 
1169811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1170811af0c4SBarry Smith           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1171577e0e04SVaclav Hapla @*/
1172*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1173*d71ae5a4SJacob Faibussowitsch {
1174577e0e04SVaclav Hapla   PetscFunctionBegin;
1175577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1176577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1177577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1178b7e81f0fSVaclav Hapla   PetscValidPointer(value, 5);
11799566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
11809566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
1181577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1182577e0e04SVaclav Hapla }
1183577e0e04SVaclav Hapla 
1184577e0e04SVaclav Hapla /*@C
1185b8ee85a1SVaclav Hapla  PetscViewerHDF5ReadAttribute - Read an attribute
1186ad0c61feSMatthew G. Knepley 
1187343a20b2SVaclav Hapla   Collective
1188343a20b2SVaclav Hapla 
1189ad0c61feSMatthew G. Knepley   Input Parameters:
1190811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
11914302210dSVaclav Hapla . parent - The parent dataset/group name
1192ad0c61feSMatthew G. Knepley . name   - The attribute name
1193a2d6be1bSVaclav Hapla . datatype - The attribute type
1194a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1195ad0c61feSMatthew G. Knepley 
1196ad0c61feSMatthew G. Knepley   Output Parameter:
1197a2d6be1bSVaclav Hapla . value    - The pointer to the read HDF5 attribute value
1198ad0c61feSMatthew G. Knepley 
1199a2d6be1bSVaclav Hapla   Notes:
1200a2d6be1bSVaclav Hapla   If defaultValue is NULL and the attribute is not found, an error occurs.
1201a2d6be1bSVaclav Hapla   If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value.
1202a2d6be1bSVaclav Hapla   The pointers defaultValue and value can be the same; for instance
1203811af0c4SBarry Smith $  flg = `PETSC_FALSE`;
1204811af0c4SBarry Smith $  PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
1205a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1206a2d6be1bSVaclav Hapla 
1207811af0c4SBarry Smith   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1208708d4cb3SBarry Smith 
1209343a20b2SVaclav 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.
1210ad0c61feSMatthew G. Knepley 
1211343a20b2SVaclav Hapla   Level: advanced
12124302210dSVaclav Hapla 
1213811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1214ad0c61feSMatthew G. Knepley @*/
1215*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1216*d71ae5a4SJacob Faibussowitsch {
12174302210dSVaclav Hapla   char     *parentAbsPath;
1218a2d6be1bSVaclav Hapla   hid_t     h5, obj, attribute, dtype;
1219969748fdSVaclav Hapla   PetscBool has;
1220ad0c61feSMatthew G. Knepley 
1221ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
12225cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12234302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1224c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
1225a2d6be1bSVaclav Hapla   if (defaultValue) PetscValidPointer(defaultValue, 5);
1226a2d6be1bSVaclav Hapla   PetscValidPointer(value, 6);
12279566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
122877717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
12299566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
12309566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1231a2d6be1bSVaclav Hapla   if (!has) {
1232a2d6be1bSVaclav Hapla     if (defaultValue) {
1233a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1234a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
12359566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1236a2d6be1bSVaclav Hapla         } else {
1237a2d6be1bSVaclav Hapla           size_t len;
1238792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
12399566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1240a2d6be1bSVaclav Hapla         }
1241a2d6be1bSVaclav Hapla       }
12429566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
1243a2d6be1bSVaclav Hapla       PetscFunctionReturn(0);
124498921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1245a2d6be1bSVaclav Hapla   }
12469566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1247792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1248792fecdfSBarry Smith   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1249f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1250f0b58479SMatthew G. Knepley     size_t len;
1251a2d6be1bSVaclav Hapla     hid_t  atype;
1252792fecdfSBarry Smith     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1253792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
12549566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1255792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1256792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1257708d4cb3SBarry Smith   } else {
1258792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1259708d4cb3SBarry Smith   }
1260792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1261e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1262792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
12639566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
1264ad0c61feSMatthew G. Knepley   PetscFunctionReturn(0);
1265ad0c61feSMatthew G. Knepley }
1266ad0c61feSMatthew G. Knepley 
1267577e0e04SVaclav Hapla /*@C
1268811af0c4SBarry Smith  PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1269577e0e04SVaclav Hapla 
1270343a20b2SVaclav Hapla   Collective
1271343a20b2SVaclav Hapla 
1272577e0e04SVaclav Hapla   Input Parameters:
1273811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1274577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1275577e0e04SVaclav Hapla . name     - The attribute name
1276577e0e04SVaclav Hapla - datatype - The attribute type
1277577e0e04SVaclav Hapla 
1278577e0e04SVaclav Hapla   Output Parameter:
1279577e0e04SVaclav Hapla . value    - The attribute value
1280577e0e04SVaclav Hapla 
1281811af0c4SBarry Smith   Note:
1282577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1283811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1284577e0e04SVaclav Hapla 
1285577e0e04SVaclav Hapla   Level: advanced
1286577e0e04SVaclav Hapla 
1287811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1288577e0e04SVaclav Hapla @*/
1289*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1290*d71ae5a4SJacob Faibussowitsch {
1291577e0e04SVaclav Hapla   PetscFunctionBegin;
1292577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1293577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1294577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
1295064a246eSJacob Faibussowitsch   PetscValidPointer(value, 6);
12969566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12979566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
1298577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1299577e0e04SVaclav Hapla }
1300577e0e04SVaclav Hapla 
1301*d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1302*d71ae5a4SJacob Faibussowitsch {
1303a75c3fd4SVaclav Hapla   htri_t exists;
1304a75c3fd4SVaclav Hapla   hid_t  group;
1305a75c3fd4SVaclav Hapla 
1306a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1307792fecdfSBarry Smith   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1308792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1309a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1310792fecdfSBarry Smith     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1311792fecdfSBarry Smith     PetscCallHDF5(H5Gclose, (group));
1312a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1313a75c3fd4SVaclav Hapla   }
1314a75c3fd4SVaclav Hapla   *exists_ = (PetscBool)exists;
1315a75c3fd4SVaclav Hapla   PetscFunctionReturn(0);
1316a75c3fd4SVaclav Hapla }
1317a75c3fd4SVaclav Hapla 
1318*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1319*d71ae5a4SJacob Faibussowitsch {
132090e03537SVaclav Hapla   const char rootGroupName[] = "/";
13215cdeb108SMatthew G. Knepley   hid_t      h5;
1322e5bf9ebcSVaclav Hapla   PetscBool  exists = PETSC_FALSE;
132315b861d2SVaclav Hapla   PetscInt   i;
132415b861d2SVaclav Hapla   int        n;
132585688ae2SVaclav Hapla   char     **hierarchy;
132685688ae2SVaclav Hapla   char       buf[PETSC_MAX_PATH_LEN] = "";
13275cdeb108SMatthew G. Knepley 
13285cdeb108SMatthew G. Knepley   PetscFunctionBegin;
13295cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
133090e03537SVaclav Hapla   if (name) PetscValidCharPointer(name, 2);
133190e03537SVaclav Hapla   else name = rootGroupName;
1332ccd66a83SVaclav Hapla   if (has) {
1333064a246eSJacob Faibussowitsch     PetscValidBoolPointer(has, 4);
13345cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1335ccd66a83SVaclav Hapla   }
1336ccd66a83SVaclav Hapla   if (otype) {
1337064a246eSJacob Faibussowitsch     PetscValidIntPointer(otype, 5);
133856cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1339ccd66a83SVaclav Hapla   }
13409566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
134185688ae2SVaclav Hapla 
134285688ae2SVaclav Hapla   /*
134385688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
134485688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
134585688ae2SVaclav Hapla      1) whether it's a valid link
134685688ae2SVaclav Hapla      2) whether this link resolves to an object
134785688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
134885688ae2SVaclav Hapla   */
13499566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
135085688ae2SVaclav Hapla   if (!n) {
135185688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1352ccd66a83SVaclav Hapla     if (has) *has = PETSC_TRUE;
1353ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
13549566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
135585688ae2SVaclav Hapla     PetscFunctionReturn(0);
135685688ae2SVaclav Hapla   }
135785688ae2SVaclav Hapla   for (i = 0; i < n; i++) {
13589566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, "/"));
13599566063dSJacob Faibussowitsch     PetscCall(PetscStrcat(buf, hierarchy[i]));
13609566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1361a75c3fd4SVaclav Hapla     if (!exists) break;
136285688ae2SVaclav Hapla   }
13639566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
136485688ae2SVaclav Hapla 
136585688ae2SVaclav Hapla   /* If the object exists, get its type */
1366ccd66a83SVaclav Hapla   if (exists && otype) {
13675cdeb108SMatthew G. Knepley     H5O_info_t info;
13685cdeb108SMatthew G. Knepley 
136976276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1370792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
137156cc0592SVaclav Hapla     *otype = info.type;
13725cdeb108SMatthew G. Knepley   }
1373ccd66a83SVaclav Hapla   if (has) *has = exists;
13745cdeb108SMatthew G. Knepley   PetscFunctionReturn(0);
13755cdeb108SMatthew G. Knepley }
13765cdeb108SMatthew G. Knepley 
13774302210dSVaclav Hapla /*@C
13780a9f38efSVaclav Hapla  PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
13790a9f38efSVaclav Hapla 
1380343a20b2SVaclav Hapla   Collective
1381343a20b2SVaclav Hapla 
13820a9f38efSVaclav Hapla   Input Parameters:
1383811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1384a0558868SVaclav Hapla - path - (Optional) The path relative to the pushed group
13850a9f38efSVaclav Hapla 
13860a9f38efSVaclav Hapla   Output Parameter:
13870a9f38efSVaclav Hapla . has - Flag for group existence
13880a9f38efSVaclav Hapla 
13890a9f38efSVaclav Hapla   Level: advanced
13900a9f38efSVaclav Hapla 
13914302210dSVaclav Hapla   Notes:
1392785c443eSVaclav 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.
1393785c443eSVaclav Hapla   So NULL or empty path means the current pushed group.
13944302210dSVaclav Hapla 
1395811af0c4SBarry Smith   If path exists but is not a group, `PETSC_FALSE` is returned.
1396811af0c4SBarry Smith 
1397811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
13980a9f38efSVaclav Hapla @*/
1399*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1400*d71ae5a4SJacob Faibussowitsch {
14010a9f38efSVaclav Hapla   H5O_type_t type;
14024302210dSVaclav Hapla   char      *abspath;
14030a9f38efSVaclav Hapla 
14040a9f38efSVaclav Hapla   PetscFunctionBegin;
14050a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14064302210dSVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
14074302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
140877717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
14099566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
14104302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
14119566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
14120a9f38efSVaclav Hapla   PetscFunctionReturn(0);
14130a9f38efSVaclav Hapla }
14140a9f38efSVaclav Hapla 
141589e0ef10SVaclav Hapla /*@C
141689e0ef10SVaclav Hapla  PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
141789e0ef10SVaclav Hapla 
1418343a20b2SVaclav Hapla   Collective
1419343a20b2SVaclav Hapla 
142089e0ef10SVaclav Hapla   Input Parameters:
1421811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
142289e0ef10SVaclav Hapla - path - The dataset path
142389e0ef10SVaclav Hapla 
142489e0ef10SVaclav Hapla   Output Parameter:
142589e0ef10SVaclav Hapla . has - Flag whether dataset exists
142689e0ef10SVaclav Hapla 
142789e0ef10SVaclav Hapla   Level: advanced
142889e0ef10SVaclav Hapla 
142989e0ef10SVaclav Hapla   Notes:
1430343a20b2SVaclav 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.
143189e0ef10SVaclav Hapla 
1432811af0c4SBarry Smith   If path is NULL or empty, has is set to `PETSC_FALSE`.
1433811af0c4SBarry Smith 
1434811af0c4SBarry Smith   If path exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1435811af0c4SBarry Smith 
1436811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
143789e0ef10SVaclav Hapla @*/
1438*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1439*d71ae5a4SJacob Faibussowitsch {
144089e0ef10SVaclav Hapla   H5O_type_t type;
144189e0ef10SVaclav Hapla   char      *abspath;
144289e0ef10SVaclav Hapla 
144389e0ef10SVaclav Hapla   PetscFunctionBegin;
144489e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
144589e0ef10SVaclav Hapla   if (path) PetscValidCharPointer(path, 2);
144689e0ef10SVaclav Hapla   PetscValidBoolPointer(has, 3);
144777717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
14489566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
144989e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
14509566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
145189e0ef10SVaclav Hapla   PetscFunctionReturn(0);
145289e0ef10SVaclav Hapla }
145389e0ef10SVaclav Hapla 
14540a9f38efSVaclav Hapla /*@
1455e3f143f7SVaclav Hapla  PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1456ecce1506SVaclav Hapla 
1457343a20b2SVaclav Hapla   Collective
1458343a20b2SVaclav Hapla 
1459ecce1506SVaclav Hapla   Input Parameters:
1460811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1461ecce1506SVaclav Hapla - obj    - The named object
1462ecce1506SVaclav Hapla 
1463ecce1506SVaclav Hapla   Output Parameter:
146489e0ef10SVaclav Hapla . has    - Flag for dataset existence
1465ecce1506SVaclav Hapla 
1466e3f143f7SVaclav Hapla   Notes:
146789e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1468811af0c4SBarry Smith 
1469811af0c4SBarry Smith   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1470e3f143f7SVaclav Hapla 
1471ecce1506SVaclav Hapla   Level: advanced
1472ecce1506SVaclav Hapla 
1473811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1474ecce1506SVaclav Hapla @*/
1475*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1476*d71ae5a4SJacob Faibussowitsch {
147789e0ef10SVaclav Hapla   size_t len;
1478ecce1506SVaclav Hapla 
1479ecce1506SVaclav Hapla   PetscFunctionBegin;
1480c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1481c57a1a9aSVaclav Hapla   PetscValidHeader(obj, 2);
14824302210dSVaclav Hapla   PetscValidBoolPointer(has, 3);
14839566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
14845f80ce2aSJacob Faibussowitsch   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
14859566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
1486ecce1506SVaclav Hapla   PetscFunctionReturn(0);
1487ecce1506SVaclav Hapla }
1488ecce1506SVaclav Hapla 
1489df863907SAlex Fikl /*@C
1490ece83b6cSVaclav Hapla  PetscViewerHDF5HasAttribute - Check whether an attribute exists
1491e7bdbf8eSMatthew G. Knepley 
1492343a20b2SVaclav Hapla   Collective
1493343a20b2SVaclav Hapla 
1494e7bdbf8eSMatthew G. Knepley   Input Parameters:
1495811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
14964302210dSVaclav Hapla . parent - The parent dataset/group name
1497e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1498e7bdbf8eSMatthew G. Knepley 
1499e7bdbf8eSMatthew G. Knepley   Output Parameter:
1500e7bdbf8eSMatthew G. Knepley . has    - Flag for attribute existence
1501e7bdbf8eSMatthew G. Knepley 
1502e7bdbf8eSMatthew G. Knepley   Level: advanced
1503e7bdbf8eSMatthew G. Knepley 
1504811af0c4SBarry Smith   Note:
1505343a20b2SVaclav 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.
15064302210dSVaclav Hapla 
1507811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1508e7bdbf8eSMatthew G. Knepley @*/
1509*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1510*d71ae5a4SJacob Faibussowitsch {
15114302210dSVaclav Hapla   char *parentAbsPath;
1512e7bdbf8eSMatthew G. Knepley 
1513e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
15145cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
15154302210dSVaclav Hapla   if (parent) PetscValidCharPointer(parent, 2);
1516c57a1a9aSVaclav Hapla   PetscValidCharPointer(name, 3);
15174302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
151877717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
15199566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
15209566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
15219566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
152206db490cSVaclav Hapla   PetscFunctionReturn(0);
152306db490cSVaclav Hapla }
152406db490cSVaclav Hapla 
1525577e0e04SVaclav Hapla /*@C
1526811af0c4SBarry Smith  PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1527577e0e04SVaclav Hapla 
1528343a20b2SVaclav Hapla   Collective
1529343a20b2SVaclav Hapla 
1530577e0e04SVaclav Hapla   Input Parameters:
1531811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1532577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1533577e0e04SVaclav Hapla - name   - The attribute name
1534577e0e04SVaclav Hapla 
1535577e0e04SVaclav Hapla   Output Parameter:
1536577e0e04SVaclav Hapla . has    - Flag for attribute existence
1537577e0e04SVaclav Hapla 
1538811af0c4SBarry Smith   Note:
1539577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1540811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1541577e0e04SVaclav Hapla 
1542577e0e04SVaclav Hapla   Level: advanced
1543577e0e04SVaclav Hapla 
1544811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1545577e0e04SVaclav Hapla @*/
1546*d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1547*d71ae5a4SJacob Faibussowitsch {
1548577e0e04SVaclav Hapla   PetscFunctionBegin;
1549577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1550577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
1551577e0e04SVaclav Hapla   PetscValidCharPointer(name, 3);
15524302210dSVaclav Hapla   PetscValidBoolPointer(has, 4);
15539566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
15549566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
1555577e0e04SVaclav Hapla   PetscFunctionReturn(0);
1556577e0e04SVaclav Hapla }
1557577e0e04SVaclav Hapla 
1558*d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1559*d71ae5a4SJacob Faibussowitsch {
156006db490cSVaclav Hapla   hid_t  h5;
156106db490cSVaclav Hapla   htri_t hhas;
156206db490cSVaclav Hapla 
156306db490cSVaclav Hapla   PetscFunctionBegin;
15649566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1565792fecdfSBarry Smith   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
15665cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
1567e7bdbf8eSMatthew G. Knepley   PetscFunctionReturn(0);
1568e7bdbf8eSMatthew G. Knepley }
1569e7bdbf8eSMatthew G. Knepley 
1570a75e6a4aSMatthew G. Knepley /*
1571a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1572a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1573a75e6a4aSMatthew G. Knepley */
1574d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1575a75e6a4aSMatthew G. Knepley 
1576a75e6a4aSMatthew G. Knepley /*@C
1577811af0c4SBarry Smith   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1578a75e6a4aSMatthew G. Knepley 
1579d083f849SBarry Smith   Collective
1580a75e6a4aSMatthew G. Knepley 
1581a75e6a4aSMatthew G. Knepley   Input Parameter:
1582811af0c4SBarry Smith . comm - the MPI communicator to share the HDF5 `PetscViewer`
1583a75e6a4aSMatthew G. Knepley 
1584a75e6a4aSMatthew G. Knepley   Level: intermediate
1585a75e6a4aSMatthew G. Knepley 
1586811af0c4SBarry Smith   Options Database Key:
158710699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1588a75e6a4aSMatthew G. Knepley 
1589811af0c4SBarry Smith   Environmental variable:
1590811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1591a75e6a4aSMatthew G. Knepley 
1592811af0c4SBarry Smith   Note:
1593811af0c4SBarry Smith   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1594811af0c4SBarry Smith   an error code.  The HDF5 `PetscViewer` is usually used in the form
1595a75e6a4aSMatthew G. Knepley $       XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1596a75e6a4aSMatthew G. Knepley 
1597811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1598a75e6a4aSMatthew G. Knepley @*/
1599*d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1600*d71ae5a4SJacob Faibussowitsch {
1601a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
1602a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1603a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1604a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1605a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1606a75e6a4aSMatthew G. Knepley 
1607a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
16089371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16099371c9d4SSatish Balay   if (ierr) {
16109371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16119371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16129371c9d4SSatish Balay   }
1613a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
1614908793a3SLisandro Dalcin     ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
16159371c9d4SSatish Balay     if (ierr) {
16169371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16179371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16189371c9d4SSatish Balay     }
1619a75e6a4aSMatthew G. Knepley   }
162047435625SJed Brown   ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg);
16219371c9d4SSatish Balay   if (ierr) {
16229371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16239371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16249371c9d4SSatish Balay   }
1625a75e6a4aSMatthew G. Knepley   if (!flg) { /* PetscViewer not yet created */
1626a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
16279371c9d4SSatish Balay     if (ierr) {
16289371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16299371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16309371c9d4SSatish Balay     }
1631a75e6a4aSMatthew G. Knepley     if (!flg) {
1632a75e6a4aSMatthew G. Knepley       ierr = PetscStrcpy(fname, "output.h5");
16339371c9d4SSatish Balay       if (ierr) {
16349371c9d4SSatish Balay         PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16359371c9d4SSatish Balay         PetscFunctionReturn(NULL);
16369371c9d4SSatish Balay       }
1637a75e6a4aSMatthew G. Knepley     }
1638a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
16399371c9d4SSatish Balay     if (ierr) {
16409371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16419371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16429371c9d4SSatish Balay     }
1643a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
16449371c9d4SSatish Balay     if (ierr) {
16459371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16469371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16479371c9d4SSatish Balay     }
164847435625SJed Brown     ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
16499371c9d4SSatish Balay     if (ierr) {
16509371c9d4SSatish Balay       PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16519371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16529371c9d4SSatish Balay     }
1653a75e6a4aSMatthew G. Knepley   }
1654a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
16559371c9d4SSatish Balay   if (ierr) {
16569371c9d4SSatish Balay     PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
16579371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16589371c9d4SSatish Balay   }
1659a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1660a75e6a4aSMatthew G. Knepley }
1661