xref: /petsc/src/sys/classes/viewer/impls/hdf5/hdf5v.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
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 *);
534e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer, const char *[]);
606db490cSVaclav Hapla 
777717648SVaclav Hapla /*@C
877717648SVaclav Hapla   PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`.
977717648SVaclav Hapla 
1020f4b53cSBarry Smith   Not Collective
1177717648SVaclav Hapla 
1277717648SVaclav Hapla   Input Parameters:
1377717648SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
1477717648SVaclav Hapla - path   - (Optional) The path relative to the pushed group
1577717648SVaclav Hapla 
1677717648SVaclav Hapla   Output Parameter:
1777717648SVaclav Hapla . abspath - The absolute HDF5 path (group)
1877717648SVaclav Hapla 
1977717648SVaclav Hapla   Level: intermediate
2077717648SVaclav Hapla 
2177717648SVaclav Hapla   Notes:
2277717648SVaclav 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.
233f423023SBarry Smith   `NULL` or empty path means the current pushed group.
2477717648SVaclav Hapla 
2577717648SVaclav Hapla   The output abspath is newly allocated so needs to be freed.
2677717648SVaclav Hapla 
27d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
2877717648SVaclav Hapla @*/
PetscViewerHDF5GetGroup(PetscViewer viewer,const char path[],const char * abspath[])29ce78bad3SBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], const char *abspath[])
30d71ae5a4SJacob Faibussowitsch {
3177717648SVaclav Hapla   size_t      len;
324302210dSVaclav Hapla   PetscBool   relative = PETSC_FALSE;
336c132bc1SVaclav Hapla   const char *group;
346c132bc1SVaclav Hapla   char        buf[PETSC_MAX_PATH_LEN] = "";
356c132bc1SVaclav Hapla 
366c132bc1SVaclav Hapla   PetscFunctionBegin;
3777717648SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
384f572ea9SToby Isaac   if (path) PetscAssertPointer(path, 2);
394f572ea9SToby Isaac   PetscAssertPointer(abspath, 3);
4077717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group));
4177717648SVaclav Hapla   PetscCall(PetscStrlen(path, &len));
4277717648SVaclav Hapla   relative = (PetscBool)(!len || path[0] != '/');
434302210dSVaclav Hapla   if (relative) {
44c6a7a370SJeremy L Thompson     PetscCall(PetscStrncpy(buf, group, sizeof(buf)));
45c6a7a370SJeremy L Thompson     if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
46c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(buf, path, sizeof(buf)));
47ce78bad3SBarry Smith     PetscCall(PetscStrallocpy(buf, (char **)abspath));
484302210dSVaclav Hapla   } else {
49ce78bad3SBarry Smith     PetscCall(PetscStrallocpy(path, (char **)abspath));
504302210dSVaclav Hapla   }
513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
526c132bc1SVaclav Hapla }
536c132bc1SVaclav Hapla 
PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer,PetscObject obj)54d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj)
55d71ae5a4SJacob Faibussowitsch {
56577e0e04SVaclav Hapla   PetscBool has;
57577e0e04SVaclav Hapla 
58577e0e04SVaclav Hapla   PetscFunctionBegin;
599566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has));
60577e0e04SVaclav Hapla   if (!has) {
61ce78bad3SBarry Smith     const char *group;
6277717648SVaclav Hapla     PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group));
6377717648SVaclav Hapla     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group);
64577e0e04SVaclav Hapla   }
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
66577e0e04SVaclav Hapla }
67577e0e04SVaclav Hapla 
PetscViewerSetFromOptions_HDF5(PetscViewer v,PetscOptionItems PetscOptionsObject)68ce78bad3SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems PetscOptionsObject)
69d71ae5a4SJacob Faibussowitsch {
70ee8b9732SVaclav Hapla   PetscBool         flg  = PETSC_FALSE, set;
7182ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
7282ea9b62SBarry Smith 
7382ea9b62SBarry Smith   PetscFunctionBegin;
74d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options");
759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL));
769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL));
779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set));
789566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg));
7919a20e4cSMatthew G. Knepley   flg = PETSC_FALSE;
809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set));
819566063dSJacob Faibussowitsch   if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg));
82ea87367fSMatthew G. Knepley   PetscCall(PetscOptionsBool("-viewer_hdf5_compress", "Enable compression", "PetscViewerHDF5SetCompress", hdf5->compress, &hdf5->compress, NULL));
83d0609cedSBarry Smith   PetscOptionsHeadEnd();
843ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8582ea9b62SBarry Smith }
8682ea9b62SBarry Smith 
PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer)87d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer)
88d71ae5a4SJacob Faibussowitsch {
891b793a25SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data;
9003fa8834SVaclav Hapla   PetscBool         flg;
911b793a25SVaclav Hapla 
921b793a25SVaclav Hapla   PetscFunctionBegin;
9348a46eb9SPierre Jolivet   if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename));
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2]));
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput]));
969566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetCollective(v, &flg));
979566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent"));
989566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping]));
99ea87367fSMatthew G. Knepley   PetscCall(PetscViewerASCIIPrintf(viewer, "Compression: %s\n", PetscBools[hdf5->compress]));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1011b793a25SVaclav Hapla }
1021b793a25SVaclav Hapla 
PetscViewerFileClose_HDF5(PetscViewer viewer)103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer)
104d71ae5a4SJacob Faibussowitsch {
1055c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1065c6c1daeSBarry Smith 
1075c6c1daeSBarry Smith   PetscFunctionBegin;
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5->filename));
109792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
1103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1115c6c1daeSBarry Smith }
1125c6c1daeSBarry Smith 
PetscViewerFlush_HDF5(PetscViewer viewer)113d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer)
114d71ae5a4SJacob Faibussowitsch {
1156226335aSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1166226335aSVaclav Hapla 
1176226335aSVaclav Hapla   PetscFunctionBegin;
118792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL));
1193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1206226335aSVaclav Hapla }
1216226335aSVaclav Hapla 
PetscViewerDestroy_HDF5(PetscViewer viewer)122d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer)
123d71ae5a4SJacob Faibussowitsch {
1245c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1255c6c1daeSBarry Smith 
1265c6c1daeSBarry Smith   PetscFunctionBegin;
127792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (hdf5->dxpl_id));
1289566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileClose_HDF5(viewer));
1295c6c1daeSBarry Smith   while (hdf5->groups) {
1307d964842SVaclav Hapla     PetscViewerHDF5GroupList *tmp = hdf5->groups->next;
1315c6c1daeSBarry Smith 
1329566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups->name));
1339566063dSJacob Faibussowitsch     PetscCall(PetscFree(hdf5->groups));
1345c6c1daeSBarry Smith     hdf5->groups = tmp;
1355c6c1daeSBarry Smith   }
1369566063dSJacob Faibussowitsch   PetscCall(PetscFree(hdf5));
1379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
1389566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
1399566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
1402e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
1419566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL));
1429566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL));
1432e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL));
1442e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL));
1459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL));
1469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL));
147ea87367fSMatthew G. Knepley   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCompress_C", NULL));
148ea87367fSMatthew G. Knepley   PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCompress_C", NULL));
1493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1505c6c1daeSBarry Smith }
1515c6c1daeSBarry Smith 
PetscViewerFileSetMode_HDF5(PetscViewer viewer,PetscFileMode type)152d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type)
153d71ae5a4SJacob Faibussowitsch {
1545c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1555c6c1daeSBarry Smith 
1565c6c1daeSBarry Smith   PetscFunctionBegin;
1575c6c1daeSBarry Smith   hdf5->btype = type;
1583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1595c6c1daeSBarry Smith }
1605c6c1daeSBarry Smith 
PetscViewerFileGetMode_HDF5(PetscViewer viewer,PetscFileMode * type)161d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type)
162d71ae5a4SJacob Faibussowitsch {
1630b62783dSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1640b62783dSVaclav Hapla 
1650b62783dSVaclav Hapla   PetscFunctionBegin;
1660b62783dSVaclav Hapla   *type = hdf5->btype;
1673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1680b62783dSVaclav Hapla }
1690b62783dSVaclav Hapla 
PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer,PetscBool flg)170d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg)
171d71ae5a4SJacob Faibussowitsch {
17282ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
17382ea9b62SBarry Smith 
17482ea9b62SBarry Smith   PetscFunctionBegin;
17582ea9b62SBarry Smith   hdf5->basedimension2 = flg;
1763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17782ea9b62SBarry Smith }
17882ea9b62SBarry Smith 
179df863907SAlex Fikl /*@
18082ea9b62SBarry Smith   PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
18182ea9b62SBarry Smith   dimension of 2.
18282ea9b62SBarry Smith 
183c3339decSBarry Smith   Logically Collective
18482ea9b62SBarry Smith 
18582ea9b62SBarry Smith   Input Parameters:
186811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored
187811af0c4SBarry 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
18882ea9b62SBarry Smith 
189811af0c4SBarry Smith   Options Database Key:
19082ea9b62SBarry 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
19182ea9b62SBarry Smith 
1923f423023SBarry Smith   Level: intermediate
1933f423023SBarry Smith 
194811af0c4SBarry Smith   Note:
19595452b02SPatrick 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
19682ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
19782ea9b62SBarry Smith 
198aec76313SJacob Faibussowitsch .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
19982ea9b62SBarry Smith @*/
PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg)200d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg)
201d71ae5a4SJacob Faibussowitsch {
20282ea9b62SBarry Smith   PetscFunctionBegin;
20382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
204cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg));
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
20682ea9b62SBarry Smith }
20782ea9b62SBarry Smith 
208df863907SAlex Fikl /*@
20982ea9b62SBarry Smith   PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a
21082ea9b62SBarry Smith   dimension of 2.
21182ea9b62SBarry Smith 
212c3339decSBarry Smith   Logically Collective
21382ea9b62SBarry Smith 
21482ea9b62SBarry Smith   Input Parameter:
215811af0c4SBarry Smith . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5`
21682ea9b62SBarry Smith 
21782ea9b62SBarry Smith   Output Parameter:
218811af0c4SBarry 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
21982ea9b62SBarry Smith 
2203f423023SBarry Smith   Level: intermediate
2213f423023SBarry Smith 
222811af0c4SBarry Smith   Note:
22395452b02SPatrick 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
22482ea9b62SBarry Smith   of one when the dimension is lower. Others think the option is crazy.
22582ea9b62SBarry Smith 
226d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`
22782ea9b62SBarry Smith @*/
PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool * flg)228d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg)
229d71ae5a4SJacob Faibussowitsch {
23082ea9b62SBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
23182ea9b62SBarry Smith 
23282ea9b62SBarry Smith   PetscFunctionBegin;
23382ea9b62SBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
23482ea9b62SBarry Smith   *flg = hdf5->basedimension2;
2353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
23682ea9b62SBarry Smith }
23782ea9b62SBarry Smith 
PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer,PetscBool flg)238d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg)
239d71ae5a4SJacob Faibussowitsch {
2409a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
2419a0502c6SHåkon Strandenes 
2429a0502c6SHåkon Strandenes   PetscFunctionBegin;
2439a0502c6SHåkon Strandenes   hdf5->spoutput = flg;
2443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2459a0502c6SHåkon Strandenes }
2469a0502c6SHåkon Strandenes 
247df863907SAlex Fikl /*@
2489a0502c6SHåkon Strandenes   PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is
249811af0c4SBarry Smith   compiled with double precision `PetscReal`.
2509a0502c6SHåkon Strandenes 
251c3339decSBarry Smith   Logically Collective
2529a0502c6SHåkon Strandenes 
2539a0502c6SHåkon Strandenes   Input Parameters:
254811af0c4SBarry Smith + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored
255811af0c4SBarry Smith - flg    - if `PETSC_TRUE` the data will be written to disk with single precision
2569a0502c6SHåkon Strandenes 
257811af0c4SBarry Smith   Options Database Key:
2589a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision
2599a0502c6SHåkon Strandenes 
2603f423023SBarry Smith   Level: intermediate
2613f423023SBarry Smith 
262811af0c4SBarry Smith   Note:
26395452b02SPatrick Sanan   Setting this option does not make any difference if PETSc is compiled with single precision
2649a0502c6SHåkon Strandenes   in the first place. It does not affect reading datasets (HDF5 handle this internally).
2659a0502c6SHåkon Strandenes 
266d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
267811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5GetSPOutput()`
2689a0502c6SHåkon Strandenes @*/
PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg)269d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg)
270d71ae5a4SJacob Faibussowitsch {
2719a0502c6SHåkon Strandenes   PetscFunctionBegin;
2729a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
273cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg));
2743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2759a0502c6SHåkon Strandenes }
2769a0502c6SHåkon Strandenes 
277df863907SAlex Fikl /*@
2789a0502c6SHåkon Strandenes   PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is
279811af0c4SBarry Smith   compiled with double precision `PetscReal`.
2809a0502c6SHåkon Strandenes 
281c3339decSBarry Smith   Logically Collective
2829a0502c6SHåkon Strandenes 
2839a0502c6SHåkon Strandenes   Input Parameter:
284811af0c4SBarry Smith . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5`
2859a0502c6SHåkon Strandenes 
2869a0502c6SHåkon Strandenes   Output Parameter:
287811af0c4SBarry Smith . flg - if `PETSC_TRUE` the data will be written to disk with single precision
2889a0502c6SHåkon Strandenes 
2899a0502c6SHåkon Strandenes   Level: intermediate
2909a0502c6SHåkon Strandenes 
291d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`,
292811af0c4SBarry Smith           `PetscReal`, `PetscViewerHDF5SetSPOutput()`
2939a0502c6SHåkon Strandenes @*/
PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool * flg)294d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg)
295d71ae5a4SJacob Faibussowitsch {
2969a0502c6SHåkon Strandenes   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
2979a0502c6SHåkon Strandenes 
2989a0502c6SHåkon Strandenes   PetscFunctionBegin;
2999a0502c6SHåkon Strandenes   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
3009a0502c6SHåkon Strandenes   *flg = hdf5->spoutput;
3013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3029a0502c6SHåkon Strandenes }
3039a0502c6SHåkon Strandenes 
PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer,PetscBool flg)304d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg)
305d71ae5a4SJacob Faibussowitsch {
306ee8b9732SVaclav Hapla   PetscFunctionBegin;
307ee8b9732SVaclav Hapla   /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions
308ee8b9732SVaclav Hapla      - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */
309c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL)
310ee8b9732SVaclav Hapla   {
311ee8b9732SVaclav Hapla     PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
312792fecdfSBarry Smith     PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT));
313ee8b9732SVaclav Hapla   }
314ee8b9732SVaclav Hapla #else
3159566063dSJacob 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"));
316ee8b9732SVaclav Hapla #endif
3173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
318ee8b9732SVaclav Hapla }
319ee8b9732SVaclav Hapla 
320ee8b9732SVaclav Hapla /*@
321ee8b9732SVaclav Hapla   PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes.
322ee8b9732SVaclav Hapla 
323ee8b9732SVaclav Hapla   Logically Collective; flg must contain common value
324ee8b9732SVaclav Hapla 
325ee8b9732SVaclav Hapla   Input Parameters:
326811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
327811af0c4SBarry Smith - flg    - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default)
328ee8b9732SVaclav Hapla 
329811af0c4SBarry Smith   Options Database Key:
330ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers
331ee8b9732SVaclav Hapla 
3323f423023SBarry Smith   Level: intermediate
3333f423023SBarry Smith 
334811af0c4SBarry Smith   Note:
335ee8b9732SVaclav Hapla   Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better.
33653effed7SBarry 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.
337ee8b9732SVaclav Hapla 
338aec76313SJacob Faibussowitsch   Developer Notes:
339811af0c4SBarry Smith   In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively.
340ee8b9732SVaclav Hapla   This in turn means use of MPI_File_{read,write}_all /  MPI_File_{read,write} in the MPI-IO layer, respectively.
341ee8b9732SVaclav Hapla   See HDF5 documentation and MPI-IO documentation for details.
342ee8b9732SVaclav Hapla 
343d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
344ee8b9732SVaclav Hapla @*/
PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg)345d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg)
346d71ae5a4SJacob Faibussowitsch {
347ee8b9732SVaclav Hapla   PetscFunctionBegin;
348ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
349ee8b9732SVaclav Hapla   PetscValidLogicalCollectiveBool(viewer, flg, 2);
350cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg));
3513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
352ee8b9732SVaclav Hapla }
353ee8b9732SVaclav Hapla 
PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer,PetscBool * flg)354d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg)
355d71ae5a4SJacob Faibussowitsch {
356c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
357ee8b9732SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
358ee8b9732SVaclav Hapla   H5FD_mpio_xfer_t  mode;
359c796909eSBarry Smith #endif
360ee8b9732SVaclav Hapla 
361ee8b9732SVaclav Hapla   PetscFunctionBegin;
362c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL)
363c796909eSBarry Smith   *flg = PETSC_FALSE;
364c796909eSBarry Smith #else
365792fecdfSBarry Smith   PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode));
366ee8b9732SVaclav Hapla   *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE;
367c796909eSBarry Smith #endif
3683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
369ee8b9732SVaclav Hapla }
370ee8b9732SVaclav Hapla 
371ee8b9732SVaclav Hapla /*@
372ee8b9732SVaclav Hapla   PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes.
373ee8b9732SVaclav Hapla 
374ee8b9732SVaclav Hapla   Not Collective
375ee8b9732SVaclav Hapla 
37620f4b53cSBarry Smith   Input Parameter:
377811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer`
378ee8b9732SVaclav Hapla 
37920f4b53cSBarry Smith   Output Parameter:
380ee8b9732SVaclav Hapla . flg - the flag
381ee8b9732SVaclav Hapla 
382ee8b9732SVaclav Hapla   Level: intermediate
383ee8b9732SVaclav Hapla 
384811af0c4SBarry Smith   Note:
385811af0c4SBarry 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.
386811af0c4SBarry Smith   For more details, see `PetscViewerHDF5SetCollective()`.
387ee8b9732SVaclav Hapla 
388d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()`
389ee8b9732SVaclav Hapla @*/
PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool * flg)390d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg)
391d71ae5a4SJacob Faibussowitsch {
392ee8b9732SVaclav Hapla   PetscFunctionBegin;
393ee8b9732SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
3944f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
395ee8b9732SVaclav Hapla 
396cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg));
3973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
398ee8b9732SVaclav Hapla }
399ee8b9732SVaclav Hapla 
PetscViewerFileSetName_HDF5(PetscViewer viewer,const char name[])400d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[])
401d71ae5a4SJacob Faibussowitsch {
4025c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
4035c6c1daeSBarry Smith   hid_t             plist_id;
4045c6c1daeSBarry Smith 
4055c6c1daeSBarry Smith   PetscFunctionBegin;
406792fecdfSBarry Smith   if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id));
4079566063dSJacob Faibussowitsch   if (hdf5->filename) PetscCall(PetscFree(hdf5->filename));
4089566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(name, &hdf5->filename));
4095c6c1daeSBarry Smith   /* Set up file access property list with parallel I/O access */
410792fecdfSBarry Smith   PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS));
411c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL)
412792fecdfSBarry Smith   PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL));
413c796909eSBarry Smith #endif
4145c6c1daeSBarry Smith   /* Create or open the file collectively */
4155c6c1daeSBarry Smith   switch (hdf5->btype) {
4165c6c1daeSBarry Smith   case FILE_MODE_READ:
4178a2871f6SBarry Smith     if (PetscDefined(USE_DEBUG)) {
4188a2871f6SBarry Smith       PetscMPIInt rank;
4198a2871f6SBarry Smith       PetscBool   flg;
4208a2871f6SBarry Smith 
4219566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
4228a2871f6SBarry Smith       if (rank == 0) {
4239566063dSJacob Faibussowitsch         PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
4245f80ce2aSJacob Faibussowitsch         PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename);
4258a2871f6SBarry Smith       }
4269566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer)));
4278a2871f6SBarry Smith     }
428792fecdfSBarry Smith     PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id));
4295c6c1daeSBarry Smith     break;
4305c6c1daeSBarry Smith   case FILE_MODE_APPEND:
4319371c9d4SSatish Balay   case FILE_MODE_UPDATE: {
4327e4fd573SVaclav Hapla     PetscBool flg;
4339566063dSJacob Faibussowitsch     PetscCall(PetscTestFile(hdf5->filename, 'r', &flg));
434792fecdfSBarry Smith     if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id));
435792fecdfSBarry Smith     else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id));
4365c6c1daeSBarry Smith     break;
4377e4fd573SVaclav Hapla   }
438d71ae5a4SJacob Faibussowitsch   case FILE_MODE_WRITE:
439d71ae5a4SJacob Faibussowitsch     PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id));
440d71ae5a4SJacob Faibussowitsch     break;
441d71ae5a4SJacob Faibussowitsch   case FILE_MODE_UNDEFINED:
442d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()");
443d71ae5a4SJacob Faibussowitsch   default:
444d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]);
4455c6c1daeSBarry Smith   }
4465f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name);
447792fecdfSBarry Smith   PetscCallHDF5(H5Pclose, (plist_id));
448671abe24SVaclav Hapla   PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer));
4493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4505c6c1daeSBarry Smith }
4515c6c1daeSBarry Smith 
PetscViewerFileGetName_HDF5(PetscViewer viewer,const char ** name)452d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name)
453d71ae5a4SJacob Faibussowitsch {
454d1232d7fSToby Isaac   PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data;
455d1232d7fSToby Isaac 
456d1232d7fSToby Isaac   PetscFunctionBegin;
457d1232d7fSToby Isaac   *name = vhdf5->filename;
4583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
459d1232d7fSToby Isaac }
460d1232d7fSToby Isaac 
PetscViewerSetUp_HDF5(PetscViewer viewer)461d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer)
462d71ae5a4SJacob Faibussowitsch {
4635dc64a97SVaclav Hapla   /*
464b723ab35SVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data;
4655dc64a97SVaclav Hapla   */
466b723ab35SVaclav Hapla 
467b723ab35SVaclav Hapla   PetscFunctionBegin;
4683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
469b723ab35SVaclav Hapla }
470b723ab35SVaclav Hapla 
PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer,PetscBool flg)471d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg)
472d71ae5a4SJacob Faibussowitsch {
47319a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
47419a20e4cSMatthew G. Knepley 
47519a20e4cSMatthew G. Knepley   PetscFunctionBegin;
47619a20e4cSMatthew G. Knepley   hdf5->defTimestepping = flg;
4773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47819a20e4cSMatthew G. Knepley }
47919a20e4cSMatthew G. Knepley 
PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer,PetscBool * flg)480d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg)
481d71ae5a4SJacob Faibussowitsch {
48219a20e4cSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
48319a20e4cSMatthew G. Knepley 
48419a20e4cSMatthew G. Knepley   PetscFunctionBegin;
48519a20e4cSMatthew G. Knepley   *flg = hdf5->defTimestepping;
4863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
48719a20e4cSMatthew G. Knepley }
48819a20e4cSMatthew G. Knepley 
48919a20e4cSMatthew G. Knepley /*@
49019a20e4cSMatthew G. Knepley   PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping
49119a20e4cSMatthew G. Knepley 
492c3339decSBarry Smith   Logically Collective
49319a20e4cSMatthew G. Knepley 
49419a20e4cSMatthew G. Knepley   Input Parameters:
495811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
496811af0c4SBarry Smith - flg    - if `PETSC_TRUE` we will assume that timestepping is on
49719a20e4cSMatthew G. Knepley 
498811af0c4SBarry Smith   Options Database Key:
49919a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping
50019a20e4cSMatthew G. Knepley 
5013f423023SBarry Smith   Level: intermediate
5023f423023SBarry Smith 
503811af0c4SBarry Smith   Note:
50419a20e4cSMatthew G. Knepley   If the timestepping attribute is not found for an object, then the default timestepping is used
50519a20e4cSMatthew G. Knepley 
506d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
50719a20e4cSMatthew G. Knepley @*/
PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer,PetscBool flg)508d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg)
509d71ae5a4SJacob Faibussowitsch {
51019a20e4cSMatthew G. Knepley   PetscFunctionBegin;
51119a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
512cac4c232SBarry Smith   PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg));
5133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51419a20e4cSMatthew G. Knepley }
51519a20e4cSMatthew G. Knepley 
51619a20e4cSMatthew G. Knepley /*@
51719a20e4cSMatthew G. Knepley   PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping
51819a20e4cSMatthew G. Knepley 
51920f4b53cSBarry Smith   Not Collective
52019a20e4cSMatthew G. Knepley 
52119a20e4cSMatthew G. Knepley   Input Parameter:
522811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
52319a20e4cSMatthew G. Knepley 
52419a20e4cSMatthew G. Knepley   Output Parameter:
525811af0c4SBarry Smith . flg - if `PETSC_TRUE` we will assume that timestepping is on
52619a20e4cSMatthew G. Knepley 
52719a20e4cSMatthew G. Knepley   Level: intermediate
52819a20e4cSMatthew G. Knepley 
529d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()`
53019a20e4cSMatthew G. Knepley @*/
PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer,PetscBool * flg)531d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg)
532d71ae5a4SJacob Faibussowitsch {
53319a20e4cSMatthew G. Knepley   PetscFunctionBegin;
53419a20e4cSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
535cac4c232SBarry Smith   PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg));
5363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
53719a20e4cSMatthew G. Knepley }
53819a20e4cSMatthew G. Knepley 
PetscViewerHDF5SetCompress_HDF5(PetscViewer viewer,PetscBool flg)539ea87367fSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5SetCompress_HDF5(PetscViewer viewer, PetscBool flg)
540ea87367fSMatthew G. Knepley {
541ea87367fSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
542ea87367fSMatthew G. Knepley 
543ea87367fSMatthew G. Knepley   PetscFunctionBegin;
544ea87367fSMatthew G. Knepley   hdf5->compress = flg;
545ea87367fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
546ea87367fSMatthew G. Knepley }
547ea87367fSMatthew G. Knepley 
PetscViewerHDF5GetCompress_HDF5(PetscViewer viewer,PetscBool * flg)548ea87367fSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5GetCompress_HDF5(PetscViewer viewer, PetscBool *flg)
549ea87367fSMatthew G. Knepley {
550ea87367fSMatthew G. Knepley   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
551ea87367fSMatthew G. Knepley 
552ea87367fSMatthew G. Knepley   PetscFunctionBegin;
553ea87367fSMatthew G. Knepley   *flg = hdf5->compress;
554ea87367fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
555ea87367fSMatthew G. Knepley }
556ea87367fSMatthew G. Knepley 
557ea87367fSMatthew G. Knepley /*@
558ea87367fSMatthew G. Knepley   PetscViewerHDF5SetCompress - Set the flag for compression
559ea87367fSMatthew G. Knepley 
560ea87367fSMatthew G. Knepley   Logically Collective
561ea87367fSMatthew G. Knepley 
562ea87367fSMatthew G. Knepley   Input Parameters:
563ea87367fSMatthew G. Knepley + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored
564ea87367fSMatthew G. Knepley - flg    - if `PETSC_TRUE` we will turn on compression
565ea87367fSMatthew G. Knepley 
566ea87367fSMatthew G. Knepley   Options Database Key:
567ea87367fSMatthew G. Knepley . -viewer_hdf5_compress - turns on (true) or off (false) compression
568ea87367fSMatthew G. Knepley 
569ea87367fSMatthew G. Knepley   Level: intermediate
570ea87367fSMatthew G. Knepley 
571ea87367fSMatthew G. Knepley .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCompress()`
572ea87367fSMatthew G. Knepley @*/
PetscViewerHDF5SetCompress(PetscViewer viewer,PetscBool flg)573ea87367fSMatthew G. Knepley PetscErrorCode PetscViewerHDF5SetCompress(PetscViewer viewer, PetscBool flg)
574ea87367fSMatthew G. Knepley {
575ea87367fSMatthew G. Knepley   PetscFunctionBegin;
576ea87367fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
577ea87367fSMatthew G. Knepley   PetscTryMethod(viewer, "PetscViewerHDF5SetCompress_C", (PetscViewer, PetscBool), (viewer, flg));
578ea87367fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
579ea87367fSMatthew G. Knepley }
580ea87367fSMatthew G. Knepley 
581ea87367fSMatthew G. Knepley /*@
582ea87367fSMatthew G. Knepley   PetscViewerHDF5GetCompress - Get the flag for compression
583ea87367fSMatthew G. Knepley 
584ea87367fSMatthew G. Knepley   Not Collective
585ea87367fSMatthew G. Knepley 
586ea87367fSMatthew G. Knepley   Input Parameter:
587ea87367fSMatthew G. Knepley . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
588ea87367fSMatthew G. Knepley 
589ea87367fSMatthew G. Knepley   Output Parameter:
590ea87367fSMatthew G. Knepley . flg - if `PETSC_TRUE` we will turn on compression
591ea87367fSMatthew G. Knepley 
592ea87367fSMatthew G. Knepley   Level: intermediate
593ea87367fSMatthew G. Knepley 
594ea87367fSMatthew G. Knepley .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCompress()`
595ea87367fSMatthew G. Knepley @*/
PetscViewerHDF5GetCompress(PetscViewer viewer,PetscBool * flg)596ea87367fSMatthew G. Knepley PetscErrorCode PetscViewerHDF5GetCompress(PetscViewer viewer, PetscBool *flg)
597ea87367fSMatthew G. Knepley {
598ea87367fSMatthew G. Knepley   PetscFunctionBegin;
599ea87367fSMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
600ea87367fSMatthew G. Knepley   PetscUseMethod(viewer, "PetscViewerHDF5GetCompress_C", (PetscViewer, PetscBool *), (viewer, flg));
601ea87367fSMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
602ea87367fSMatthew G. Knepley }
603ea87367fSMatthew G. Knepley 
6048556b5ebSBarry Smith /*MC
6058556b5ebSBarry Smith    PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file
6068556b5ebSBarry Smith 
607811af0c4SBarry Smith   Level: beginner
608811af0c4SBarry Smith 
609d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
610db781477SPatrick Sanan           `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
611db781477SPatrick Sanan           `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
612db781477SPatrick Sanan           `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
6138556b5ebSBarry Smith M*/
614d1232d7fSToby Isaac 
PetscViewerCreate_HDF5(PetscViewer v)615d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v)
616d71ae5a4SJacob Faibussowitsch {
6175c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5;
6185c6c1daeSBarry Smith 
6195c6c1daeSBarry Smith   PetscFunctionBegin;
62099335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL)
62199335c34SVaclav Hapla   {
62299335c34SVaclav Hapla     PetscMPIInt size;
6239566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
6245f80ce2aSJacob 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)");
62599335c34SVaclav Hapla   }
62699335c34SVaclav Hapla #endif
62799335c34SVaclav Hapla 
6284dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&hdf5));
6295c6c1daeSBarry Smith 
6305c6c1daeSBarry Smith   v->data                = (void *)hdf5;
6315c6c1daeSBarry Smith   v->ops->destroy        = PetscViewerDestroy_HDF5;
63282ea9b62SBarry Smith   v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5;
633b723ab35SVaclav Hapla   v->ops->setup          = PetscViewerSetUp_HDF5;
6341b793a25SVaclav Hapla   v->ops->view           = PetscViewerView_HDF5;
6356226335aSVaclav Hapla   v->ops->flush          = PetscViewerFlush_HDF5;
6367e4fd573SVaclav Hapla   hdf5->btype            = FILE_MODE_UNDEFINED;
637908793a3SLisandro Dalcin   hdf5->filename         = NULL;
6385c6c1daeSBarry Smith   hdf5->timestep         = -1;
6390298fd71SBarry Smith   hdf5->groups           = NULL;
6405c6c1daeSBarry Smith 
641792fecdfSBarry Smith   PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER));
6429c5072fbSVaclav Hapla 
6439566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5));
6449566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5));
6459566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5));
6469566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5));
6479566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5));
6489566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5));
6499566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5));
6509566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5));
6519566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5));
6529566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5));
653ea87367fSMatthew G. Knepley   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCompress_C", PetscViewerHDF5GetCompress_HDF5));
654ea87367fSMatthew G. Knepley   PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCompress_C", PetscViewerHDF5SetCompress_HDF5));
6553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6565c6c1daeSBarry Smith }
6575c6c1daeSBarry Smith 
6585c6c1daeSBarry Smith /*@C
659811af0c4SBarry Smith   PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer`
6605c6c1daeSBarry Smith 
661d083f849SBarry Smith   Collective
6625c6c1daeSBarry Smith 
6635c6c1daeSBarry Smith   Input Parameters:
6645c6c1daeSBarry Smith + comm - MPI communicator
6655c6c1daeSBarry Smith . name - name of file
6665c6c1daeSBarry Smith - type - type of file
6675c6c1daeSBarry Smith 
6685c6c1daeSBarry Smith   Output Parameter:
669811af0c4SBarry Smith . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file
6705c6c1daeSBarry Smith 
671811af0c4SBarry Smith   Options Database Keys:
672a2b725a8SWilliam 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
673a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output       - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal
67482ea9b62SBarry Smith 
6755c6c1daeSBarry Smith   Level: beginner
6765c6c1daeSBarry Smith 
6777e4fd573SVaclav Hapla   Notes:
6787e4fd573SVaclav Hapla   Reading is always available, regardless of the mode. Available modes are
679811af0c4SBarry Smith .vb
680811af0c4SBarry Smith   FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY]
681811af0c4SBarry Smith   FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC]
682811af0c4SBarry 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]
683811af0c4SBarry Smith   FILE_MODE_UPDATE - same as FILE_MODE_APPEND
684811af0c4SBarry Smith .ve
6857e4fd573SVaclav Hapla 
686da81f932SPierre Jolivet   In case of `FILE_MODE_APPEND` / `FILE_MODE_UPDATE`, any stored object (dataset, attribute) can be selectively overwritten if the same fully qualified name (/group/path/to/object) is specified.
6877e4fd573SVaclav Hapla 
6885c6c1daeSBarry Smith   This PetscViewer should be destroyed with PetscViewerDestroy().
6895c6c1daeSBarry Smith 
690d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`,
691db781477SPatrick Sanan           `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`,
692db781477SPatrick Sanan           `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()`
6935c6c1daeSBarry Smith @*/
PetscViewerHDF5Open(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer * hdf5v)694d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v)
695d71ae5a4SJacob Faibussowitsch {
6965c6c1daeSBarry Smith   PetscFunctionBegin;
6979566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, hdf5v));
6989566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5));
6999566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(*hdf5v, type));
7009566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(*hdf5v, name));
7019566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetFromOptions(*hdf5v));
7023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7035c6c1daeSBarry Smith }
7045c6c1daeSBarry Smith 
7055c6c1daeSBarry Smith /*@C
7065c6c1daeSBarry Smith   PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls
7075c6c1daeSBarry Smith 
70820f4b53cSBarry Smith   Not Collective
7095c6c1daeSBarry Smith 
7105c6c1daeSBarry Smith   Input Parameter:
711811af0c4SBarry Smith . viewer - the `PetscViewer`
7125c6c1daeSBarry Smith 
7135c6c1daeSBarry Smith   Output Parameter:
7145c6c1daeSBarry Smith . file_id - The file id
7155c6c1daeSBarry Smith 
7165c6c1daeSBarry Smith   Level: intermediate
7175c6c1daeSBarry Smith 
718d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`
7195c6c1daeSBarry Smith @*/
PetscViewerHDF5GetFileId(PetscViewer viewer,hid_t * file_id)720d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id)
721d71ae5a4SJacob Faibussowitsch {
7225c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7235c6c1daeSBarry Smith 
7245c6c1daeSBarry Smith   PetscFunctionBegin;
7255c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7265c6c1daeSBarry Smith   if (file_id) *file_id = hdf5->file_id;
7273ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7285c6c1daeSBarry Smith }
7295c6c1daeSBarry Smith 
7305c6c1daeSBarry Smith /*@C
7315c6c1daeSBarry Smith   PetscViewerHDF5PushGroup - Set the current HDF5 group for output
7325c6c1daeSBarry Smith 
73320f4b53cSBarry Smith   Not Collective
7345c6c1daeSBarry Smith 
7355c6c1daeSBarry Smith   Input Parameters:
736811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
7375c6c1daeSBarry Smith - name   - The group name
7385c6c1daeSBarry Smith 
7395c6c1daeSBarry Smith   Level: intermediate
7405c6c1daeSBarry Smith 
7412e361344SVaclav Hapla   Notes:
7422e361344SVaclav Hapla   This is designed to mnemonically resemble the Unix cd command.
743811af0c4SBarry Smith .vb
744811af0c4SBarry 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.
7453f423023SBarry Smith   `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/".
746811af0c4SBarry Smith   "." means the current group is pushed again.
747811af0c4SBarry Smith .ve
7482e361344SVaclav Hapla 
7492e361344SVaclav Hapla   Example:
7502e361344SVaclav Hapla   Suppose the current group is "/a".
7513f423023SBarry Smith .vb
7523f423023SBarry Smith   If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/".
7533f423023SBarry Smith   If name is ".", then the new group will be "/a".
7543f423023SBarry Smith   If name is "b", then the new group will be "/a/b".
7553f423023SBarry Smith   If name is "/b", then the new group will be "/b".
7563f423023SBarry Smith .ve
7572e361344SVaclav Hapla 
758aec76313SJacob Faibussowitsch   Developer Notes:
7593f423023SBarry Smith   The root group "/" is internally stored as `NULL`.
760820fc2d1SVaclav Hapla 
761d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
7625c6c1daeSBarry Smith @*/
PetscViewerHDF5PushGroup(PetscViewer viewer,const char name[])763d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[])
764d71ae5a4SJacob Faibussowitsch {
7655c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
7667d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
7672e361344SVaclav Hapla   size_t                    i, len;
7682e361344SVaclav Hapla   char                      buf[PETSC_MAX_PATH_LEN];
7692e361344SVaclav Hapla   const char               *gname;
7705c6c1daeSBarry Smith 
7715c6c1daeSBarry Smith   PetscFunctionBegin;
7725c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
7734f572ea9SToby Isaac   if (name) PetscAssertPointer(name, 2);
7749566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(name, &len));
7752e361344SVaclav Hapla   gname = NULL;
7762e361344SVaclav Hapla   if (len) {
7772e361344SVaclav Hapla     if (len == 1 && name[0] == '.') {
7782e361344SVaclav Hapla       /* use current name */
7792e361344SVaclav Hapla       gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL;
7802e361344SVaclav Hapla     } else if (name[0] == '/') {
7812e361344SVaclav Hapla       /* absolute */
7822e361344SVaclav Hapla       for (i = 1; i < len; i++) {
7832e361344SVaclav Hapla         if (name[i] != '/') {
7842e361344SVaclav Hapla           gname = name;
7852e361344SVaclav Hapla           break;
7862e361344SVaclav Hapla         }
7872e361344SVaclav Hapla       }
7882e361344SVaclav Hapla     } else {
7892e361344SVaclav Hapla       /* relative */
7902e361344SVaclav Hapla       const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : "";
7919566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name));
7922e361344SVaclav Hapla       gname = buf;
7932e361344SVaclav Hapla     }
7942e361344SVaclav Hapla   }
7959566063dSJacob Faibussowitsch   PetscCall(PetscNew(&groupNode));
7969566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name));
7975c6c1daeSBarry Smith   groupNode->next = hdf5->groups;
7985c6c1daeSBarry Smith   hdf5->groups    = groupNode;
7993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8005c6c1daeSBarry Smith }
8015c6c1daeSBarry Smith 
8023ef9c667SSatish Balay /*@
8035c6c1daeSBarry Smith   PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value
8045c6c1daeSBarry Smith 
80520f4b53cSBarry Smith   Not Collective
8065c6c1daeSBarry Smith 
8075c6c1daeSBarry Smith   Input Parameter:
808811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
8095c6c1daeSBarry Smith 
8105c6c1daeSBarry Smith   Level: intermediate
8115c6c1daeSBarry Smith 
812d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()`
8135c6c1daeSBarry Smith @*/
PetscViewerHDF5PopGroup(PetscViewer viewer)814d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer)
815d71ae5a4SJacob Faibussowitsch {
8165c6c1daeSBarry Smith   PetscViewer_HDF5         *hdf5 = (PetscViewer_HDF5 *)viewer->data;
8177d964842SVaclav Hapla   PetscViewerHDF5GroupList *groupNode;
8185c6c1daeSBarry Smith 
8195c6c1daeSBarry Smith   PetscFunctionBegin;
8205c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8215f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop");
8225c6c1daeSBarry Smith   groupNode    = hdf5->groups;
8235c6c1daeSBarry Smith   hdf5->groups = hdf5->groups->next;
8249566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode->name));
8259566063dSJacob Faibussowitsch   PetscCall(PetscFree(groupNode));
8263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8275c6c1daeSBarry Smith }
8285c6c1daeSBarry Smith 
PetscViewerHDF5GetGroup_Internal(PetscViewer viewer,const char * name[])82934e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[])
830d71ae5a4SJacob Faibussowitsch {
8315c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
8325c6c1daeSBarry Smith 
8335c6c1daeSBarry Smith   PetscFunctionBegin;
8345c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
8354f572ea9SToby Isaac   PetscAssertPointer(name, 2);
836a297a907SKarl Rupp   if (hdf5->groups) *name = hdf5->groups->name;
8370298fd71SBarry Smith   else *name = NULL;
8383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8395c6c1daeSBarry Smith }
8405c6c1daeSBarry Smith 
8413014b61aSVaclav Hapla /*@C
842811af0c4SBarry Smith   PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`,
843874270d9SVaclav Hapla   and return this group's ID and file ID.
844811af0c4SBarry Smith   If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID.
845874270d9SVaclav Hapla 
84620f4b53cSBarry Smith   Not Collective
847874270d9SVaclav Hapla 
8483014b61aSVaclav Hapla   Input Parameters:
8493014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
8503014b61aSVaclav Hapla - path   - (Optional) The path relative to the pushed group
851874270d9SVaclav Hapla 
852d8d19677SJose E. Roman   Output Parameters:
853874270d9SVaclav Hapla + fileId  - The HDF5 file ID
854874270d9SVaclav Hapla - groupId - The HDF5 group ID
855874270d9SVaclav Hapla 
8563f423023SBarry Smith   Level: intermediate
8573f423023SBarry Smith 
858811af0c4SBarry Smith   Note:
8593014b61aSVaclav 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.
8603f423023SBarry Smith   `NULL` or empty path means the current pushed group.
8613014b61aSVaclav Hapla 
862e74616beSVaclav Hapla   If the viewer is writable, the group is created if it doesn't exist yet.
863e74616beSVaclav Hapla 
864d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()`
865874270d9SVaclav Hapla @*/
PetscViewerHDF5OpenGroup(PetscViewer viewer,const char path[],hid_t * fileId,hid_t * groupId)866d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId)
867d71ae5a4SJacob Faibussowitsch {
86890e03537SVaclav Hapla   hid_t       file_id;
86990e03537SVaclav Hapla   H5O_type_t  type;
8703014b61aSVaclav Hapla   const char *fileName  = NULL;
871ce78bad3SBarry Smith   const char *groupName = NULL;
87276d59af2SVaclav Hapla   PetscBool   writable, has;
87354dbf706SVaclav Hapla 
87454dbf706SVaclav Hapla   PetscFunctionBegin;
8759566063dSJacob Faibussowitsch   PetscCall(PetscViewerWritable(viewer, &writable));
8769566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id));
8779566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileGetName(viewer, &fileName));
8783014b61aSVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName));
8799566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type));
88076d59af2SVaclav Hapla   if (!has) {
8815f80ce2aSJacob 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);
882f7d195e4SLawrence Mitchell     SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName);
88376d59af2SVaclav Hapla   }
8845f80ce2aSJacob 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);
8853014b61aSVaclav Hapla   PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT));
8863014b61aSVaclav Hapla   PetscCall(PetscFree(groupName));
88754dbf706SVaclav Hapla   *fileId = file_id;
8883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
88954dbf706SVaclav Hapla }
89054dbf706SVaclav Hapla 
89163cb69f5SVaclav Hapla /*@C
89263cb69f5SVaclav Hapla   PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file
89363cb69f5SVaclav Hapla 
89420f4b53cSBarry Smith   Not Collective
89563cb69f5SVaclav Hapla 
89663cb69f5SVaclav Hapla   Input Parameters:
89763cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
89863cb69f5SVaclav Hapla - path   - (Optional) The path relative to the pushed group
89963cb69f5SVaclav Hapla 
9003f423023SBarry Smith   Level: intermediate
9013f423023SBarry Smith 
90263cb69f5SVaclav Hapla   Note:
90363cb69f5SVaclav 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.
9043f423023SBarry Smith   `NULL` or empty path means the current pushed group.
90563cb69f5SVaclav Hapla 
90663cb69f5SVaclav Hapla   This will fail if the viewer is not writable.
90763cb69f5SVaclav Hapla 
908d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
90963cb69f5SVaclav Hapla @*/
PetscViewerHDF5WriteGroup(PetscViewer viewer,const char path[])910d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[])
911d71ae5a4SJacob Faibussowitsch {
91263cb69f5SVaclav Hapla   hid_t fileId, groupId;
91363cb69f5SVaclav Hapla 
91463cb69f5SVaclav Hapla   PetscFunctionBegin;
91563cb69f5SVaclav Hapla   PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created
91663cb69f5SVaclav Hapla   PetscCallHDF5(H5Gclose, (groupId));
9173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91863cb69f5SVaclav Hapla }
91963cb69f5SVaclav Hapla 
9203ef9c667SSatish Balay /*@
921d7dd068bSVaclav Hapla   PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing.
9225c6c1daeSBarry Smith 
92320f4b53cSBarry Smith   Not Collective
9245c6c1daeSBarry Smith 
9255c6c1daeSBarry Smith   Input Parameter:
926811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
9275c6c1daeSBarry Smith 
9285c6c1daeSBarry Smith   Level: intermediate
9295c6c1daeSBarry Smith 
930d7dd068bSVaclav Hapla   Notes:
931811af0c4SBarry Smith   On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0.
932811af0c4SBarry Smith   Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`.
933811af0c4SBarry 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()`].
934811af0c4SBarry Smith   Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory.
935811af0c4SBarry Smith   Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`.
936d7dd068bSVaclav Hapla 
937d7dd068bSVaclav Hapla   If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again.
938d7dd068bSVaclav Hapla   Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error.
939d7dd068bSVaclav Hapla 
940aec76313SJacob Faibussowitsch   Developer Notes:
941d7dd068bSVaclav Hapla   Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true.
942d7dd068bSVaclav Hapla 
943d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
944d7dd068bSVaclav Hapla @*/
PetscViewerHDF5PushTimestepping(PetscViewer viewer)945d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer)
946d71ae5a4SJacob Faibussowitsch {
947d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
948d7dd068bSVaclav Hapla 
949d7dd068bSVaclav Hapla   PetscFunctionBegin;
950d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9515f80ce2aSJacob Faibussowitsch   PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed");
952d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_TRUE;
953d7dd068bSVaclav Hapla   if (hdf5->timestep < 0) hdf5->timestep = 0;
9543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
955d7dd068bSVaclav Hapla }
956d7dd068bSVaclav Hapla 
957d7dd068bSVaclav Hapla /*@
958d7dd068bSVaclav Hapla   PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing.
959d7dd068bSVaclav Hapla 
96020f4b53cSBarry Smith   Not Collective
961d7dd068bSVaclav Hapla 
962d7dd068bSVaclav Hapla   Input Parameter:
963811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
964d7dd068bSVaclav Hapla 
965d7dd068bSVaclav Hapla   Level: intermediate
966d7dd068bSVaclav Hapla 
967811af0c4SBarry Smith   Note:
968811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
969d7dd068bSVaclav Hapla 
970d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
971d7dd068bSVaclav Hapla @*/
PetscViewerHDF5PopTimestepping(PetscViewer viewer)972d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer)
973d71ae5a4SJacob Faibussowitsch {
974d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
975d7dd068bSVaclav Hapla 
976d7dd068bSVaclav Hapla   PetscFunctionBegin;
977d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
9785f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
979d7dd068bSVaclav Hapla   hdf5->timestepping = PETSC_FALSE;
9803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
981d7dd068bSVaclav Hapla }
982d7dd068bSVaclav Hapla 
983d7dd068bSVaclav Hapla /*@
984d7dd068bSVaclav Hapla   PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently.
985d7dd068bSVaclav Hapla 
98620f4b53cSBarry Smith   Not Collective
987d7dd068bSVaclav Hapla 
988d7dd068bSVaclav Hapla   Input Parameter:
989811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
990d7dd068bSVaclav Hapla 
991d7dd068bSVaclav Hapla   Output Parameter:
992d7dd068bSVaclav Hapla . flg - is timestepping active?
993d7dd068bSVaclav Hapla 
994d7dd068bSVaclav Hapla   Level: intermediate
995d7dd068bSVaclav Hapla 
996811af0c4SBarry Smith   Note:
997811af0c4SBarry Smith   See `PetscViewerHDF5PushTimestepping()` for details.
998d7dd068bSVaclav Hapla 
999d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
1000d7dd068bSVaclav Hapla @*/
PetscViewerHDF5IsTimestepping(PetscViewer viewer,PetscBool * flg)1001d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg)
1002d71ae5a4SJacob Faibussowitsch {
1003d7dd068bSVaclav Hapla   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
1004d7dd068bSVaclav Hapla 
1005d7dd068bSVaclav Hapla   PetscFunctionBegin;
1006d7dd068bSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1007d7dd068bSVaclav Hapla   *flg = hdf5->timestepping;
10083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1009d7dd068bSVaclav Hapla }
1010d7dd068bSVaclav Hapla 
1011d7dd068bSVaclav Hapla /*@
1012d7dd068bSVaclav Hapla   PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time.
1013d7dd068bSVaclav Hapla 
101420f4b53cSBarry Smith   Not Collective
1015d7dd068bSVaclav Hapla 
1016d7dd068bSVaclav Hapla   Input Parameter:
1017811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
1018d7dd068bSVaclav Hapla 
1019d7dd068bSVaclav Hapla   Level: intermediate
1020d7dd068bSVaclav Hapla 
1021811af0c4SBarry Smith   Note:
1022811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1023d7dd068bSVaclav Hapla 
1024d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()`
10255c6c1daeSBarry Smith @*/
PetscViewerHDF5IncrementTimestep(PetscViewer viewer)1026d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer)
1027d71ae5a4SJacob Faibussowitsch {
10285c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
10295c6c1daeSBarry Smith 
10305c6c1daeSBarry Smith   PetscFunctionBegin;
10315c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
10325f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
10335c6c1daeSBarry Smith   ++hdf5->timestep;
10343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10355c6c1daeSBarry Smith }
10365c6c1daeSBarry Smith 
10373ef9c667SSatish Balay /*@
1038d7dd068bSVaclav Hapla   PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time.
10395c6c1daeSBarry Smith 
104020f4b53cSBarry Smith   Logically Collective
10415c6c1daeSBarry Smith 
10425c6c1daeSBarry Smith   Input Parameters:
1043811af0c4SBarry Smith + viewer   - the `PetscViewer` of type `PETSCVIEWERHDF5`
1044d7dd068bSVaclav Hapla - timestep - The timestep
10455c6c1daeSBarry Smith 
10465c6c1daeSBarry Smith   Level: intermediate
10475c6c1daeSBarry Smith 
1048811af0c4SBarry Smith   Note:
1049811af0c4SBarry Smith   This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1050d7dd068bSVaclav Hapla 
1051d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()`
10525c6c1daeSBarry Smith @*/
PetscViewerHDF5SetTimestep(PetscViewer viewer,PetscInt timestep)1053d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep)
1054d71ae5a4SJacob Faibussowitsch {
10555c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
10565c6c1daeSBarry Smith 
10575c6c1daeSBarry Smith   PetscFunctionBegin;
10585c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1059d7dd068bSVaclav Hapla   PetscValidLogicalCollectiveInt(viewer, timestep, 2);
10605f80ce2aSJacob Faibussowitsch   PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep);
10615f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
10625c6c1daeSBarry Smith   hdf5->timestep = timestep;
10633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10645c6c1daeSBarry Smith }
10655c6c1daeSBarry Smith 
10663ef9c667SSatish Balay /*@
10675c6c1daeSBarry Smith   PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time.
10685c6c1daeSBarry Smith 
106920f4b53cSBarry Smith   Not Collective
10705c6c1daeSBarry Smith 
10715c6c1daeSBarry Smith   Input Parameter:
1072811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5`
10735c6c1daeSBarry Smith 
10745c6c1daeSBarry Smith   Output Parameter:
1075d7dd068bSVaclav Hapla . timestep - The timestep
10765c6c1daeSBarry Smith 
10775c6c1daeSBarry Smith   Level: intermediate
10785c6c1daeSBarry Smith 
1079811af0c4SBarry Smith   Note:
1080811af0c4SBarry Smith   This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details.
1081d7dd068bSVaclav Hapla 
1082d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()`
10835c6c1daeSBarry Smith @*/
PetscViewerHDF5GetTimestep(PetscViewer viewer,PetscInt * timestep)1084d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep)
1085d71ae5a4SJacob Faibussowitsch {
10865c6c1daeSBarry Smith   PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data;
10875c6c1daeSBarry Smith 
10885c6c1daeSBarry Smith   PetscFunctionBegin;
10895c6c1daeSBarry Smith   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
10904f572ea9SToby Isaac   PetscAssertPointer(timestep, 2);
10915f80ce2aSJacob Faibussowitsch   PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first");
10925c6c1daeSBarry Smith   *timestep = hdf5->timestep;
10933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10945c6c1daeSBarry Smith }
10955c6c1daeSBarry Smith 
109636bce6e8SMatthew G. Knepley /*@C
109736bce6e8SMatthew G. Knepley   PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name.
109836bce6e8SMatthew G. Knepley 
109920f4b53cSBarry Smith   Not Collective
110036bce6e8SMatthew G. Knepley 
110136bce6e8SMatthew G. Knepley   Input Parameter:
1102811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
110336bce6e8SMatthew G. Knepley 
110436bce6e8SMatthew G. Knepley   Output Parameter:
1105aec76313SJacob Faibussowitsch . htype - the HDF5  datatype
110636bce6e8SMatthew G. Knepley 
110736bce6e8SMatthew G. Knepley   Level: advanced
110836bce6e8SMatthew G. Knepley 
1109d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()`
111036bce6e8SMatthew G. Knepley @*/
PetscDataTypeToHDF5DataType(PetscDataType ptype,hid_t * htype)1111d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype)
1112d71ae5a4SJacob Faibussowitsch {
111336bce6e8SMatthew G. Knepley   PetscFunctionBegin;
1114*fc2fb351SPierre Jolivet   if (ptype == PETSC_INT) *htype = PetscDefined(USE_64BIT_INDICES) ? H5T_NATIVE_LLONG : H5T_NATIVE_INT;
111536bce6e8SMatthew G. Knepley   else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE;
111636bce6e8SMatthew G. Knepley   else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG;
111736bce6e8SMatthew G. Knepley   else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT;
1118c9450970SBarry Smith   else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT;
11191dc74096SMartin Diehl   else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_HBOOL;
112036bce6e8SMatthew G. Knepley   else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT;
112136bce6e8SMatthew G. Knepley   else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR;
112236bce6e8SMatthew G. Knepley   else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR;
11237e97c476SMatthew G. Knepley   else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1);
112436bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype");
11253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
112636bce6e8SMatthew G. Knepley }
112736bce6e8SMatthew G. Knepley 
112836bce6e8SMatthew G. Knepley /*@C
112936bce6e8SMatthew G. Knepley   PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name
113036bce6e8SMatthew G. Knepley 
113120f4b53cSBarry Smith   Not Collective
113236bce6e8SMatthew G. Knepley 
113336bce6e8SMatthew G. Knepley   Input Parameter:
1134811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...)
113536bce6e8SMatthew G. Knepley 
113636bce6e8SMatthew G. Knepley   Output Parameter:
1137811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`)
113836bce6e8SMatthew G. Knepley 
113936bce6e8SMatthew G. Knepley   Level: advanced
114036bce6e8SMatthew G. Knepley 
114142747ad1SJacob Faibussowitsch .seealso: [](sec_viewers), `PetscDataType`
114236bce6e8SMatthew G. Knepley @*/
PetscHDF5DataTypeToPetscDataType(hid_t htype,PetscDataType * ptype)1143d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype)
1144d71ae5a4SJacob Faibussowitsch {
114536bce6e8SMatthew G. Knepley   PetscFunctionBegin;
1146*fc2fb351SPierre Jolivet   if (htype == H5T_NATIVE_INT) *ptype = PetscDefined(USE_64BIT_INDICES) ? PETSC_LONG : PETSC_INT;
114736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES)
114836bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT;
114936bce6e8SMatthew G. Knepley #endif
115036bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE;
115136bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG;
115236bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT;
11531dc74096SMartin Diehl   else if (htype == H5T_NATIVE_HBOOL) *ptype = PETSC_BOOL;
115436bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT;
115536bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR;
115636bce6e8SMatthew G. Knepley   else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR;
11577e97c476SMatthew G. Knepley   else if (htype == H5T_C_S1) *ptype = PETSC_STRING;
115836bce6e8SMatthew G. Knepley   else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype");
11593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
116036bce6e8SMatthew G. Knepley }
116136bce6e8SMatthew G. Knepley 
1162df863907SAlex Fikl /*@C
1163b8ee85a1SVaclav Hapla   PetscViewerHDF5WriteAttribute - Write an attribute
116436bce6e8SMatthew G. Knepley 
1165343a20b2SVaclav Hapla   Collective
1166343a20b2SVaclav Hapla 
116736bce6e8SMatthew G. Knepley   Input Parameters:
1168811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
11694302210dSVaclav Hapla . parent   - The parent dataset/group name
117036bce6e8SMatthew G. Knepley . name     - The attribute name
117136bce6e8SMatthew G. Knepley . datatype - The attribute type
117236bce6e8SMatthew G. Knepley - value    - The attribute value
117336bce6e8SMatthew G. Knepley 
117436bce6e8SMatthew G. Knepley   Level: advanced
117536bce6e8SMatthew G. Knepley 
1176811af0c4SBarry Smith   Note:
1177343a20b2SVaclav 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.
11784302210dSVaclav Hapla 
1179d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`,
1180811af0c4SBarry Smith           `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
118136bce6e8SMatthew G. Knepley @*/
PetscViewerHDF5WriteAttribute(PetscViewer viewer,const char parent[],const char name[],PetscDataType datatype,const void * value)1182d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value)
1183d71ae5a4SJacob Faibussowitsch {
1184ce78bad3SBarry Smith   const char *parentAbsPath;
118560568592SMatthew G. Knepley   hid_t       h5, dataspace, obj, attribute, dtype;
1186080f144cSVaclav Hapla   PetscBool   has;
118736bce6e8SMatthew G. Knepley 
118836bce6e8SMatthew G. Knepley   PetscFunctionBegin;
11895cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
11904f572ea9SToby Isaac   if (parent) PetscAssertPointer(parent, 2);
11914f572ea9SToby Isaac   PetscAssertPointer(name, 3);
11924302210dSVaclav Hapla   PetscValidLogicalCollectiveEnum(viewer, datatype, 4);
11934f572ea9SToby Isaac   PetscAssertPointer(value, 5);
119477717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
11959566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL));
11969566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
11979566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
11987e97c476SMatthew G. Knepley   if (datatype == PETSC_STRING) {
11997e97c476SMatthew G. Knepley     size_t len;
12009566063dSJacob Faibussowitsch     PetscCall(PetscStrlen((const char *)value, &len));
1201792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
12027e97c476SMatthew G. Knepley   }
12039566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1204792fecdfSBarry Smith   PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR));
1205792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1206080f144cSVaclav Hapla   if (has) {
1207792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1208080f144cSVaclav Hapla   } else {
1209792fecdfSBarry Smith     PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT));
1210080f144cSVaclav Hapla   }
1211792fecdfSBarry Smith   PetscCallHDF5(H5Awrite, (attribute, dtype, value));
1212792fecdfSBarry Smith   if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype));
1213792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1214792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
1215792fecdfSBarry Smith   PetscCallHDF5(H5Sclose, (dataspace));
12169566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
12173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121836bce6e8SMatthew G. Knepley }
121936bce6e8SMatthew G. Knepley 
1220df863907SAlex Fikl /*@C
1221811af0c4SBarry Smith   PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name
1222577e0e04SVaclav Hapla 
1223343a20b2SVaclav Hapla   Collective
1224343a20b2SVaclav Hapla 
1225577e0e04SVaclav Hapla   Input Parameters:
1226811af0c4SBarry Smith + viewer   - The `PETSCVIEWERHDF5` viewer
1227577e0e04SVaclav Hapla . obj      - The object whose name is used to lookup the parent dataset, relative to the current group.
1228577e0e04SVaclav Hapla . name     - The attribute name
1229577e0e04SVaclav Hapla . datatype - The attribute type
1230577e0e04SVaclav Hapla - value    - The attribute value
1231577e0e04SVaclav Hapla 
12323f423023SBarry Smith   Level: advanced
12333f423023SBarry Smith 
1234811af0c4SBarry Smith   Note:
1235343a20b2SVaclav 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).
1236811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1237577e0e04SVaclav Hapla 
1238d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`,
1239811af0c4SBarry Smith           `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1240577e0e04SVaclav Hapla @*/
PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer,PetscObject obj,const char name[],PetscDataType datatype,const void * value)1241d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value)
1242d71ae5a4SJacob Faibussowitsch {
1243577e0e04SVaclav Hapla   PetscFunctionBegin;
1244577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1245577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
12464f572ea9SToby Isaac   PetscAssertPointer(name, 3);
12474f572ea9SToby Isaac   PetscAssertPointer(value, 5);
12489566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
12499566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value));
12503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1251577e0e04SVaclav Hapla }
1252577e0e04SVaclav Hapla 
1253577e0e04SVaclav Hapla /*@C
1254b8ee85a1SVaclav Hapla   PetscViewerHDF5ReadAttribute - Read an attribute
1255ad0c61feSMatthew G. Knepley 
1256343a20b2SVaclav Hapla   Collective
1257343a20b2SVaclav Hapla 
1258ad0c61feSMatthew G. Knepley   Input Parameters:
1259811af0c4SBarry Smith + viewer       - The `PETSCVIEWERHDF5` viewer
12604302210dSVaclav Hapla . parent       - The parent dataset/group name
1261ad0c61feSMatthew G. Knepley . name         - The attribute name
1262a2d6be1bSVaclav Hapla . datatype     - The attribute type
1263a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value
1264ad0c61feSMatthew G. Knepley 
1265ad0c61feSMatthew G. Knepley   Output Parameter:
1266a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value
1267ad0c61feSMatthew G. Knepley 
12683f423023SBarry Smith   Level: advanced
12693f423023SBarry Smith 
1270a2d6be1bSVaclav Hapla   Notes:
12713f423023SBarry Smith   If defaultValue is `NULL` and the attribute is not found, an error occurs.
12723f423023SBarry Smith 
12733f423023SBarry Smith   If defaultValue is not `NULL` and the attribute is not found, `defaultValue` is copied to value.
12743f423023SBarry Smith 
12753f423023SBarry Smith   The pointers `defaultValue` and value can be the same; for instance
12763f423023SBarry Smith .vb
12773f423023SBarry Smith   flg = PETSC_FALSE;
12783f423023SBarry Smith   PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg));
12793f423023SBarry Smith .ve
1280a2d6be1bSVaclav Hapla   is valid, but make sure the default value is initialized.
1281a2d6be1bSVaclav Hapla 
1282811af0c4SBarry Smith   If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed.
1283708d4cb3SBarry Smith 
12843f423023SBarry Smith   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.
12854302210dSVaclav Hapla 
1286d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1287ad0c61feSMatthew G. Knepley @*/
PetscViewerHDF5ReadAttribute(PetscViewer viewer,const char parent[],const char name[],PetscDataType datatype,const void * defaultValue,void * value)1288d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value)
1289d71ae5a4SJacob Faibussowitsch {
1290ce78bad3SBarry Smith   const char *parentAbsPath;
1291a2d6be1bSVaclav Hapla   hid_t       h5, obj, attribute, dtype;
1292969748fdSVaclav Hapla   PetscBool   has;
1293ad0c61feSMatthew G. Knepley 
1294ad0c61feSMatthew G. Knepley   PetscFunctionBegin;
12955cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
12964f572ea9SToby Isaac   if (parent) PetscAssertPointer(parent, 2);
12974f572ea9SToby Isaac   PetscAssertPointer(name, 3);
12984f572ea9SToby Isaac   if (defaultValue) PetscAssertPointer(defaultValue, 5);
12994f572ea9SToby Isaac   PetscAssertPointer(value, 6);
13009566063dSJacob Faibussowitsch   PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype));
130177717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
13029566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL));
13039566063dSJacob Faibussowitsch   if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has));
1304a2d6be1bSVaclav Hapla   if (!has) {
1305a2d6be1bSVaclav Hapla     if (defaultValue) {
1306a2d6be1bSVaclav Hapla       if (defaultValue != value) {
1307a2d6be1bSVaclav Hapla         if (datatype == PETSC_STRING) {
13089566063dSJacob Faibussowitsch           PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value));
1309a2d6be1bSVaclav Hapla         } else {
1310a2d6be1bSVaclav Hapla           size_t len;
1311792fecdfSBarry Smith           PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype));
13129566063dSJacob Faibussowitsch           PetscCall(PetscMemcpy(value, defaultValue, len));
1313a2d6be1bSVaclav Hapla         }
1314a2d6be1bSVaclav Hapla       }
13159566063dSJacob Faibussowitsch       PetscCall(PetscFree(parentAbsPath));
13163ba16761SJacob Faibussowitsch       PetscFunctionReturn(PETSC_SUCCESS);
131798921bdaSJacob Faibussowitsch     } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name);
1318a2d6be1bSVaclav Hapla   }
13199566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1320792fecdfSBarry Smith   PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT));
1321792fecdfSBarry Smith   PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name));
1322f0b58479SMatthew G. Knepley   if (datatype == PETSC_STRING) {
1323f0b58479SMatthew G. Knepley     size_t len;
1324a2d6be1bSVaclav Hapla     hid_t  atype;
1325792fecdfSBarry Smith     PetscCallHDF5Return(atype, H5Aget_type, (attribute));
1326792fecdfSBarry Smith     PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype));
13279566063dSJacob Faibussowitsch     PetscCall(PetscMalloc((len + 1) * sizeof(char), value));
1328792fecdfSBarry Smith     PetscCallHDF5(H5Tset_size, (dtype, len + 1));
1329792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value));
1330708d4cb3SBarry Smith   } else {
1331792fecdfSBarry Smith     PetscCallHDF5(H5Aread, (attribute, dtype, value));
1332708d4cb3SBarry Smith   }
1333792fecdfSBarry Smith   PetscCallHDF5(H5Aclose, (attribute));
1334e90c5075SMarek Pecha   /* H5Oclose can be used to close groups, datasets, or committed datatypes */
1335792fecdfSBarry Smith   PetscCallHDF5(H5Oclose, (obj));
13369566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
13373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1338ad0c61feSMatthew G. Knepley }
1339ad0c61feSMatthew G. Knepley 
1340577e0e04SVaclav Hapla /*@C
1341811af0c4SBarry Smith   PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name
1342577e0e04SVaclav Hapla 
1343343a20b2SVaclav Hapla   Collective
1344343a20b2SVaclav Hapla 
1345577e0e04SVaclav Hapla   Input Parameters:
1346811af0c4SBarry Smith + viewer       - The `PETSCVIEWERHDF5` viewer
1347577e0e04SVaclav Hapla . obj          - The object whose name is used to lookup the parent dataset, relative to the current group.
1348577e0e04SVaclav Hapla . name         - The attribute name
134910450e9eSJacob Faibussowitsch . datatype     - The attribute type
135010450e9eSJacob Faibussowitsch - defaultValue - The default attribute value
1351577e0e04SVaclav Hapla 
1352577e0e04SVaclav Hapla   Output Parameter:
1353577e0e04SVaclav Hapla . value - The attribute value
1354577e0e04SVaclav Hapla 
13553f423023SBarry Smith   Level: advanced
13563f423023SBarry Smith 
1357811af0c4SBarry Smith   Note:
1358577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1359811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1360577e0e04SVaclav Hapla 
1361d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1362577e0e04SVaclav Hapla @*/
PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer,PetscObject obj,const char name[],PetscDataType datatype,void * defaultValue,void * value)1363d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value)
1364d71ae5a4SJacob Faibussowitsch {
1365577e0e04SVaclav Hapla   PetscFunctionBegin;
1366577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1367577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
13684f572ea9SToby Isaac   PetscAssertPointer(name, 3);
13694f572ea9SToby Isaac   PetscAssertPointer(value, 6);
13709566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
13719566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value));
13723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1373577e0e04SVaclav Hapla }
1374577e0e04SVaclav Hapla 
PetscViewerHDF5Traverse_Inner_Internal(hid_t h5,const char name[],PetscBool createGroup,PetscBool * exists_)137534e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_)
1376d71ae5a4SJacob Faibussowitsch {
1377a75c3fd4SVaclav Hapla   htri_t exists;
1378a75c3fd4SVaclav Hapla   hid_t  group;
1379a75c3fd4SVaclav Hapla 
1380a75c3fd4SVaclav Hapla   PetscFunctionBegin;
1381792fecdfSBarry Smith   PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT));
1382792fecdfSBarry Smith   if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT));
1383a75c3fd4SVaclav Hapla   if (!exists && createGroup) {
1384792fecdfSBarry Smith     PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT));
1385792fecdfSBarry Smith     PetscCallHDF5(H5Gclose, (group));
1386a75c3fd4SVaclav Hapla     exists = PETSC_TRUE;
1387a75c3fd4SVaclav Hapla   }
1388a75c3fd4SVaclav Hapla   *exists_ = (PetscBool)exists;
13893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1390a75c3fd4SVaclav Hapla }
1391a75c3fd4SVaclav Hapla 
PetscViewerHDF5Traverse_Internal(PetscViewer viewer,const char name[],PetscBool createGroup,PetscBool * has,H5O_type_t * otype)1392d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype)
1393d71ae5a4SJacob Faibussowitsch {
139490e03537SVaclav Hapla   const char rootGroupName[] = "/";
13955cdeb108SMatthew G. Knepley   hid_t      h5;
1396e5bf9ebcSVaclav Hapla   PetscBool  exists = PETSC_FALSE;
139715b861d2SVaclav Hapla   PetscInt   i;
139815b861d2SVaclav Hapla   int        n;
139985688ae2SVaclav Hapla   char     **hierarchy;
140085688ae2SVaclav Hapla   char       buf[PETSC_MAX_PATH_LEN] = "";
14015cdeb108SMatthew G. Knepley 
14025cdeb108SMatthew G. Knepley   PetscFunctionBegin;
14035cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14044f572ea9SToby Isaac   if (name) PetscAssertPointer(name, 2);
140590e03537SVaclav Hapla   else name = rootGroupName;
1406ccd66a83SVaclav Hapla   if (has) {
14074f572ea9SToby Isaac     PetscAssertPointer(has, 4);
14085cdeb108SMatthew G. Knepley     *has = PETSC_FALSE;
1409ccd66a83SVaclav Hapla   }
1410ccd66a83SVaclav Hapla   if (otype) {
14114f572ea9SToby Isaac     PetscAssertPointer(otype, 5);
141256cc0592SVaclav Hapla     *otype = H5O_TYPE_UNKNOWN;
1413ccd66a83SVaclav Hapla   }
14149566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
141585688ae2SVaclav Hapla 
141685688ae2SVaclav Hapla   /*
141785688ae2SVaclav Hapla      Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing.
141885688ae2SVaclav Hapla      Hence, each of them needs to be tested separately:
141985688ae2SVaclav Hapla      1) whether it's a valid link
142085688ae2SVaclav Hapla      2) whether this link resolves to an object
142185688ae2SVaclav Hapla      See H5Oexists_by_name() documentation.
142285688ae2SVaclav Hapla   */
14239566063dSJacob Faibussowitsch   PetscCall(PetscStrToArray(name, '/', &n, &hierarchy));
142485688ae2SVaclav Hapla   if (!n) {
142585688ae2SVaclav Hapla     /*  Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */
1426ccd66a83SVaclav Hapla     if (has) *has = PETSC_TRUE;
1427ccd66a83SVaclav Hapla     if (otype) *otype = H5O_TYPE_GROUP;
14289566063dSJacob Faibussowitsch     PetscCall(PetscStrToArrayDestroy(n, hierarchy));
14293ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
143085688ae2SVaclav Hapla   }
143185688ae2SVaclav Hapla   for (i = 0; i < n; i++) {
1432c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(buf, "/", sizeof(buf)));
1433c6a7a370SJeremy L Thompson     PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf)));
14349566063dSJacob Faibussowitsch     PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists));
1435a75c3fd4SVaclav Hapla     if (!exists) break;
143685688ae2SVaclav Hapla   }
14379566063dSJacob Faibussowitsch   PetscCall(PetscStrToArrayDestroy(n, hierarchy));
143885688ae2SVaclav Hapla 
143985688ae2SVaclav Hapla   /* If the object exists, get its type */
1440ccd66a83SVaclav Hapla   if (exists && otype) {
14415cdeb108SMatthew G. Knepley     H5O_info_t info;
14425cdeb108SMatthew G. Knepley 
144376276748SVaclav Hapla     /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */
1444792fecdfSBarry Smith     PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT));
144556cc0592SVaclav Hapla     *otype = info.type;
14465cdeb108SMatthew G. Knepley   }
1447ccd66a83SVaclav Hapla   if (has) *has = exists;
14483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14495cdeb108SMatthew G. Knepley }
14505cdeb108SMatthew G. Knepley 
14514302210dSVaclav Hapla /*@C
14520a9f38efSVaclav Hapla   PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file
14530a9f38efSVaclav Hapla 
1454343a20b2SVaclav Hapla   Collective
1455343a20b2SVaclav Hapla 
14560a9f38efSVaclav Hapla   Input Parameters:
1457811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1458a0558868SVaclav Hapla - path   - (Optional) The path relative to the pushed group
14590a9f38efSVaclav Hapla 
14600a9f38efSVaclav Hapla   Output Parameter:
14610a9f38efSVaclav Hapla . has - Flag for group existence
14620a9f38efSVaclav Hapla 
14630a9f38efSVaclav Hapla   Level: advanced
14640a9f38efSVaclav Hapla 
14654302210dSVaclav Hapla   Notes:
1466785c443eSVaclav 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.
14673f423023SBarry Smith   `NULL` or empty path means the current pushed group.
14684302210dSVaclav Hapla 
1469811af0c4SBarry Smith   If path exists but is not a group, `PETSC_FALSE` is returned.
1470811af0c4SBarry Smith 
1471d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`
14720a9f38efSVaclav Hapla @*/
PetscViewerHDF5HasGroup(PetscViewer viewer,const char path[],PetscBool * has)1473d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has)
1474d71ae5a4SJacob Faibussowitsch {
14750a9f38efSVaclav Hapla   H5O_type_t  type;
1476ce78bad3SBarry Smith   const char *abspath;
14770a9f38efSVaclav Hapla 
14780a9f38efSVaclav Hapla   PetscFunctionBegin;
14790a9f38efSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
14804f572ea9SToby Isaac   if (path) PetscAssertPointer(path, 2);
14814f572ea9SToby Isaac   PetscAssertPointer(has, 3);
148277717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
14839566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
14844302210dSVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_GROUP);
14859566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
14863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
14870a9f38efSVaclav Hapla }
14880a9f38efSVaclav Hapla 
148989e0ef10SVaclav Hapla /*@C
149089e0ef10SVaclav Hapla   PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file
149189e0ef10SVaclav Hapla 
1492343a20b2SVaclav Hapla   Collective
1493343a20b2SVaclav Hapla 
149489e0ef10SVaclav Hapla   Input Parameters:
1495811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
149689e0ef10SVaclav Hapla - path   - The dataset path
149789e0ef10SVaclav Hapla 
149889e0ef10SVaclav Hapla   Output Parameter:
149989e0ef10SVaclav Hapla . has - Flag whether dataset exists
150089e0ef10SVaclav Hapla 
150189e0ef10SVaclav Hapla   Level: advanced
150289e0ef10SVaclav Hapla 
150389e0ef10SVaclav Hapla   Notes:
1504343a20b2SVaclav 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.
150589e0ef10SVaclav Hapla 
15063f423023SBarry Smith   If `path` is `NULL` or empty, has is set to `PETSC_FALSE`.
1507811af0c4SBarry Smith 
15083f423023SBarry Smith   If `path` exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1509811af0c4SBarry Smith 
1510d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
151189e0ef10SVaclav Hapla @*/
PetscViewerHDF5HasDataset(PetscViewer viewer,const char path[],PetscBool * has)1512d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has)
1513d71ae5a4SJacob Faibussowitsch {
151489e0ef10SVaclav Hapla   H5O_type_t  type;
1515ce78bad3SBarry Smith   const char *abspath;
151689e0ef10SVaclav Hapla 
151789e0ef10SVaclav Hapla   PetscFunctionBegin;
151889e0ef10SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
15194f572ea9SToby Isaac   if (path) PetscAssertPointer(path, 2);
15204f572ea9SToby Isaac   PetscAssertPointer(has, 3);
152177717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath));
15229566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type));
152389e0ef10SVaclav Hapla   *has = (PetscBool)(type == H5O_TYPE_DATASET);
15249566063dSJacob Faibussowitsch   PetscCall(PetscFree(abspath));
15253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
152689e0ef10SVaclav Hapla }
152789e0ef10SVaclav Hapla 
15280a9f38efSVaclav Hapla /*@
1529e3f143f7SVaclav Hapla   PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group
1530ecce1506SVaclav Hapla 
1531343a20b2SVaclav Hapla   Collective
1532343a20b2SVaclav Hapla 
1533ecce1506SVaclav Hapla   Input Parameters:
1534811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1535ecce1506SVaclav Hapla - obj    - The named object
1536ecce1506SVaclav Hapla 
1537ecce1506SVaclav Hapla   Output Parameter:
153889e0ef10SVaclav Hapla . has - Flag for dataset existence
1539ecce1506SVaclav Hapla 
15403f423023SBarry Smith   Level: advanced
15413f423023SBarry Smith 
1542e3f143f7SVaclav Hapla   Notes:
154389e0ef10SVaclav Hapla   If the object is unnamed, an error occurs.
1544811af0c4SBarry Smith 
1545811af0c4SBarry Smith   If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well.
1546e3f143f7SVaclav Hapla 
1547d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1548ecce1506SVaclav Hapla @*/
PetscViewerHDF5HasObject(PetscViewer viewer,PetscObject obj,PetscBool * has)1549d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has)
1550d71ae5a4SJacob Faibussowitsch {
155189e0ef10SVaclav Hapla   size_t len;
1552ecce1506SVaclav Hapla 
1553ecce1506SVaclav Hapla   PetscFunctionBegin;
1554c57a1a9aSVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1555c57a1a9aSVaclav Hapla   PetscValidHeader(obj, 2);
15564f572ea9SToby Isaac   PetscAssertPointer(has, 3);
15579566063dSJacob Faibussowitsch   PetscCall(PetscStrlen(obj->name, &len));
15585f80ce2aSJacob Faibussowitsch   PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named");
15599566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has));
15603ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1561ecce1506SVaclav Hapla }
1562ecce1506SVaclav Hapla 
1563df863907SAlex Fikl /*@C
1564ece83b6cSVaclav Hapla   PetscViewerHDF5HasAttribute - Check whether an attribute exists
1565e7bdbf8eSMatthew G. Knepley 
1566343a20b2SVaclav Hapla   Collective
1567343a20b2SVaclav Hapla 
1568e7bdbf8eSMatthew G. Knepley   Input Parameters:
1569811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
15704302210dSVaclav Hapla . parent - The parent dataset/group name
1571e7bdbf8eSMatthew G. Knepley - name   - The attribute name
1572e7bdbf8eSMatthew G. Knepley 
1573e7bdbf8eSMatthew G. Knepley   Output Parameter:
1574e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence
1575e7bdbf8eSMatthew G. Knepley 
1576e7bdbf8eSMatthew G. Knepley   Level: advanced
1577e7bdbf8eSMatthew G. Knepley 
1578811af0c4SBarry Smith   Note:
15793f423023SBarry Smith   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.
15804302210dSVaclav Hapla 
1581d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1582e7bdbf8eSMatthew G. Knepley @*/
PetscViewerHDF5HasAttribute(PetscViewer viewer,const char parent[],const char name[],PetscBool * has)1583d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1584d71ae5a4SJacob Faibussowitsch {
1585ce78bad3SBarry Smith   const char *parentAbsPath;
1586e7bdbf8eSMatthew G. Knepley 
1587e7bdbf8eSMatthew G. Knepley   PetscFunctionBegin;
15885cdeb108SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
15894f572ea9SToby Isaac   if (parent) PetscAssertPointer(parent, 2);
15904f572ea9SToby Isaac   PetscAssertPointer(name, 3);
15914f572ea9SToby Isaac   PetscAssertPointer(has, 4);
159277717648SVaclav Hapla   PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath));
15939566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL));
15949566063dSJacob Faibussowitsch   if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has));
15959566063dSJacob Faibussowitsch   PetscCall(PetscFree(parentAbsPath));
15963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
159706db490cSVaclav Hapla }
159806db490cSVaclav Hapla 
1599577e0e04SVaclav Hapla /*@C
1600811af0c4SBarry Smith   PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name
1601577e0e04SVaclav Hapla 
1602343a20b2SVaclav Hapla   Collective
1603343a20b2SVaclav Hapla 
1604577e0e04SVaclav Hapla   Input Parameters:
1605811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer
1606577e0e04SVaclav Hapla . obj    - The object whose name is used to lookup the parent dataset, relative to the current group.
1607577e0e04SVaclav Hapla - name   - The attribute name
1608577e0e04SVaclav Hapla 
1609577e0e04SVaclav Hapla   Output Parameter:
1610577e0e04SVaclav Hapla . has - Flag for attribute existence
1611577e0e04SVaclav Hapla 
16123f423023SBarry Smith   Level: advanced
16133f423023SBarry Smith 
1614811af0c4SBarry Smith   Note:
1615577e0e04SVaclav Hapla   This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset).
1616811af0c4SBarry Smith   You might want to check first if it does using `PetscViewerHDF5HasObject()`.
1617577e0e04SVaclav Hapla 
1618d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`
1619577e0e04SVaclav Hapla @*/
PetscViewerHDF5HasObjectAttribute(PetscViewer viewer,PetscObject obj,const char name[],PetscBool * has)1620d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has)
1621d71ae5a4SJacob Faibussowitsch {
1622577e0e04SVaclav Hapla   PetscFunctionBegin;
1623577e0e04SVaclav Hapla   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1);
1624577e0e04SVaclav Hapla   PetscValidHeader(obj, 2);
16254f572ea9SToby Isaac   PetscAssertPointer(name, 3);
16264f572ea9SToby Isaac   PetscAssertPointer(has, 4);
16279566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj));
16289566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has));
16293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1630577e0e04SVaclav Hapla }
1631577e0e04SVaclav Hapla 
PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer,const char parent[],const char name[],PetscBool * has)1632d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has)
1633d71ae5a4SJacob Faibussowitsch {
163406db490cSVaclav Hapla   hid_t  h5;
163506db490cSVaclav Hapla   htri_t hhas;
163606db490cSVaclav Hapla 
163706db490cSVaclav Hapla   PetscFunctionBegin;
16389566063dSJacob Faibussowitsch   PetscCall(PetscViewerHDF5GetFileId(viewer, &h5));
1639792fecdfSBarry Smith   PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT));
16405cdeb108SMatthew G. Knepley   *has = hhas ? PETSC_TRUE : PETSC_FALSE;
16413ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1642e7bdbf8eSMatthew G. Knepley }
1643e7bdbf8eSMatthew G. Knepley 
1644a75e6a4aSMatthew G. Knepley /*
1645a75e6a4aSMatthew G. Knepley   The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that
1646a75e6a4aSMatthew G. Knepley   is attached to a communicator, in this case the attribute is a PetscViewer.
1647a75e6a4aSMatthew G. Knepley */
1648d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID;
1649a75e6a4aSMatthew G. Knepley 
1650a75e6a4aSMatthew G. Knepley /*@C
1651811af0c4SBarry Smith   PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator.
1652a75e6a4aSMatthew G. Knepley 
1653d083f849SBarry Smith   Collective
1654a75e6a4aSMatthew G. Knepley 
1655a75e6a4aSMatthew G. Knepley   Input Parameter:
16563f423023SBarry Smith . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer`
1657a75e6a4aSMatthew G. Knepley 
1658811af0c4SBarry Smith   Options Database Key:
165910699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file
1660a75e6a4aSMatthew G. Knepley 
1661811af0c4SBarry Smith   Environmental variable:
1662811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file
1663a75e6a4aSMatthew G. Knepley 
166420f4b53cSBarry Smith   Level: intermediate
166520f4b53cSBarry Smith 
1666811af0c4SBarry Smith   Note:
1667811af0c4SBarry Smith   Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return
1668811af0c4SBarry Smith   an error code.  The HDF5 `PetscViewer` is usually used in the form
1669b44f4de4SBarry Smith .vb
1670b44f4de4SBarry Smith   XXXView(XXX object, PETSC_VIEWER_HDF5_(comm));
1671b44f4de4SBarry Smith .ve
1672a75e6a4aSMatthew G. Knepley 
1673d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()`
1674a75e6a4aSMatthew G. Knepley @*/
PETSC_VIEWER_HDF5_(MPI_Comm comm)1675d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm)
1676d71ae5a4SJacob Faibussowitsch {
1677a75e6a4aSMatthew G. Knepley   PetscErrorCode ierr;
16783ba16761SJacob Faibussowitsch   PetscMPIInt    mpi_ierr;
1679a75e6a4aSMatthew G. Knepley   PetscBool      flg;
1680b8b5be36SMartin Diehl   PetscMPIInt    iflg;
1681a75e6a4aSMatthew G. Knepley   PetscViewer    viewer;
1682a75e6a4aSMatthew G. Knepley   char           fname[PETSC_MAX_PATH_LEN];
1683a75e6a4aSMatthew G. Knepley   MPI_Comm       ncomm;
1684a75e6a4aSMatthew G. Knepley 
1685a75e6a4aSMatthew G. Knepley   PetscFunctionBegin;
16869371c9d4SSatish Balay   ierr = PetscCommDuplicate(comm, &ncomm, NULL);
16879371c9d4SSatish Balay   if (ierr) {
16883ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16899371c9d4SSatish Balay     PetscFunctionReturn(NULL);
16909371c9d4SSatish Balay   }
1691a75e6a4aSMatthew G. Knepley   if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) {
16923ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL);
16933ba16761SJacob Faibussowitsch     if (mpi_ierr) {
16943ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
16959371c9d4SSatish Balay       PetscFunctionReturn(NULL);
16969371c9d4SSatish Balay     }
1697a75e6a4aSMatthew G. Knepley   }
1698b8b5be36SMartin Diehl   mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, &iflg);
16993ba16761SJacob Faibussowitsch   if (mpi_ierr) {
17003ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
17019371c9d4SSatish Balay     PetscFunctionReturn(NULL);
17029371c9d4SSatish Balay   }
1703b8b5be36SMartin Diehl   if (!iflg) { /* PetscViewer not yet created */
1704a75e6a4aSMatthew G. Knepley     ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg);
17059371c9d4SSatish Balay     if (ierr) {
17063ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
17079371c9d4SSatish Balay       PetscFunctionReturn(NULL);
17089371c9d4SSatish Balay     }
1709a75e6a4aSMatthew G. Knepley     if (!flg) {
1710c6a7a370SJeremy L Thompson       ierr = PetscStrncpy(fname, "output.h5", sizeof(fname));
17119371c9d4SSatish Balay       if (ierr) {
17123ba16761SJacob Faibussowitsch         ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
17139371c9d4SSatish Balay         PetscFunctionReturn(NULL);
17149371c9d4SSatish Balay       }
1715a75e6a4aSMatthew G. Knepley     }
1716a75e6a4aSMatthew G. Knepley     ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer);
17179371c9d4SSatish Balay     if (ierr) {
17183ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
17199371c9d4SSatish Balay       PetscFunctionReturn(NULL);
17209371c9d4SSatish Balay     }
1721a75e6a4aSMatthew G. Knepley     ierr = PetscObjectRegisterDestroy((PetscObject)viewer);
17229371c9d4SSatish Balay     if (ierr) {
17233ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
17249371c9d4SSatish Balay       PetscFunctionReturn(NULL);
17259371c9d4SSatish Balay     }
17263ba16761SJacob Faibussowitsch     mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer);
17273ba16761SJacob Faibussowitsch     if (mpi_ierr) {
17283ba16761SJacob Faibussowitsch       ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " ");
17299371c9d4SSatish Balay       PetscFunctionReturn(NULL);
17309371c9d4SSatish Balay     }
1731a75e6a4aSMatthew G. Knepley   }
1732a75e6a4aSMatthew G. Knepley   ierr = PetscCommDestroy(&ncomm);
17339371c9d4SSatish Balay   if (ierr) {
17343ba16761SJacob Faibussowitsch     ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " ");
17359371c9d4SSatish Balay     PetscFunctionReturn(NULL);
17369371c9d4SSatish Balay   }
1737a75e6a4aSMatthew G. Knepley   PetscFunctionReturn(viewer);
1738a75e6a4aSMatthew G. Knepley }
1739