16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 25c6c1daeSBarry Smith 3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *); 406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *); 506db490cSVaclav Hapla 677717648SVaclav Hapla /*@C 777717648SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`. 877717648SVaclav Hapla 9*20f4b53cSBarry Smith Not Collective 1077717648SVaclav Hapla 1177717648SVaclav Hapla Input Parameters: 1277717648SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 1377717648SVaclav Hapla - path - (Optional) The path relative to the pushed group 1477717648SVaclav Hapla 1577717648SVaclav Hapla Output Parameter: 1677717648SVaclav Hapla . abspath - The absolute HDF5 path (group) 1777717648SVaclav Hapla 1877717648SVaclav Hapla Level: intermediate 1977717648SVaclav Hapla 2077717648SVaclav Hapla Notes: 2177717648SVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 223f423023SBarry Smith `NULL` or empty path means the current pushed group. 2377717648SVaclav Hapla 2477717648SVaclav Hapla The output abspath is newly allocated so needs to be freed. 2577717648SVaclav Hapla 26d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 2777717648SVaclav Hapla @*/ 28d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[]) 29d71ae5a4SJacob Faibussowitsch { 3077717648SVaclav Hapla size_t len; 314302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 326c132bc1SVaclav Hapla const char *group; 336c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 346c132bc1SVaclav Hapla 356c132bc1SVaclav Hapla PetscFunctionBegin; 3677717648SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 3777717648SVaclav Hapla if (path) PetscValidCharPointer(path, 2); 3877717648SVaclav Hapla PetscValidPointer(abspath, 3); 3977717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group)); 4077717648SVaclav Hapla PetscCall(PetscStrlen(path, &len)); 4177717648SVaclav Hapla relative = (PetscBool)(!len || path[0] != '/'); 424302210dSVaclav Hapla if (relative) { 43c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(buf, group, sizeof(buf))); 44c6a7a370SJeremy L Thompson if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf))); 45c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, path, sizeof(buf))); 469566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(buf, abspath)); 474302210dSVaclav Hapla } else { 489566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(path, abspath)); 494302210dSVaclav Hapla } 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516c132bc1SVaclav Hapla } 526c132bc1SVaclav Hapla 53d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 54d71ae5a4SJacob Faibussowitsch { 55577e0e04SVaclav Hapla PetscBool has; 56577e0e04SVaclav Hapla 57577e0e04SVaclav Hapla PetscFunctionBegin; 589566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has)); 59577e0e04SVaclav Hapla if (!has) { 6077717648SVaclav Hapla char *group; 6177717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 6277717648SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group); 63577e0e04SVaclav Hapla } 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65577e0e04SVaclav Hapla } 66577e0e04SVaclav Hapla 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject) 68d71ae5a4SJacob Faibussowitsch { 69ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 7082ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 7182ea9b62SBarry Smith 7282ea9b62SBarry Smith PetscFunctionBegin; 73d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options"); 749566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL)); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL)); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set)); 779566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg)); 7819a20e4cSMatthew G. Knepley flg = PETSC_FALSE; 799566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set)); 809566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg)); 81d0609cedSBarry Smith PetscOptionsHeadEnd(); 823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8382ea9b62SBarry Smith } 8482ea9b62SBarry Smith 85d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer) 86d71ae5a4SJacob Faibussowitsch { 871b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 8803fa8834SVaclav Hapla PetscBool flg; 891b793a25SVaclav Hapla 901b793a25SVaclav Hapla PetscFunctionBegin; 9148a46eb9SPierre Jolivet if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename)); 929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2])); 939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput])); 949566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetCollective(v, &flg)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent")); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping])); 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 981b793a25SVaclav Hapla } 991b793a25SVaclav Hapla 100d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 101d71ae5a4SJacob Faibussowitsch { 1025c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1035c6c1daeSBarry Smith 1045c6c1daeSBarry Smith PetscFunctionBegin; 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->filename)); 106792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1085c6c1daeSBarry Smith } 1095c6c1daeSBarry Smith 110d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) 111d71ae5a4SJacob Faibussowitsch { 1126226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1136226335aSVaclav Hapla 1146226335aSVaclav Hapla PetscFunctionBegin; 115792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL)); 1163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1176226335aSVaclav Hapla } 1186226335aSVaclav Hapla 119d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 120d71ae5a4SJacob Faibussowitsch { 1215c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1225c6c1daeSBarry Smith 1235c6c1daeSBarry Smith PetscFunctionBegin; 124792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (hdf5->dxpl_id)); 1259566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_HDF5(viewer)); 1265c6c1daeSBarry Smith while (hdf5->groups) { 1277d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1285c6c1daeSBarry Smith 1299566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups->name)); 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups)); 1315c6c1daeSBarry Smith hdf5->groups = tmp; 1325c6c1daeSBarry Smith } 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5)); 1349566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 1372e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 1389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL)); 1402e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL)); 1412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL)); 1443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1455c6c1daeSBarry Smith } 1465c6c1daeSBarry Smith 147d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 148d71ae5a4SJacob Faibussowitsch { 1495c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1505c6c1daeSBarry Smith 1515c6c1daeSBarry Smith PetscFunctionBegin; 1525c6c1daeSBarry Smith hdf5->btype = type; 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1545c6c1daeSBarry Smith } 1555c6c1daeSBarry Smith 156d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 157d71ae5a4SJacob Faibussowitsch { 1580b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1590b62783dSVaclav Hapla 1600b62783dSVaclav Hapla PetscFunctionBegin; 1610b62783dSVaclav Hapla *type = hdf5->btype; 1623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1630b62783dSVaclav Hapla } 1640b62783dSVaclav Hapla 165d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 166d71ae5a4SJacob Faibussowitsch { 16782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 16882ea9b62SBarry Smith 16982ea9b62SBarry Smith PetscFunctionBegin; 17082ea9b62SBarry Smith hdf5->basedimension2 = flg; 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17282ea9b62SBarry Smith } 17382ea9b62SBarry Smith 174df863907SAlex Fikl /*@ 17582ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17682ea9b62SBarry Smith dimension of 2. 17782ea9b62SBarry Smith 178c3339decSBarry Smith Logically Collective 17982ea9b62SBarry Smith 18082ea9b62SBarry Smith Input Parameters: 181811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored 182811af0c4SBarry Smith - flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1 18382ea9b62SBarry Smith 184811af0c4SBarry Smith Options Database Key: 18582ea9b62SBarry Smith . -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1 18682ea9b62SBarry Smith 1873f423023SBarry Smith Level: intermediate 1883f423023SBarry Smith 189811af0c4SBarry Smith Note: 19095452b02SPatrick 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 19182ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19282ea9b62SBarry Smith 193d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 19482ea9b62SBarry Smith @*/ 195d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg) 196d71ae5a4SJacob Faibussowitsch { 19782ea9b62SBarry Smith PetscFunctionBegin; 19882ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 199cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20182ea9b62SBarry Smith } 20282ea9b62SBarry Smith 203df863907SAlex Fikl /*@ 20482ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 20582ea9b62SBarry Smith dimension of 2. 20682ea9b62SBarry Smith 207c3339decSBarry Smith Logically Collective 20882ea9b62SBarry Smith 20982ea9b62SBarry Smith Input Parameter: 210811af0c4SBarry Smith . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5` 21182ea9b62SBarry Smith 21282ea9b62SBarry Smith Output Parameter: 213811af0c4SBarry Smith . flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1 21482ea9b62SBarry Smith 2153f423023SBarry Smith Level: intermediate 2163f423023SBarry Smith 217811af0c4SBarry Smith Note: 21895452b02SPatrick 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 21982ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 22082ea9b62SBarry Smith 221d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 22282ea9b62SBarry Smith @*/ 223d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg) 224d71ae5a4SJacob Faibussowitsch { 22582ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 22682ea9b62SBarry Smith 22782ea9b62SBarry Smith PetscFunctionBegin; 22882ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 22982ea9b62SBarry Smith *flg = hdf5->basedimension2; 2303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23182ea9b62SBarry Smith } 23282ea9b62SBarry Smith 233d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 234d71ae5a4SJacob Faibussowitsch { 2359a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2369a0502c6SHåkon Strandenes 2379a0502c6SHåkon Strandenes PetscFunctionBegin; 2389a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2409a0502c6SHåkon Strandenes } 2419a0502c6SHåkon Strandenes 242df863907SAlex Fikl /*@ 2439a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 244811af0c4SBarry Smith compiled with double precision `PetscReal`. 2459a0502c6SHåkon Strandenes 246c3339decSBarry Smith Logically Collective 2479a0502c6SHåkon Strandenes 2489a0502c6SHåkon Strandenes Input Parameters: 249811af0c4SBarry Smith + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored 250811af0c4SBarry Smith - flg - if `PETSC_TRUE` the data will be written to disk with single precision 2519a0502c6SHåkon Strandenes 252811af0c4SBarry Smith Options Database Key: 2539a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2549a0502c6SHåkon Strandenes 2553f423023SBarry Smith Level: intermediate 2563f423023SBarry Smith 257811af0c4SBarry Smith Note: 25895452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2599a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2609a0502c6SHåkon Strandenes 261d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 262811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5GetSPOutput()` 2639a0502c6SHåkon Strandenes @*/ 264d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg) 265d71ae5a4SJacob Faibussowitsch { 2669a0502c6SHåkon Strandenes PetscFunctionBegin; 2679a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 268cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg)); 2693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2709a0502c6SHåkon Strandenes } 2719a0502c6SHåkon Strandenes 272df863907SAlex Fikl /*@ 2739a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 274811af0c4SBarry Smith compiled with double precision `PetscReal`. 2759a0502c6SHåkon Strandenes 276c3339decSBarry Smith Logically Collective 2779a0502c6SHåkon Strandenes 2789a0502c6SHåkon Strandenes Input Parameter: 279811af0c4SBarry Smith . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5` 2809a0502c6SHåkon Strandenes 2819a0502c6SHåkon Strandenes Output Parameter: 282811af0c4SBarry Smith . flg - if `PETSC_TRUE` the data will be written to disk with single precision 2839a0502c6SHåkon Strandenes 2849a0502c6SHåkon Strandenes Level: intermediate 2859a0502c6SHåkon Strandenes 286d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 287811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5SetSPOutput()` 2889a0502c6SHåkon Strandenes @*/ 289d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg) 290d71ae5a4SJacob Faibussowitsch { 2919a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2929a0502c6SHåkon Strandenes 2939a0502c6SHåkon Strandenes PetscFunctionBegin; 2949a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 2959a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2979a0502c6SHåkon Strandenes } 2989a0502c6SHåkon Strandenes 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 300d71ae5a4SJacob Faibussowitsch { 301ee8b9732SVaclav Hapla PetscFunctionBegin; 302ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 303ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 304c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL) 305ee8b9732SVaclav Hapla { 306ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 307792fecdfSBarry Smith PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 308ee8b9732SVaclav Hapla } 309ee8b9732SVaclav Hapla #else 3109566063dSJacob Faibussowitsch if (flg) PetscCall(PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n")); 311ee8b9732SVaclav Hapla #endif 3123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 313ee8b9732SVaclav Hapla } 314ee8b9732SVaclav Hapla 315ee8b9732SVaclav Hapla /*@ 316ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 317ee8b9732SVaclav Hapla 318ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 319ee8b9732SVaclav Hapla 320ee8b9732SVaclav Hapla Input Parameters: 321811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 322811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default) 323ee8b9732SVaclav Hapla 324811af0c4SBarry Smith Options Database Key: 325ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 326ee8b9732SVaclav Hapla 3273f423023SBarry Smith Level: intermediate 3283f423023SBarry Smith 329811af0c4SBarry Smith Note: 330ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 33153effed7SBarry 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. 332ee8b9732SVaclav Hapla 333811af0c4SBarry Smith Developer Note: 334811af0c4SBarry Smith In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively. 335ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 336ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 337ee8b9732SVaclav Hapla 338d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 339ee8b9732SVaclav Hapla @*/ 340d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg) 341d71ae5a4SJacob Faibussowitsch { 342ee8b9732SVaclav Hapla PetscFunctionBegin; 343ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 344ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer, flg, 2); 345cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg)); 3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 347ee8b9732SVaclav Hapla } 348ee8b9732SVaclav Hapla 349d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 350d71ae5a4SJacob Faibussowitsch { 351c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 352ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 353ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 354c796909eSBarry Smith #endif 355ee8b9732SVaclav Hapla 356ee8b9732SVaclav Hapla PetscFunctionBegin; 357c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 358c796909eSBarry Smith *flg = PETSC_FALSE; 359c796909eSBarry Smith #else 360792fecdfSBarry Smith PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode)); 361ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 362c796909eSBarry Smith #endif 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 364ee8b9732SVaclav Hapla } 365ee8b9732SVaclav Hapla 366ee8b9732SVaclav Hapla /*@ 367ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 368ee8b9732SVaclav Hapla 369ee8b9732SVaclav Hapla Not Collective 370ee8b9732SVaclav Hapla 371*20f4b53cSBarry Smith Input Parameter: 372811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer` 373ee8b9732SVaclav Hapla 374*20f4b53cSBarry Smith Output Parameter: 375ee8b9732SVaclav Hapla . flg - the flag 376ee8b9732SVaclav Hapla 377ee8b9732SVaclav Hapla Level: intermediate 378ee8b9732SVaclav Hapla 379811af0c4SBarry Smith Note: 380811af0c4SBarry Smith This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, `PETSC_FALSE` will be always returned. 381811af0c4SBarry Smith For more details, see `PetscViewerHDF5SetCollective()`. 382ee8b9732SVaclav Hapla 383d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 384ee8b9732SVaclav Hapla @*/ 385d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg) 386d71ae5a4SJacob Faibussowitsch { 387ee8b9732SVaclav Hapla PetscFunctionBegin; 388ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 389534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 390ee8b9732SVaclav Hapla 391cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg)); 3923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 393ee8b9732SVaclav Hapla } 394ee8b9732SVaclav Hapla 395d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 396d71ae5a4SJacob Faibussowitsch { 3975c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 3985c6c1daeSBarry Smith hid_t plist_id; 3995c6c1daeSBarry Smith 4005c6c1daeSBarry Smith PetscFunctionBegin; 401792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 4029566063dSJacob Faibussowitsch if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 4039566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &hdf5->filename)); 4045c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 405792fecdfSBarry Smith PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS)); 406c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 407792fecdfSBarry Smith PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 408c796909eSBarry Smith #endif 4095c6c1daeSBarry Smith /* Create or open the file collectively */ 4105c6c1daeSBarry Smith switch (hdf5->btype) { 4115c6c1daeSBarry Smith case FILE_MODE_READ: 4128a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 4138a2871f6SBarry Smith PetscMPIInt rank; 4148a2871f6SBarry Smith PetscBool flg; 4158a2871f6SBarry Smith 4169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 4178a2871f6SBarry Smith if (rank == 0) { 4189566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4195f80ce2aSJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename); 4208a2871f6SBarry Smith } 4219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 4228a2871f6SBarry Smith } 423792fecdfSBarry Smith PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id)); 4245c6c1daeSBarry Smith break; 4255c6c1daeSBarry Smith case FILE_MODE_APPEND: 4269371c9d4SSatish Balay case FILE_MODE_UPDATE: { 4277e4fd573SVaclav Hapla PetscBool flg; 4289566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 429792fecdfSBarry Smith if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id)); 430792fecdfSBarry Smith else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4315c6c1daeSBarry Smith break; 4327e4fd573SVaclav Hapla } 433d71ae5a4SJacob Faibussowitsch case FILE_MODE_WRITE: 434d71ae5a4SJacob Faibussowitsch PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 435d71ae5a4SJacob Faibussowitsch break; 436d71ae5a4SJacob Faibussowitsch case FILE_MODE_UNDEFINED: 437d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 438d71ae5a4SJacob Faibussowitsch default: 439d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]); 4405c6c1daeSBarry Smith } 4415f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 442792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (plist_id)); 443671abe24SVaclav Hapla PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer)); 4443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4455c6c1daeSBarry Smith } 4465c6c1daeSBarry Smith 447d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name) 448d71ae5a4SJacob Faibussowitsch { 449d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data; 450d1232d7fSToby Isaac 451d1232d7fSToby Isaac PetscFunctionBegin; 452d1232d7fSToby Isaac *name = vhdf5->filename; 4533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 454d1232d7fSToby Isaac } 455d1232d7fSToby Isaac 456d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 457d71ae5a4SJacob Faibussowitsch { 4585dc64a97SVaclav Hapla /* 459b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4605dc64a97SVaclav Hapla */ 461b723ab35SVaclav Hapla 462b723ab35SVaclav Hapla PetscFunctionBegin; 4633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 464b723ab35SVaclav Hapla } 465b723ab35SVaclav Hapla 466d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 467d71ae5a4SJacob Faibussowitsch { 46819a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 46919a20e4cSMatthew G. Knepley 47019a20e4cSMatthew G. Knepley PetscFunctionBegin; 47119a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 4723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47319a20e4cSMatthew G. Knepley } 47419a20e4cSMatthew G. Knepley 475d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 476d71ae5a4SJacob Faibussowitsch { 47719a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 47819a20e4cSMatthew G. Knepley 47919a20e4cSMatthew G. Knepley PetscFunctionBegin; 48019a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 4813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48219a20e4cSMatthew G. Knepley } 48319a20e4cSMatthew G. Knepley 48419a20e4cSMatthew G. Knepley /*@ 48519a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 48619a20e4cSMatthew G. Knepley 487c3339decSBarry Smith Logically Collective 48819a20e4cSMatthew G. Knepley 48919a20e4cSMatthew G. Knepley Input Parameters: 490811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 491811af0c4SBarry Smith - flg - if `PETSC_TRUE` we will assume that timestepping is on 49219a20e4cSMatthew G. Knepley 493811af0c4SBarry Smith Options Database Key: 49419a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 49519a20e4cSMatthew G. Knepley 4963f423023SBarry Smith Level: intermediate 4973f423023SBarry Smith 498811af0c4SBarry Smith Note: 49919a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 50019a20e4cSMatthew G. Knepley 501d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 50219a20e4cSMatthew G. Knepley @*/ 503d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 504d71ae5a4SJacob Faibussowitsch { 50519a20e4cSMatthew G. Knepley PetscFunctionBegin; 50619a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 507cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg)); 5083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 50919a20e4cSMatthew G. Knepley } 51019a20e4cSMatthew G. Knepley 51119a20e4cSMatthew G. Knepley /*@ 51219a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 51319a20e4cSMatthew G. Knepley 514*20f4b53cSBarry Smith Not Collective 51519a20e4cSMatthew G. Knepley 51619a20e4cSMatthew G. Knepley Input Parameter: 517811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 51819a20e4cSMatthew G. Knepley 51919a20e4cSMatthew G. Knepley Output Parameter: 520811af0c4SBarry Smith . flg - if `PETSC_TRUE` we will assume that timestepping is on 52119a20e4cSMatthew G. Knepley 52219a20e4cSMatthew G. Knepley Level: intermediate 52319a20e4cSMatthew G. Knepley 524d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 52519a20e4cSMatthew G. Knepley @*/ 526d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 527d71ae5a4SJacob Faibussowitsch { 52819a20e4cSMatthew G. Knepley PetscFunctionBegin; 52919a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 530cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg)); 5313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53219a20e4cSMatthew G. Knepley } 53319a20e4cSMatthew G. Knepley 5348556b5ebSBarry Smith /*MC 5358556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5368556b5ebSBarry Smith 537811af0c4SBarry Smith Level: beginner 538811af0c4SBarry Smith 539d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 540db781477SPatrick Sanan `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 541db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 542db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 5438556b5ebSBarry Smith M*/ 544d1232d7fSToby Isaac 545d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 546d71ae5a4SJacob Faibussowitsch { 5475c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5485c6c1daeSBarry Smith 5495c6c1daeSBarry Smith PetscFunctionBegin; 55099335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 55199335c34SVaclav Hapla { 55299335c34SVaclav Hapla PetscMPIInt size; 5539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 5545f80ce2aSJacob Faibussowitsch PetscCheck(size <= 1, PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)"); 55599335c34SVaclav Hapla } 55699335c34SVaclav Hapla #endif 55799335c34SVaclav Hapla 5584dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&hdf5)); 5595c6c1daeSBarry Smith 5605c6c1daeSBarry Smith v->data = (void *)hdf5; 5615c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 56282ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 563b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5641b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5656226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5667e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 567908793a3SLisandro Dalcin hdf5->filename = NULL; 5685c6c1daeSBarry Smith hdf5->timestep = -1; 5690298fd71SBarry Smith hdf5->groups = NULL; 5705c6c1daeSBarry Smith 571792fecdfSBarry Smith PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER)); 5729c5072fbSVaclav Hapla 5739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5)); 5749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5)); 5759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5)); 5769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5)); 5779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5)); 5789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5)); 5799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5)); 5809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5)); 5819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5)); 5829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5)); 5833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5845c6c1daeSBarry Smith } 5855c6c1daeSBarry Smith 5865c6c1daeSBarry Smith /*@C 587811af0c4SBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer` 5885c6c1daeSBarry Smith 589d083f849SBarry Smith Collective 5905c6c1daeSBarry Smith 5915c6c1daeSBarry Smith Input Parameters: 5925c6c1daeSBarry Smith + comm - MPI communicator 5935c6c1daeSBarry Smith . name - name of file 5945c6c1daeSBarry Smith - type - type of file 5955c6c1daeSBarry Smith 5965c6c1daeSBarry Smith Output Parameter: 597811af0c4SBarry Smith . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file 5985c6c1daeSBarry Smith 599811af0c4SBarry Smith Options Database Keys: 600a2b725a8SWilliam Gropp + -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1 601a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 60282ea9b62SBarry Smith 6035c6c1daeSBarry Smith Level: beginner 6045c6c1daeSBarry Smith 6057e4fd573SVaclav Hapla Notes: 6067e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 607811af0c4SBarry Smith .vb 608811af0c4SBarry Smith FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 609811af0c4SBarry Smith FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 610811af0c4SBarry Smith FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL] 611811af0c4SBarry Smith FILE_MODE_UPDATE - same as FILE_MODE_APPEND 612811af0c4SBarry Smith .ve 6137e4fd573SVaclav Hapla 614da81f932SPierre 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. 6157e4fd573SVaclav Hapla 6165c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 6175c6c1daeSBarry Smith 618d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 619db781477SPatrick Sanan `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 620db781477SPatrick Sanan `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 6215c6c1daeSBarry Smith @*/ 622d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 623d71ae5a4SJacob Faibussowitsch { 6245c6c1daeSBarry Smith PetscFunctionBegin; 6259566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, hdf5v)); 6269566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 6279566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 6289566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*hdf5v, name)); 6299566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*hdf5v)); 6303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6315c6c1daeSBarry Smith } 6325c6c1daeSBarry Smith 6335c6c1daeSBarry Smith /*@C 6345c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6355c6c1daeSBarry Smith 636*20f4b53cSBarry Smith Not Collective 6375c6c1daeSBarry Smith 6385c6c1daeSBarry Smith Input Parameter: 639811af0c4SBarry Smith . viewer - the `PetscViewer` 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith Output Parameter: 6425c6c1daeSBarry Smith . file_id - The file id 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith Level: intermediate 6455c6c1daeSBarry Smith 646d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()` 6475c6c1daeSBarry Smith @*/ 648d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 649d71ae5a4SJacob Faibussowitsch { 6505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6515c6c1daeSBarry Smith 6525c6c1daeSBarry Smith PetscFunctionBegin; 6535c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 6545c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6565c6c1daeSBarry Smith } 6575c6c1daeSBarry Smith 6585c6c1daeSBarry Smith /*@C 6595c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6605c6c1daeSBarry Smith 661*20f4b53cSBarry Smith Not Collective 6625c6c1daeSBarry Smith 6635c6c1daeSBarry Smith Input Parameters: 664811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 6655c6c1daeSBarry Smith - name - The group name 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith Level: intermediate 6685c6c1daeSBarry Smith 6692e361344SVaclav Hapla Notes: 6702e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 671811af0c4SBarry Smith .vb 672811af0c4SBarry Smith If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group. 6733f423023SBarry Smith `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 674811af0c4SBarry Smith "." means the current group is pushed again. 675811af0c4SBarry Smith .ve 6762e361344SVaclav Hapla 6772e361344SVaclav Hapla Example: 6782e361344SVaclav Hapla Suppose the current group is "/a". 6793f423023SBarry Smith .vb 6803f423023SBarry Smith If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6813f423023SBarry Smith If name is ".", then the new group will be "/a". 6823f423023SBarry Smith If name is "b", then the new group will be "/a/b". 6833f423023SBarry Smith If name is "/b", then the new group will be "/b". 6843f423023SBarry Smith .ve 6852e361344SVaclav Hapla 686811af0c4SBarry Smith Developer Note: 6873f423023SBarry Smith The root group "/" is internally stored as `NULL`. 688820fc2d1SVaclav Hapla 689d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 6905c6c1daeSBarry Smith @*/ 691d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 692d71ae5a4SJacob Faibussowitsch { 6935c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6947d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6952e361344SVaclav Hapla size_t i, len; 6962e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6972e361344SVaclav Hapla const char *gname; 6985c6c1daeSBarry Smith 6995c6c1daeSBarry Smith PetscFunctionBegin; 7005c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 701820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 7029566063dSJacob Faibussowitsch PetscCall(PetscStrlen(name, &len)); 7032e361344SVaclav Hapla gname = NULL; 7042e361344SVaclav Hapla if (len) { 7052e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 7062e361344SVaclav Hapla /* use current name */ 7072e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 7082e361344SVaclav Hapla } else if (name[0] == '/') { 7092e361344SVaclav Hapla /* absolute */ 7102e361344SVaclav Hapla for (i = 1; i < len; i++) { 7112e361344SVaclav Hapla if (name[i] != '/') { 7122e361344SVaclav Hapla gname = name; 7132e361344SVaclav Hapla break; 7142e361344SVaclav Hapla } 7152e361344SVaclav Hapla } 7162e361344SVaclav Hapla } else { 7172e361344SVaclav Hapla /* relative */ 7182e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 7199566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 7202e361344SVaclav Hapla gname = buf; 7212e361344SVaclav Hapla } 7222e361344SVaclav Hapla } 7239566063dSJacob Faibussowitsch PetscCall(PetscNew(&groupNode)); 7249566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name)); 7255c6c1daeSBarry Smith groupNode->next = hdf5->groups; 7265c6c1daeSBarry Smith hdf5->groups = groupNode; 7273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7285c6c1daeSBarry Smith } 7295c6c1daeSBarry Smith 7303ef9c667SSatish Balay /*@ 7315c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 7325c6c1daeSBarry Smith 733*20f4b53cSBarry Smith Not Collective 7345c6c1daeSBarry Smith 7355c6c1daeSBarry Smith Input Parameter: 736811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7375c6c1daeSBarry Smith 7385c6c1daeSBarry Smith Level: intermediate 7395c6c1daeSBarry Smith 740d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 7415c6c1daeSBarry Smith @*/ 742d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 743d71ae5a4SJacob Faibussowitsch { 7445c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7457d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7465c6c1daeSBarry Smith 7475c6c1daeSBarry Smith PetscFunctionBegin; 7485c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 7495f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7505c6c1daeSBarry Smith groupNode = hdf5->groups; 7515c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7529566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode->name)); 7539566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode)); 7543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7555c6c1daeSBarry Smith } 7565c6c1daeSBarry Smith 757d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[]) 758d71ae5a4SJacob Faibussowitsch { 7595c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7605c6c1daeSBarry Smith 7615c6c1daeSBarry Smith PetscFunctionBegin; 7625c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 763c959eef4SJed Brown PetscValidPointer(name, 2); 764a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 7650298fd71SBarry Smith else *name = NULL; 7663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7675c6c1daeSBarry Smith } 7685c6c1daeSBarry Smith 7693014b61aSVaclav Hapla /*@C 770811af0c4SBarry Smith PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`, 771874270d9SVaclav Hapla and return this group's ID and file ID. 772811af0c4SBarry Smith If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID. 773874270d9SVaclav Hapla 774*20f4b53cSBarry Smith Not Collective 775874270d9SVaclav Hapla 7763014b61aSVaclav Hapla Input Parameters: 7773014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7783014b61aSVaclav Hapla - path - (Optional) The path relative to the pushed group 779874270d9SVaclav Hapla 780d8d19677SJose E. Roman Output Parameters: 781874270d9SVaclav Hapla + fileId - The HDF5 file ID 782874270d9SVaclav Hapla - groupId - The HDF5 group ID 783874270d9SVaclav Hapla 7843f423023SBarry Smith Level: intermediate 7853f423023SBarry Smith 786811af0c4SBarry Smith Note: 7873014b61aSVaclav 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. 7883f423023SBarry Smith `NULL` or empty path means the current pushed group. 7893014b61aSVaclav Hapla 790e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 791e74616beSVaclav Hapla 792d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()` 793874270d9SVaclav Hapla @*/ 794d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId) 795d71ae5a4SJacob Faibussowitsch { 79690e03537SVaclav Hapla hid_t file_id; 79790e03537SVaclav Hapla H5O_type_t type; 7983014b61aSVaclav Hapla const char *fileName = NULL; 7993014b61aSVaclav Hapla char *groupName = NULL; 80076d59af2SVaclav Hapla PetscBool writable, has; 80154dbf706SVaclav Hapla 80254dbf706SVaclav Hapla PetscFunctionBegin; 8039566063dSJacob Faibussowitsch PetscCall(PetscViewerWritable(viewer, &writable)); 8049566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 8059566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetName(viewer, &fileName)); 8063014b61aSVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName)); 8079566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 80876d59af2SVaclav Hapla if (!has) { 8095f80ce2aSJacob Faibussowitsch PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName); 810f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 81176d59af2SVaclav Hapla } 8125f80ce2aSJacob Faibussowitsch PetscCheck(type == H5O_TYPE_GROUP, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName); 8133014b61aSVaclav Hapla PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT)); 8143014b61aSVaclav Hapla PetscCall(PetscFree(groupName)); 81554dbf706SVaclav Hapla *fileId = file_id; 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81754dbf706SVaclav Hapla } 81854dbf706SVaclav Hapla 81963cb69f5SVaclav Hapla /*@C 82063cb69f5SVaclav Hapla PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file 82163cb69f5SVaclav Hapla 822*20f4b53cSBarry Smith Not Collective 82363cb69f5SVaclav Hapla 82463cb69f5SVaclav Hapla Input Parameters: 82563cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 82663cb69f5SVaclav Hapla - path - (Optional) The path relative to the pushed group 82763cb69f5SVaclav Hapla 8283f423023SBarry Smith Level: intermediate 8293f423023SBarry Smith 83063cb69f5SVaclav Hapla Note: 83163cb69f5SVaclav 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. 8323f423023SBarry Smith `NULL` or empty path means the current pushed group. 83363cb69f5SVaclav Hapla 83463cb69f5SVaclav Hapla This will fail if the viewer is not writable. 83563cb69f5SVaclav Hapla 836d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 83763cb69f5SVaclav Hapla @*/ 838d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[]) 839d71ae5a4SJacob Faibussowitsch { 84063cb69f5SVaclav Hapla hid_t fileId, groupId; 84163cb69f5SVaclav Hapla 84263cb69f5SVaclav Hapla PetscFunctionBegin; 84363cb69f5SVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created 84463cb69f5SVaclav Hapla PetscCallHDF5(H5Gclose, (groupId)); 8453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84663cb69f5SVaclav Hapla } 84763cb69f5SVaclav Hapla 8483ef9c667SSatish Balay /*@ 849d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8505c6c1daeSBarry Smith 851*20f4b53cSBarry Smith Not Collective 8525c6c1daeSBarry Smith 8535c6c1daeSBarry Smith Input Parameter: 854811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 8555c6c1daeSBarry Smith 8565c6c1daeSBarry Smith Level: intermediate 8575c6c1daeSBarry Smith 858d7dd068bSVaclav Hapla Notes: 859811af0c4SBarry Smith On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0. 860811af0c4SBarry Smith Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`. 861811af0c4SBarry 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()`]. 862811af0c4SBarry Smith Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 863811af0c4SBarry Smith Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`. 864d7dd068bSVaclav Hapla 865d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 866d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 867d7dd068bSVaclav Hapla 868811af0c4SBarry Smith Developer note: 869d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 870d7dd068bSVaclav Hapla 871d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 872d7dd068bSVaclav Hapla @*/ 873d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 874d71ae5a4SJacob Faibussowitsch { 875d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 876d7dd068bSVaclav Hapla 877d7dd068bSVaclav Hapla PetscFunctionBegin; 878d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 8795f80ce2aSJacob Faibussowitsch PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 880d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 881d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 8823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 883d7dd068bSVaclav Hapla } 884d7dd068bSVaclav Hapla 885d7dd068bSVaclav Hapla /*@ 886d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 887d7dd068bSVaclav Hapla 888*20f4b53cSBarry Smith Not Collective 889d7dd068bSVaclav Hapla 890d7dd068bSVaclav Hapla Input Parameter: 891811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 892d7dd068bSVaclav Hapla 893d7dd068bSVaclav Hapla Level: intermediate 894d7dd068bSVaclav Hapla 895811af0c4SBarry Smith Note: 896811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 897d7dd068bSVaclav Hapla 898d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 899d7dd068bSVaclav Hapla @*/ 900d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 901d71ae5a4SJacob Faibussowitsch { 902d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 903d7dd068bSVaclav Hapla 904d7dd068bSVaclav Hapla PetscFunctionBegin; 905d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 9065f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 907d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 9083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 909d7dd068bSVaclav Hapla } 910d7dd068bSVaclav Hapla 911d7dd068bSVaclav Hapla /*@ 912d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 913d7dd068bSVaclav Hapla 914*20f4b53cSBarry Smith Not Collective 915d7dd068bSVaclav Hapla 916d7dd068bSVaclav Hapla Input Parameter: 917811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 918d7dd068bSVaclav Hapla 919d7dd068bSVaclav Hapla Output Parameter: 920d7dd068bSVaclav Hapla . flg - is timestepping active? 921d7dd068bSVaclav Hapla 922d7dd068bSVaclav Hapla Level: intermediate 923d7dd068bSVaclav Hapla 924811af0c4SBarry Smith Note: 925811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 926d7dd068bSVaclav Hapla 927d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 928d7dd068bSVaclav Hapla @*/ 929d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 930d71ae5a4SJacob Faibussowitsch { 931d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 932d7dd068bSVaclav Hapla 933d7dd068bSVaclav Hapla PetscFunctionBegin; 934d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 935d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 9363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 937d7dd068bSVaclav Hapla } 938d7dd068bSVaclav Hapla 939d7dd068bSVaclav Hapla /*@ 940d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 941d7dd068bSVaclav Hapla 942*20f4b53cSBarry Smith Not Collective 943d7dd068bSVaclav Hapla 944d7dd068bSVaclav Hapla Input Parameter: 945811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 946d7dd068bSVaclav Hapla 947d7dd068bSVaclav Hapla Level: intermediate 948d7dd068bSVaclav Hapla 949811af0c4SBarry Smith Note: 950811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 951d7dd068bSVaclav Hapla 952d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 9535c6c1daeSBarry Smith @*/ 954d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 955d71ae5a4SJacob Faibussowitsch { 9565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9575c6c1daeSBarry Smith 9585c6c1daeSBarry Smith PetscFunctionBegin; 9595c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 9605f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9615c6c1daeSBarry Smith ++hdf5->timestep; 9623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9635c6c1daeSBarry Smith } 9645c6c1daeSBarry Smith 9653ef9c667SSatish Balay /*@ 966d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9675c6c1daeSBarry Smith 968*20f4b53cSBarry Smith Logically Collective 9695c6c1daeSBarry Smith 9705c6c1daeSBarry Smith Input Parameters: 971811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 972d7dd068bSVaclav Hapla - timestep - The timestep 9735c6c1daeSBarry Smith 9745c6c1daeSBarry Smith Level: intermediate 9755c6c1daeSBarry Smith 976811af0c4SBarry Smith Note: 977811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 978d7dd068bSVaclav Hapla 979d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 9805c6c1daeSBarry Smith @*/ 981d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 982d71ae5a4SJacob Faibussowitsch { 9835c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9845c6c1daeSBarry Smith 9855c6c1daeSBarry Smith PetscFunctionBegin; 9865c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 987d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 9885f80ce2aSJacob Faibussowitsch PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 9895f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9905c6c1daeSBarry Smith hdf5->timestep = timestep; 9913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9925c6c1daeSBarry Smith } 9935c6c1daeSBarry Smith 9943ef9c667SSatish Balay /*@ 9955c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 9965c6c1daeSBarry Smith 997*20f4b53cSBarry Smith Not Collective 9985c6c1daeSBarry Smith 9995c6c1daeSBarry Smith Input Parameter: 1000811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 10015c6c1daeSBarry Smith 10025c6c1daeSBarry Smith Output Parameter: 1003d7dd068bSVaclav Hapla . timestep - The timestep 10045c6c1daeSBarry Smith 10055c6c1daeSBarry Smith Level: intermediate 10065c6c1daeSBarry Smith 1007811af0c4SBarry Smith Note: 1008811af0c4SBarry Smith This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 1009d7dd068bSVaclav Hapla 1010d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 10115c6c1daeSBarry Smith @*/ 1012d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 1013d71ae5a4SJacob Faibussowitsch { 10145c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 10155c6c1daeSBarry Smith 10165c6c1daeSBarry Smith PetscFunctionBegin; 10175c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1018d7dd068bSVaclav Hapla PetscValidIntPointer(timestep, 2); 10195f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 10205c6c1daeSBarry Smith *timestep = hdf5->timestep; 10213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10225c6c1daeSBarry Smith } 10235c6c1daeSBarry Smith 102436bce6e8SMatthew G. Knepley /*@C 102536bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 102636bce6e8SMatthew G. Knepley 1027*20f4b53cSBarry Smith Not Collective 102836bce6e8SMatthew G. Knepley 102936bce6e8SMatthew G. Knepley Input Parameter: 1030811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 103136bce6e8SMatthew G. Knepley 103236bce6e8SMatthew G. Knepley Output Parameter: 1033811af0c4SBarry Smith . mtype - the HDF5 datatype 103436bce6e8SMatthew G. Knepley 103536bce6e8SMatthew G. Knepley Level: advanced 103636bce6e8SMatthew G. Knepley 1037d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 103836bce6e8SMatthew G. Knepley @*/ 1039d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 1040d71ae5a4SJacob Faibussowitsch { 104136bce6e8SMatthew G. Knepley PetscFunctionBegin; 104236bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 104336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 104436bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 104536bce6e8SMatthew G. Knepley #else 104636bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 104736bce6e8SMatthew G. Knepley #endif 104836bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 104936bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 105036bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1051c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1052de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 105336bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 105436bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 105536bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10567e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 105736bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 10583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 105936bce6e8SMatthew G. Knepley } 106036bce6e8SMatthew G. Knepley 106136bce6e8SMatthew G. Knepley /*@C 106236bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 106336bce6e8SMatthew G. Knepley 1064*20f4b53cSBarry Smith Not Collective 106536bce6e8SMatthew G. Knepley 106636bce6e8SMatthew G. Knepley Input Parameter: 1067811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...) 106836bce6e8SMatthew G. Knepley 106936bce6e8SMatthew G. Knepley Output Parameter: 1070811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 107136bce6e8SMatthew G. Knepley 107236bce6e8SMatthew G. Knepley Level: advanced 107336bce6e8SMatthew G. Knepley 1074d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 107536bce6e8SMatthew G. Knepley @*/ 1076d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 1077d71ae5a4SJacob Faibussowitsch { 107836bce6e8SMatthew G. Knepley PetscFunctionBegin; 107936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 108036bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 108136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 108236bce6e8SMatthew G. Knepley #else 108336bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 108436bce6e8SMatthew G. Knepley #endif 108536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 108636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 108736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 108836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 108936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 109036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 10917e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 109236bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 10933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109436bce6e8SMatthew G. Knepley } 109536bce6e8SMatthew G. Knepley 1096df863907SAlex Fikl /*@C 1097b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 109836bce6e8SMatthew G. Knepley 1099343a20b2SVaclav Hapla Collective 1100343a20b2SVaclav Hapla 110136bce6e8SMatthew G. Knepley Input Parameters: 1102811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 11034302210dSVaclav Hapla . parent - The parent dataset/group name 110436bce6e8SMatthew G. Knepley . name - The attribute name 110536bce6e8SMatthew G. Knepley . datatype - The attribute type 110636bce6e8SMatthew G. Knepley - value - The attribute value 110736bce6e8SMatthew G. Knepley 110836bce6e8SMatthew G. Knepley Level: advanced 110936bce6e8SMatthew G. Knepley 1110811af0c4SBarry Smith Note: 1111343a20b2SVaclav 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. 11124302210dSVaclav Hapla 1113d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, 1114811af0c4SBarry Smith `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 111536bce6e8SMatthew G. Knepley @*/ 1116d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 1117d71ae5a4SJacob Faibussowitsch { 11184302210dSVaclav Hapla char *parentAbsPath; 111960568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1120080f144cSVaclav Hapla PetscBool has; 112136bce6e8SMatthew G. Knepley 112236bce6e8SMatthew G. Knepley PetscFunctionBegin; 11235cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 11244302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1125c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 11264302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer, datatype, 4); 1127b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 112877717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 11299566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 11309566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 11319566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 11327e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 11337e97c476SMatthew G. Knepley size_t len; 11349566063dSJacob Faibussowitsch PetscCall(PetscStrlen((const char *)value, &len)); 1135792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 11367e97c476SMatthew G. Knepley } 11379566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1138792fecdfSBarry Smith PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR)); 1139792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1140080f144cSVaclav Hapla if (has) { 1141792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1142080f144cSVaclav Hapla } else { 1143792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1144080f144cSVaclav Hapla } 1145792fecdfSBarry Smith PetscCallHDF5(H5Awrite, (attribute, dtype, value)); 1146792fecdfSBarry Smith if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype)); 1147792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1148792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 1149792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (dataspace)); 11509566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 11513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115236bce6e8SMatthew G. Knepley } 115336bce6e8SMatthew G. Knepley 1154df863907SAlex Fikl /*@C 1155811af0c4SBarry Smith PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name 1156577e0e04SVaclav Hapla 1157343a20b2SVaclav Hapla Collective 1158343a20b2SVaclav Hapla 1159577e0e04SVaclav Hapla Input Parameters: 1160811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1161577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1162577e0e04SVaclav Hapla . name - The attribute name 1163577e0e04SVaclav Hapla . datatype - The attribute type 1164577e0e04SVaclav Hapla - value - The attribute value 1165577e0e04SVaclav Hapla 11663f423023SBarry Smith Level: advanced 11673f423023SBarry Smith 1168811af0c4SBarry Smith Note: 1169343a20b2SVaclav 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). 1170811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1171577e0e04SVaclav Hapla 1172d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, 1173811af0c4SBarry Smith `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1174577e0e04SVaclav Hapla @*/ 1175d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1176d71ae5a4SJacob Faibussowitsch { 1177577e0e04SVaclav Hapla PetscFunctionBegin; 1178577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1179577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1180577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 1181b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 11829566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 11839566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 11843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1185577e0e04SVaclav Hapla } 1186577e0e04SVaclav Hapla 1187577e0e04SVaclav Hapla /*@C 1188b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1189ad0c61feSMatthew G. Knepley 1190343a20b2SVaclav Hapla Collective 1191343a20b2SVaclav Hapla 1192ad0c61feSMatthew G. Knepley Input Parameters: 1193811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 11944302210dSVaclav Hapla . parent - The parent dataset/group name 1195ad0c61feSMatthew G. Knepley . name - The attribute name 1196a2d6be1bSVaclav Hapla . datatype - The attribute type 1197a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1198ad0c61feSMatthew G. Knepley 1199ad0c61feSMatthew G. Knepley Output Parameter: 1200a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1201ad0c61feSMatthew G. Knepley 12023f423023SBarry Smith Level: advanced 12033f423023SBarry Smith 1204a2d6be1bSVaclav Hapla Notes: 12053f423023SBarry Smith If defaultValue is `NULL` and the attribute is not found, an error occurs. 12063f423023SBarry Smith 12073f423023SBarry Smith If defaultValue is not `NULL` and the attribute is not found, `defaultValue` is copied to value. 12083f423023SBarry Smith 12093f423023SBarry Smith The pointers `defaultValue` and value can be the same; for instance 12103f423023SBarry Smith .vb 12113f423023SBarry Smith flg = PETSC_FALSE; 12123f423023SBarry Smith PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 12133f423023SBarry Smith .ve 1214a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1215a2d6be1bSVaclav Hapla 1216811af0c4SBarry Smith If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed. 1217708d4cb3SBarry Smith 12183f423023SBarry 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. 12194302210dSVaclav Hapla 1220d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1221ad0c61feSMatthew G. Knepley @*/ 1222d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1223d71ae5a4SJacob Faibussowitsch { 12244302210dSVaclav Hapla char *parentAbsPath; 1225a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1226969748fdSVaclav Hapla PetscBool has; 1227ad0c61feSMatthew G. Knepley 1228ad0c61feSMatthew G. Knepley PetscFunctionBegin; 12295cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 12304302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1231c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1232a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1233a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 12349566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 123577717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 12369566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 12379566063dSJacob Faibussowitsch if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1238a2d6be1bSVaclav Hapla if (!has) { 1239a2d6be1bSVaclav Hapla if (defaultValue) { 1240a2d6be1bSVaclav Hapla if (defaultValue != value) { 1241a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 12429566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value)); 1243a2d6be1bSVaclav Hapla } else { 1244a2d6be1bSVaclav Hapla size_t len; 1245792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype)); 12469566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(value, defaultValue, len)); 1247a2d6be1bSVaclav Hapla } 1248a2d6be1bSVaclav Hapla } 12499566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 12503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125198921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1252a2d6be1bSVaclav Hapla } 12539566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1254792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1255792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1256f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1257f0b58479SMatthew G. Knepley size_t len; 1258a2d6be1bSVaclav Hapla hid_t atype; 1259792fecdfSBarry Smith PetscCallHDF5Return(atype, H5Aget_type, (attribute)); 1260792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype)); 12619566063dSJacob Faibussowitsch PetscCall(PetscMalloc((len + 1) * sizeof(char), value)); 1262792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 1263792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value)); 1264708d4cb3SBarry Smith } else { 1265792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, value)); 1266708d4cb3SBarry Smith } 1267792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1268e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1269792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 12709566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 12713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1272ad0c61feSMatthew G. Knepley } 1273ad0c61feSMatthew G. Knepley 1274577e0e04SVaclav Hapla /*@C 1275811af0c4SBarry Smith PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name 1276577e0e04SVaclav Hapla 1277343a20b2SVaclav Hapla Collective 1278343a20b2SVaclav Hapla 1279577e0e04SVaclav Hapla Input Parameters: 1280811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1281577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1282577e0e04SVaclav Hapla . name - The attribute name 1283577e0e04SVaclav Hapla - datatype - The attribute type 1284577e0e04SVaclav Hapla 1285577e0e04SVaclav Hapla Output Parameter: 1286577e0e04SVaclav Hapla . value - The attribute value 1287577e0e04SVaclav Hapla 12883f423023SBarry Smith Level: advanced 12893f423023SBarry Smith 1290811af0c4SBarry Smith Note: 1291577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1292811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1293577e0e04SVaclav Hapla 1294d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1295577e0e04SVaclav Hapla @*/ 1296d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1297d71ae5a4SJacob Faibussowitsch { 1298577e0e04SVaclav Hapla PetscFunctionBegin; 1299577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1300577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1301577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 1302064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 13039566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 13049566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 13053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1306577e0e04SVaclav Hapla } 1307577e0e04SVaclav Hapla 1308d71ae5a4SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1309d71ae5a4SJacob Faibussowitsch { 1310a75c3fd4SVaclav Hapla htri_t exists; 1311a75c3fd4SVaclav Hapla hid_t group; 1312a75c3fd4SVaclav Hapla 1313a75c3fd4SVaclav Hapla PetscFunctionBegin; 1314792fecdfSBarry Smith PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT)); 1315792fecdfSBarry Smith if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT)); 1316a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1317792fecdfSBarry Smith PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1318792fecdfSBarry Smith PetscCallHDF5(H5Gclose, (group)); 1319a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1320a75c3fd4SVaclav Hapla } 1321a75c3fd4SVaclav Hapla *exists_ = (PetscBool)exists; 13223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1323a75c3fd4SVaclav Hapla } 1324a75c3fd4SVaclav Hapla 1325d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1326d71ae5a4SJacob Faibussowitsch { 132790e03537SVaclav Hapla const char rootGroupName[] = "/"; 13285cdeb108SMatthew G. Knepley hid_t h5; 1329e5bf9ebcSVaclav Hapla PetscBool exists = PETSC_FALSE; 133015b861d2SVaclav Hapla PetscInt i; 133115b861d2SVaclav Hapla int n; 133285688ae2SVaclav Hapla char **hierarchy; 133385688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 13345cdeb108SMatthew G. Knepley 13355cdeb108SMatthew G. Knepley PetscFunctionBegin; 13365cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 133790e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 133890e03537SVaclav Hapla else name = rootGroupName; 1339ccd66a83SVaclav Hapla if (has) { 1340064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 13415cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1342ccd66a83SVaclav Hapla } 1343ccd66a83SVaclav Hapla if (otype) { 1344064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 134556cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1346ccd66a83SVaclav Hapla } 13479566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 134885688ae2SVaclav Hapla 134985688ae2SVaclav Hapla /* 135085688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 135185688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 135285688ae2SVaclav Hapla 1) whether it's a valid link 135385688ae2SVaclav Hapla 2) whether this link resolves to an object 135485688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 135585688ae2SVaclav Hapla */ 13569566063dSJacob Faibussowitsch PetscCall(PetscStrToArray(name, '/', &n, &hierarchy)); 135785688ae2SVaclav Hapla if (!n) { 135885688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1359ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1360ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 13619566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 13623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136385688ae2SVaclav Hapla } 136485688ae2SVaclav Hapla for (i = 0; i < n; i++) { 1365c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, "/", sizeof(buf))); 1366c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf))); 13679566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1368a75c3fd4SVaclav Hapla if (!exists) break; 136985688ae2SVaclav Hapla } 13709566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 137185688ae2SVaclav Hapla 137285688ae2SVaclav Hapla /* If the object exists, get its type */ 1373ccd66a83SVaclav Hapla if (exists && otype) { 13745cdeb108SMatthew G. Knepley H5O_info_t info; 13755cdeb108SMatthew G. Knepley 137676276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1377792fecdfSBarry Smith PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT)); 137856cc0592SVaclav Hapla *otype = info.type; 13795cdeb108SMatthew G. Knepley } 1380ccd66a83SVaclav Hapla if (has) *has = exists; 13813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13825cdeb108SMatthew G. Knepley } 13835cdeb108SMatthew G. Knepley 13844302210dSVaclav Hapla /*@C 13850a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13860a9f38efSVaclav Hapla 1387343a20b2SVaclav Hapla Collective 1388343a20b2SVaclav Hapla 13890a9f38efSVaclav Hapla Input Parameters: 1390811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1391a0558868SVaclav Hapla - path - (Optional) The path relative to the pushed group 13920a9f38efSVaclav Hapla 13930a9f38efSVaclav Hapla Output Parameter: 13940a9f38efSVaclav Hapla . has - Flag for group existence 13950a9f38efSVaclav Hapla 13960a9f38efSVaclav Hapla Level: advanced 13970a9f38efSVaclav Hapla 13984302210dSVaclav Hapla Notes: 1399785c443eSVaclav 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. 14003f423023SBarry Smith `NULL` or empty path means the current pushed group. 14014302210dSVaclav Hapla 1402811af0c4SBarry Smith If path exists but is not a group, `PETSC_FALSE` is returned. 1403811af0c4SBarry Smith 1404d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 14050a9f38efSVaclav Hapla @*/ 1406d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 1407d71ae5a4SJacob Faibussowitsch { 14080a9f38efSVaclav Hapla H5O_type_t type; 14094302210dSVaclav Hapla char *abspath; 14100a9f38efSVaclav Hapla 14110a9f38efSVaclav Hapla PetscFunctionBegin; 14120a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 14134302210dSVaclav Hapla if (path) PetscValidCharPointer(path, 2); 14144302210dSVaclav Hapla PetscValidBoolPointer(has, 3); 141577717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 14169566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 14174302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 14189566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 14193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14200a9f38efSVaclav Hapla } 14210a9f38efSVaclav Hapla 142289e0ef10SVaclav Hapla /*@C 142389e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 142489e0ef10SVaclav Hapla 1425343a20b2SVaclav Hapla Collective 1426343a20b2SVaclav Hapla 142789e0ef10SVaclav Hapla Input Parameters: 1428811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 142989e0ef10SVaclav Hapla - path - The dataset path 143089e0ef10SVaclav Hapla 143189e0ef10SVaclav Hapla Output Parameter: 143289e0ef10SVaclav Hapla . has - Flag whether dataset exists 143389e0ef10SVaclav Hapla 143489e0ef10SVaclav Hapla Level: advanced 143589e0ef10SVaclav Hapla 143689e0ef10SVaclav Hapla Notes: 1437343a20b2SVaclav 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. 143889e0ef10SVaclav Hapla 14393f423023SBarry Smith If `path` is `NULL` or empty, has is set to `PETSC_FALSE`. 1440811af0c4SBarry Smith 14413f423023SBarry Smith If `path` exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1442811af0c4SBarry Smith 1443d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 144489e0ef10SVaclav Hapla @*/ 1445d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 1446d71ae5a4SJacob Faibussowitsch { 144789e0ef10SVaclav Hapla H5O_type_t type; 144889e0ef10SVaclav Hapla char *abspath; 144989e0ef10SVaclav Hapla 145089e0ef10SVaclav Hapla PetscFunctionBegin; 145189e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 145289e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path, 2); 145389e0ef10SVaclav Hapla PetscValidBoolPointer(has, 3); 145477717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 14559566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 145689e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 14579566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 14583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 145989e0ef10SVaclav Hapla } 146089e0ef10SVaclav Hapla 14610a9f38efSVaclav Hapla /*@ 1462e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1463ecce1506SVaclav Hapla 1464343a20b2SVaclav Hapla Collective 1465343a20b2SVaclav Hapla 1466ecce1506SVaclav Hapla Input Parameters: 1467811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1468ecce1506SVaclav Hapla - obj - The named object 1469ecce1506SVaclav Hapla 1470ecce1506SVaclav Hapla Output Parameter: 147189e0ef10SVaclav Hapla . has - Flag for dataset existence 1472ecce1506SVaclav Hapla 14733f423023SBarry Smith Level: advanced 14743f423023SBarry Smith 1475e3f143f7SVaclav Hapla Notes: 147689e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1477811af0c4SBarry Smith 1478811af0c4SBarry Smith If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1479e3f143f7SVaclav Hapla 1480d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1481ecce1506SVaclav Hapla @*/ 1482d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1483d71ae5a4SJacob Faibussowitsch { 148489e0ef10SVaclav Hapla size_t len; 1485ecce1506SVaclav Hapla 1486ecce1506SVaclav Hapla PetscFunctionBegin; 1487c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1488c57a1a9aSVaclav Hapla PetscValidHeader(obj, 2); 14894302210dSVaclav Hapla PetscValidBoolPointer(has, 3); 14909566063dSJacob Faibussowitsch PetscCall(PetscStrlen(obj->name, &len)); 14915f80ce2aSJacob Faibussowitsch PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 14929566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 14933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1494ecce1506SVaclav Hapla } 1495ecce1506SVaclav Hapla 1496df863907SAlex Fikl /*@C 1497ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1498e7bdbf8eSMatthew G. Knepley 1499343a20b2SVaclav Hapla Collective 1500343a20b2SVaclav Hapla 1501e7bdbf8eSMatthew G. Knepley Input Parameters: 1502811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 15034302210dSVaclav Hapla . parent - The parent dataset/group name 1504e7bdbf8eSMatthew G. Knepley - name - The attribute name 1505e7bdbf8eSMatthew G. Knepley 1506e7bdbf8eSMatthew G. Knepley Output Parameter: 1507e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1508e7bdbf8eSMatthew G. Knepley 1509e7bdbf8eSMatthew G. Knepley Level: advanced 1510e7bdbf8eSMatthew G. Knepley 1511811af0c4SBarry Smith Note: 15123f423023SBarry 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. 15134302210dSVaclav Hapla 1514d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1515e7bdbf8eSMatthew G. Knepley @*/ 1516d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1517d71ae5a4SJacob Faibussowitsch { 15184302210dSVaclav Hapla char *parentAbsPath; 1519e7bdbf8eSMatthew G. Knepley 1520e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 15215cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 15224302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1523c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 15244302210dSVaclav Hapla PetscValidBoolPointer(has, 4); 152577717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 15269566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 15279566063dSJacob Faibussowitsch if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 15289566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 15293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153006db490cSVaclav Hapla } 153106db490cSVaclav Hapla 1532577e0e04SVaclav Hapla /*@C 1533811af0c4SBarry Smith PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name 1534577e0e04SVaclav Hapla 1535343a20b2SVaclav Hapla Collective 1536343a20b2SVaclav Hapla 1537577e0e04SVaclav Hapla Input Parameters: 1538811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1539577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1540577e0e04SVaclav Hapla - name - The attribute name 1541577e0e04SVaclav Hapla 1542577e0e04SVaclav Hapla Output Parameter: 1543577e0e04SVaclav Hapla . has - Flag for attribute existence 1544577e0e04SVaclav Hapla 15453f423023SBarry Smith Level: advanced 15463f423023SBarry Smith 1547811af0c4SBarry Smith Note: 1548577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1549811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1550577e0e04SVaclav Hapla 1551d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1552577e0e04SVaclav Hapla @*/ 1553d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1554d71ae5a4SJacob Faibussowitsch { 1555577e0e04SVaclav Hapla PetscFunctionBegin; 1556577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1557577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1558577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 15594302210dSVaclav Hapla PetscValidBoolPointer(has, 4); 15609566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 15619566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 15623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1563577e0e04SVaclav Hapla } 1564577e0e04SVaclav Hapla 1565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1566d71ae5a4SJacob Faibussowitsch { 156706db490cSVaclav Hapla hid_t h5; 156806db490cSVaclav Hapla htri_t hhas; 156906db490cSVaclav Hapla 157006db490cSVaclav Hapla PetscFunctionBegin; 15719566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1572792fecdfSBarry Smith PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT)); 15735cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 15743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1575e7bdbf8eSMatthew G. Knepley } 1576e7bdbf8eSMatthew G. Knepley 1577a75e6a4aSMatthew G. Knepley /* 1578a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1579a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1580a75e6a4aSMatthew G. Knepley */ 1581d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1582a75e6a4aSMatthew G. Knepley 1583a75e6a4aSMatthew G. Knepley /*@C 1584811af0c4SBarry Smith PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator. 1585a75e6a4aSMatthew G. Knepley 1586d083f849SBarry Smith Collective 1587a75e6a4aSMatthew G. Knepley 1588a75e6a4aSMatthew G. Knepley Input Parameter: 15893f423023SBarry Smith . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer` 1590a75e6a4aSMatthew G. Knepley 1591811af0c4SBarry Smith Options Database Key: 159210699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1593a75e6a4aSMatthew G. Knepley 1594811af0c4SBarry Smith Environmental variable: 1595811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file 1596a75e6a4aSMatthew G. Knepley 1597*20f4b53cSBarry Smith Level: intermediate 1598*20f4b53cSBarry Smith 1599811af0c4SBarry Smith Note: 1600811af0c4SBarry Smith Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return 1601811af0c4SBarry Smith an error code. The HDF5 `PetscViewer` is usually used in the form 1602a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1603a75e6a4aSMatthew G. Knepley 1604d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1605a75e6a4aSMatthew G. Knepley @*/ 1606d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1607d71ae5a4SJacob Faibussowitsch { 1608a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 16093ba16761SJacob Faibussowitsch PetscMPIInt mpi_ierr; 1610a75e6a4aSMatthew G. Knepley PetscBool flg; 1611a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1612a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1613a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1614a75e6a4aSMatthew G. Knepley 1615a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 16169371c9d4SSatish Balay ierr = PetscCommDuplicate(comm, &ncomm, NULL); 16179371c9d4SSatish Balay if (ierr) { 16183ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16199371c9d4SSatish Balay PetscFunctionReturn(NULL); 16209371c9d4SSatish Balay } 1621a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 16223ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL); 16233ba16761SJacob Faibussowitsch if (mpi_ierr) { 16243ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16259371c9d4SSatish Balay PetscFunctionReturn(NULL); 16269371c9d4SSatish Balay } 1627a75e6a4aSMatthew G. Knepley } 16283ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg); 16293ba16761SJacob Faibussowitsch if (mpi_ierr) { 16303ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16319371c9d4SSatish Balay PetscFunctionReturn(NULL); 16329371c9d4SSatish Balay } 1633a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1634a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg); 16359371c9d4SSatish Balay if (ierr) { 16363ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16379371c9d4SSatish Balay PetscFunctionReturn(NULL); 16389371c9d4SSatish Balay } 1639a75e6a4aSMatthew G. Knepley if (!flg) { 1640c6a7a370SJeremy L Thompson ierr = PetscStrncpy(fname, "output.h5", sizeof(fname)); 16419371c9d4SSatish Balay if (ierr) { 16423ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16439371c9d4SSatish Balay PetscFunctionReturn(NULL); 16449371c9d4SSatish Balay } 1645a75e6a4aSMatthew G. Knepley } 1646a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer); 16479371c9d4SSatish Balay if (ierr) { 16483ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16499371c9d4SSatish Balay PetscFunctionReturn(NULL); 16509371c9d4SSatish Balay } 1651a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16529371c9d4SSatish Balay if (ierr) { 16533ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16549371c9d4SSatish Balay PetscFunctionReturn(NULL); 16559371c9d4SSatish Balay } 16563ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer); 16573ba16761SJacob Faibussowitsch if (mpi_ierr) { 16583ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16599371c9d4SSatish Balay PetscFunctionReturn(NULL); 16609371c9d4SSatish Balay } 1661a75e6a4aSMatthew G. Knepley } 1662a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 16639371c9d4SSatish Balay if (ierr) { 16643ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16659371c9d4SSatish Balay PetscFunctionReturn(NULL); 16669371c9d4SSatish Balay } 1667a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1668a75e6a4aSMatthew G. Knepley } 1669