16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 25c6c1daeSBarry Smith 3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *); 406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *); 534e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer, const char *[]); 606db490cSVaclav Hapla 777717648SVaclav Hapla /*@C 877717648SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`. 977717648SVaclav Hapla 1020f4b53cSBarry Smith Not Collective 1177717648SVaclav Hapla 1277717648SVaclav Hapla Input Parameters: 1377717648SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 1477717648SVaclav Hapla - path - (Optional) The path relative to the pushed group 1577717648SVaclav Hapla 1677717648SVaclav Hapla Output Parameter: 1777717648SVaclav Hapla . abspath - The absolute HDF5 path (group) 1877717648SVaclav Hapla 1977717648SVaclav Hapla Level: intermediate 2077717648SVaclav Hapla 2177717648SVaclav Hapla Notes: 2277717648SVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 233f423023SBarry Smith `NULL` or empty path means the current pushed group. 2477717648SVaclav Hapla 2577717648SVaclav Hapla The output abspath is newly allocated so needs to be freed. 2677717648SVaclav Hapla 27d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 2877717648SVaclav Hapla @*/ 29ce78bad3SBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], const char *abspath[]) 30d71ae5a4SJacob Faibussowitsch { 3177717648SVaclav Hapla size_t len; 324302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 336c132bc1SVaclav Hapla const char *group; 346c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 356c132bc1SVaclav Hapla 366c132bc1SVaclav Hapla PetscFunctionBegin; 3777717648SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 384f572ea9SToby Isaac if (path) PetscAssertPointer(path, 2); 394f572ea9SToby Isaac PetscAssertPointer(abspath, 3); 4077717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group)); 4177717648SVaclav Hapla PetscCall(PetscStrlen(path, &len)); 4277717648SVaclav Hapla relative = (PetscBool)(!len || path[0] != '/'); 434302210dSVaclav Hapla if (relative) { 44c6a7a370SJeremy L Thompson PetscCall(PetscStrncpy(buf, group, sizeof(buf))); 45c6a7a370SJeremy L Thompson if (!group || len) PetscCall(PetscStrlcat(buf, "/", sizeof(buf))); 46c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, path, sizeof(buf))); 47ce78bad3SBarry Smith PetscCall(PetscStrallocpy(buf, (char **)abspath)); 484302210dSVaclav Hapla } else { 49ce78bad3SBarry Smith PetscCall(PetscStrallocpy(path, (char **)abspath)); 504302210dSVaclav Hapla } 513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 526c132bc1SVaclav Hapla } 536c132bc1SVaclav Hapla 54d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 55d71ae5a4SJacob Faibussowitsch { 56577e0e04SVaclav Hapla PetscBool has; 57577e0e04SVaclav Hapla 58577e0e04SVaclav Hapla PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has)); 60577e0e04SVaclav Hapla if (!has) { 61ce78bad3SBarry Smith const char *group; 6277717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 6377717648SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group); 64577e0e04SVaclav Hapla } 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66577e0e04SVaclav Hapla } 67577e0e04SVaclav Hapla 68ce78bad3SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems PetscOptionsObject) 69d71ae5a4SJacob Faibussowitsch { 70ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 7182ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 7282ea9b62SBarry Smith 7382ea9b62SBarry Smith PetscFunctionBegin; 74d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options"); 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL)); 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL)); 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set)); 789566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg)); 7919a20e4cSMatthew G. Knepley flg = PETSC_FALSE; 809566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set)); 819566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg)); 82d0609cedSBarry Smith PetscOptionsHeadEnd(); 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8482ea9b62SBarry Smith } 8582ea9b62SBarry Smith 86d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer) 87d71ae5a4SJacob Faibussowitsch { 881b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 8903fa8834SVaclav Hapla PetscBool flg; 901b793a25SVaclav Hapla 911b793a25SVaclav Hapla PetscFunctionBegin; 9248a46eb9SPierre Jolivet if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename)); 939566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2])); 949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput])); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetCollective(v, &flg)); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent")); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping])); 983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 991b793a25SVaclav Hapla } 1001b793a25SVaclav Hapla 101d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 102d71ae5a4SJacob Faibussowitsch { 1035c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1045c6c1daeSBarry Smith 1055c6c1daeSBarry Smith PetscFunctionBegin; 1069566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->filename)); 107792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1095c6c1daeSBarry Smith } 1105c6c1daeSBarry Smith 111d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) 112d71ae5a4SJacob Faibussowitsch { 1136226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1146226335aSVaclav Hapla 1156226335aSVaclav Hapla PetscFunctionBegin; 116792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL)); 1173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1186226335aSVaclav Hapla } 1196226335aSVaclav Hapla 120d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 121d71ae5a4SJacob Faibussowitsch { 1225c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1235c6c1daeSBarry Smith 1245c6c1daeSBarry Smith PetscFunctionBegin; 125792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (hdf5->dxpl_id)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_HDF5(viewer)); 1275c6c1daeSBarry Smith while (hdf5->groups) { 1287d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1295c6c1daeSBarry Smith 1309566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups->name)); 1319566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups)); 1325c6c1daeSBarry Smith hdf5->groups = tmp; 1335c6c1daeSBarry Smith } 1349566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5)); 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 1379566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 1382e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 1399566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL)); 1409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL)); 1412e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL)); 1422e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL)); 1449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL)); 1453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1465c6c1daeSBarry Smith } 1475c6c1daeSBarry Smith 148d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 149d71ae5a4SJacob Faibussowitsch { 1505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1515c6c1daeSBarry Smith 1525c6c1daeSBarry Smith PetscFunctionBegin; 1535c6c1daeSBarry Smith hdf5->btype = type; 1543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1555c6c1daeSBarry Smith } 1565c6c1daeSBarry Smith 157d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 158d71ae5a4SJacob Faibussowitsch { 1590b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1600b62783dSVaclav Hapla 1610b62783dSVaclav Hapla PetscFunctionBegin; 1620b62783dSVaclav Hapla *type = hdf5->btype; 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1640b62783dSVaclav Hapla } 1650b62783dSVaclav Hapla 166d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 167d71ae5a4SJacob Faibussowitsch { 16882ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 16982ea9b62SBarry Smith 17082ea9b62SBarry Smith PetscFunctionBegin; 17182ea9b62SBarry Smith hdf5->basedimension2 = flg; 1723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 17382ea9b62SBarry Smith } 17482ea9b62SBarry Smith 175df863907SAlex Fikl /*@ 17682ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17782ea9b62SBarry Smith dimension of 2. 17882ea9b62SBarry Smith 179c3339decSBarry Smith Logically Collective 18082ea9b62SBarry Smith 18182ea9b62SBarry Smith Input Parameters: 182811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored 183811af0c4SBarry 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 18482ea9b62SBarry Smith 185811af0c4SBarry Smith Options Database Key: 18682ea9b62SBarry 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 18782ea9b62SBarry Smith 1883f423023SBarry Smith Level: intermediate 1893f423023SBarry Smith 190811af0c4SBarry Smith Note: 19195452b02SPatrick 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 19282ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19382ea9b62SBarry Smith 194aec76313SJacob Faibussowitsch .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 19582ea9b62SBarry Smith @*/ 196d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg) 197d71ae5a4SJacob Faibussowitsch { 19882ea9b62SBarry Smith PetscFunctionBegin; 19982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 200cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg)); 2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20282ea9b62SBarry Smith } 20382ea9b62SBarry Smith 204df863907SAlex Fikl /*@ 20582ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 20682ea9b62SBarry Smith dimension of 2. 20782ea9b62SBarry Smith 208c3339decSBarry Smith Logically Collective 20982ea9b62SBarry Smith 21082ea9b62SBarry Smith Input Parameter: 211811af0c4SBarry Smith . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5` 21282ea9b62SBarry Smith 21382ea9b62SBarry Smith Output Parameter: 214811af0c4SBarry 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 21582ea9b62SBarry Smith 2163f423023SBarry Smith Level: intermediate 2173f423023SBarry Smith 218811af0c4SBarry Smith Note: 21995452b02SPatrick 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 22082ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 22182ea9b62SBarry Smith 222d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 22382ea9b62SBarry Smith @*/ 224d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg) 225d71ae5a4SJacob Faibussowitsch { 22682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 22782ea9b62SBarry Smith 22882ea9b62SBarry Smith PetscFunctionBegin; 22982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 23082ea9b62SBarry Smith *flg = hdf5->basedimension2; 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23282ea9b62SBarry Smith } 23382ea9b62SBarry Smith 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 235d71ae5a4SJacob Faibussowitsch { 2369a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2379a0502c6SHåkon Strandenes 2389a0502c6SHåkon Strandenes PetscFunctionBegin; 2399a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2419a0502c6SHåkon Strandenes } 2429a0502c6SHåkon Strandenes 243df863907SAlex Fikl /*@ 2449a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 245811af0c4SBarry Smith compiled with double precision `PetscReal`. 2469a0502c6SHåkon Strandenes 247c3339decSBarry Smith Logically Collective 2489a0502c6SHåkon Strandenes 2499a0502c6SHåkon Strandenes Input Parameters: 250811af0c4SBarry Smith + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored 251811af0c4SBarry Smith - flg - if `PETSC_TRUE` the data will be written to disk with single precision 2529a0502c6SHåkon Strandenes 253811af0c4SBarry Smith Options Database Key: 2549a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2559a0502c6SHåkon Strandenes 2563f423023SBarry Smith Level: intermediate 2573f423023SBarry Smith 258811af0c4SBarry Smith Note: 25995452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2609a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2619a0502c6SHåkon Strandenes 262d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 263811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5GetSPOutput()` 2649a0502c6SHåkon Strandenes @*/ 265d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg) 266d71ae5a4SJacob Faibussowitsch { 2679a0502c6SHåkon Strandenes PetscFunctionBegin; 2689a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 269cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg)); 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2719a0502c6SHåkon Strandenes } 2729a0502c6SHåkon Strandenes 273df863907SAlex Fikl /*@ 2749a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 275811af0c4SBarry Smith compiled with double precision `PetscReal`. 2769a0502c6SHåkon Strandenes 277c3339decSBarry Smith Logically Collective 2789a0502c6SHåkon Strandenes 2799a0502c6SHåkon Strandenes Input Parameter: 280811af0c4SBarry Smith . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5` 2819a0502c6SHåkon Strandenes 2829a0502c6SHåkon Strandenes Output Parameter: 283811af0c4SBarry Smith . flg - if `PETSC_TRUE` the data will be written to disk with single precision 2849a0502c6SHåkon Strandenes 2859a0502c6SHåkon Strandenes Level: intermediate 2869a0502c6SHåkon Strandenes 287d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 288811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5SetSPOutput()` 2899a0502c6SHåkon Strandenes @*/ 290d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg) 291d71ae5a4SJacob Faibussowitsch { 2929a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2939a0502c6SHåkon Strandenes 2949a0502c6SHåkon Strandenes PetscFunctionBegin; 2959a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 2969a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2989a0502c6SHåkon Strandenes } 2999a0502c6SHåkon Strandenes 300d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 301d71ae5a4SJacob Faibussowitsch { 302ee8b9732SVaclav Hapla PetscFunctionBegin; 303ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 304ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 305c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL) 306ee8b9732SVaclav Hapla { 307ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 308792fecdfSBarry Smith PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 309ee8b9732SVaclav Hapla } 310ee8b9732SVaclav Hapla #else 3119566063dSJacob 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")); 312ee8b9732SVaclav Hapla #endif 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 314ee8b9732SVaclav Hapla } 315ee8b9732SVaclav Hapla 316ee8b9732SVaclav Hapla /*@ 317ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 318ee8b9732SVaclav Hapla 319ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 320ee8b9732SVaclav Hapla 321ee8b9732SVaclav Hapla Input Parameters: 322811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 323811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default) 324ee8b9732SVaclav Hapla 325811af0c4SBarry Smith Options Database Key: 326ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 327ee8b9732SVaclav Hapla 3283f423023SBarry Smith Level: intermediate 3293f423023SBarry Smith 330811af0c4SBarry Smith Note: 331ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 33253effed7SBarry 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. 333ee8b9732SVaclav Hapla 334aec76313SJacob Faibussowitsch Developer Notes: 335811af0c4SBarry Smith In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively. 336ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 337ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 338ee8b9732SVaclav Hapla 339d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 340ee8b9732SVaclav Hapla @*/ 341d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg) 342d71ae5a4SJacob Faibussowitsch { 343ee8b9732SVaclav Hapla PetscFunctionBegin; 344ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 345ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer, flg, 2); 346cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg)); 3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 348ee8b9732SVaclav Hapla } 349ee8b9732SVaclav Hapla 350d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 351d71ae5a4SJacob Faibussowitsch { 352c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 353ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 354ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 355c796909eSBarry Smith #endif 356ee8b9732SVaclav Hapla 357ee8b9732SVaclav Hapla PetscFunctionBegin; 358c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 359c796909eSBarry Smith *flg = PETSC_FALSE; 360c796909eSBarry Smith #else 361792fecdfSBarry Smith PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode)); 362ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 363c796909eSBarry Smith #endif 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 365ee8b9732SVaclav Hapla } 366ee8b9732SVaclav Hapla 367ee8b9732SVaclav Hapla /*@ 368ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 369ee8b9732SVaclav Hapla 370ee8b9732SVaclav Hapla Not Collective 371ee8b9732SVaclav Hapla 37220f4b53cSBarry Smith Input Parameter: 373811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer` 374ee8b9732SVaclav Hapla 37520f4b53cSBarry Smith Output Parameter: 376ee8b9732SVaclav Hapla . flg - the flag 377ee8b9732SVaclav Hapla 378ee8b9732SVaclav Hapla Level: intermediate 379ee8b9732SVaclav Hapla 380811af0c4SBarry Smith Note: 381811af0c4SBarry 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. 382811af0c4SBarry Smith For more details, see `PetscViewerHDF5SetCollective()`. 383ee8b9732SVaclav Hapla 384d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 385ee8b9732SVaclav Hapla @*/ 386d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg) 387d71ae5a4SJacob Faibussowitsch { 388ee8b9732SVaclav Hapla PetscFunctionBegin; 389ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 3904f572ea9SToby Isaac PetscAssertPointer(flg, 2); 391ee8b9732SVaclav Hapla 392cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg)); 3933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 394ee8b9732SVaclav Hapla } 395ee8b9732SVaclav Hapla 396d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 397d71ae5a4SJacob Faibussowitsch { 3985c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 3995c6c1daeSBarry Smith hid_t plist_id; 4005c6c1daeSBarry Smith 4015c6c1daeSBarry Smith PetscFunctionBegin; 402792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 4039566063dSJacob Faibussowitsch if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 4049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &hdf5->filename)); 4055c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 406792fecdfSBarry Smith PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS)); 407c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 408792fecdfSBarry Smith PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 409c796909eSBarry Smith #endif 4105c6c1daeSBarry Smith /* Create or open the file collectively */ 4115c6c1daeSBarry Smith switch (hdf5->btype) { 4125c6c1daeSBarry Smith case FILE_MODE_READ: 4138a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 4148a2871f6SBarry Smith PetscMPIInt rank; 4158a2871f6SBarry Smith PetscBool flg; 4168a2871f6SBarry Smith 4179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 4188a2871f6SBarry Smith if (rank == 0) { 4199566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4205f80ce2aSJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename); 4218a2871f6SBarry Smith } 4229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 4238a2871f6SBarry Smith } 424792fecdfSBarry Smith PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id)); 4255c6c1daeSBarry Smith break; 4265c6c1daeSBarry Smith case FILE_MODE_APPEND: 4279371c9d4SSatish Balay case FILE_MODE_UPDATE: { 4287e4fd573SVaclav Hapla PetscBool flg; 4299566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 430792fecdfSBarry Smith if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id)); 431792fecdfSBarry Smith else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4325c6c1daeSBarry Smith break; 4337e4fd573SVaclav Hapla } 434d71ae5a4SJacob Faibussowitsch case FILE_MODE_WRITE: 435d71ae5a4SJacob Faibussowitsch PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 436d71ae5a4SJacob Faibussowitsch break; 437d71ae5a4SJacob Faibussowitsch case FILE_MODE_UNDEFINED: 438d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 439d71ae5a4SJacob Faibussowitsch default: 440d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]); 4415c6c1daeSBarry Smith } 4425f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 443792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (plist_id)); 444671abe24SVaclav Hapla PetscCall(PetscViewerHDF5ResetAttachedDMPlexStorageVersion(viewer)); 4453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4465c6c1daeSBarry Smith } 4475c6c1daeSBarry Smith 448d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name) 449d71ae5a4SJacob Faibussowitsch { 450d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data; 451d1232d7fSToby Isaac 452d1232d7fSToby Isaac PetscFunctionBegin; 453d1232d7fSToby Isaac *name = vhdf5->filename; 4543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 455d1232d7fSToby Isaac } 456d1232d7fSToby Isaac 457d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 458d71ae5a4SJacob Faibussowitsch { 4595dc64a97SVaclav Hapla /* 460b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4615dc64a97SVaclav Hapla */ 462b723ab35SVaclav Hapla 463b723ab35SVaclav Hapla PetscFunctionBegin; 4643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 465b723ab35SVaclav Hapla } 466b723ab35SVaclav Hapla 467d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 468d71ae5a4SJacob Faibussowitsch { 46919a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 47019a20e4cSMatthew G. Knepley 47119a20e4cSMatthew G. Knepley PetscFunctionBegin; 47219a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 4733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47419a20e4cSMatthew G. Knepley } 47519a20e4cSMatthew G. Knepley 476d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 477d71ae5a4SJacob Faibussowitsch { 47819a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 47919a20e4cSMatthew G. Knepley 48019a20e4cSMatthew G. Knepley PetscFunctionBegin; 48119a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48319a20e4cSMatthew G. Knepley } 48419a20e4cSMatthew G. Knepley 48519a20e4cSMatthew G. Knepley /*@ 48619a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 48719a20e4cSMatthew G. Knepley 488c3339decSBarry Smith Logically Collective 48919a20e4cSMatthew G. Knepley 49019a20e4cSMatthew G. Knepley Input Parameters: 491811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 492811af0c4SBarry Smith - flg - if `PETSC_TRUE` we will assume that timestepping is on 49319a20e4cSMatthew G. Knepley 494811af0c4SBarry Smith Options Database Key: 49519a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 49619a20e4cSMatthew G. Knepley 4973f423023SBarry Smith Level: intermediate 4983f423023SBarry Smith 499811af0c4SBarry Smith Note: 50019a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 50119a20e4cSMatthew G. Knepley 502d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 50319a20e4cSMatthew G. Knepley @*/ 504d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 505d71ae5a4SJacob Faibussowitsch { 50619a20e4cSMatthew G. Knepley PetscFunctionBegin; 50719a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 508cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg)); 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51019a20e4cSMatthew G. Knepley } 51119a20e4cSMatthew G. Knepley 51219a20e4cSMatthew G. Knepley /*@ 51319a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 51419a20e4cSMatthew G. Knepley 51520f4b53cSBarry Smith Not Collective 51619a20e4cSMatthew G. Knepley 51719a20e4cSMatthew G. Knepley Input Parameter: 518811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 51919a20e4cSMatthew G. Knepley 52019a20e4cSMatthew G. Knepley Output Parameter: 521811af0c4SBarry Smith . flg - if `PETSC_TRUE` we will assume that timestepping is on 52219a20e4cSMatthew G. Knepley 52319a20e4cSMatthew G. Knepley Level: intermediate 52419a20e4cSMatthew G. Knepley 525d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 52619a20e4cSMatthew G. Knepley @*/ 527d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 528d71ae5a4SJacob Faibussowitsch { 52919a20e4cSMatthew G. Knepley PetscFunctionBegin; 53019a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 531cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg)); 5323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53319a20e4cSMatthew G. Knepley } 53419a20e4cSMatthew G. Knepley 5358556b5ebSBarry Smith /*MC 5368556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5378556b5ebSBarry Smith 538811af0c4SBarry Smith Level: beginner 539811af0c4SBarry Smith 540d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 541db781477SPatrick Sanan `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 542db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 543db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 5448556b5ebSBarry Smith M*/ 545d1232d7fSToby Isaac 546d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 547d71ae5a4SJacob Faibussowitsch { 5485c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5495c6c1daeSBarry Smith 5505c6c1daeSBarry Smith PetscFunctionBegin; 55199335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 55299335c34SVaclav Hapla { 55399335c34SVaclav Hapla PetscMPIInt size; 5549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 5555f80ce2aSJacob 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)"); 55699335c34SVaclav Hapla } 55799335c34SVaclav Hapla #endif 55899335c34SVaclav Hapla 5594dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&hdf5)); 5605c6c1daeSBarry Smith 5615c6c1daeSBarry Smith v->data = (void *)hdf5; 5625c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 56382ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 564b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5651b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5666226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5677e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 568908793a3SLisandro Dalcin hdf5->filename = NULL; 5695c6c1daeSBarry Smith hdf5->timestep = -1; 5700298fd71SBarry Smith hdf5->groups = NULL; 5715c6c1daeSBarry Smith 572792fecdfSBarry Smith PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER)); 5739c5072fbSVaclav Hapla 5749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5)); 5759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5)); 5769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5)); 5779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5)); 5789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5)); 5799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5)); 5809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5)); 5819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5)); 5829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5)); 5839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5)); 5843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5855c6c1daeSBarry Smith } 5865c6c1daeSBarry Smith 5875c6c1daeSBarry Smith /*@C 588811af0c4SBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer` 5895c6c1daeSBarry Smith 590d083f849SBarry Smith Collective 5915c6c1daeSBarry Smith 5925c6c1daeSBarry Smith Input Parameters: 5935c6c1daeSBarry Smith + comm - MPI communicator 5945c6c1daeSBarry Smith . name - name of file 5955c6c1daeSBarry Smith - type - type of file 5965c6c1daeSBarry Smith 5975c6c1daeSBarry Smith Output Parameter: 598811af0c4SBarry Smith . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file 5995c6c1daeSBarry Smith 600811af0c4SBarry Smith Options Database Keys: 601a2b725a8SWilliam 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 602a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 60382ea9b62SBarry Smith 6045c6c1daeSBarry Smith Level: beginner 6055c6c1daeSBarry Smith 6067e4fd573SVaclav Hapla Notes: 6077e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 608811af0c4SBarry Smith .vb 609811af0c4SBarry Smith FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 610811af0c4SBarry Smith FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 611811af0c4SBarry 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] 612811af0c4SBarry Smith FILE_MODE_UPDATE - same as FILE_MODE_APPEND 613811af0c4SBarry Smith .ve 6147e4fd573SVaclav Hapla 615da81f932SPierre 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. 6167e4fd573SVaclav Hapla 6175c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 6185c6c1daeSBarry Smith 619d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 620db781477SPatrick Sanan `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 621db781477SPatrick Sanan `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 6225c6c1daeSBarry Smith @*/ 623d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 624d71ae5a4SJacob Faibussowitsch { 6255c6c1daeSBarry Smith PetscFunctionBegin; 6269566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, hdf5v)); 6279566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 6289566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 6299566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*hdf5v, name)); 6309566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*hdf5v)); 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6325c6c1daeSBarry Smith } 6335c6c1daeSBarry Smith 6345c6c1daeSBarry Smith /*@C 6355c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6365c6c1daeSBarry Smith 63720f4b53cSBarry Smith Not Collective 6385c6c1daeSBarry Smith 6395c6c1daeSBarry Smith Input Parameter: 640811af0c4SBarry Smith . viewer - the `PetscViewer` 6415c6c1daeSBarry Smith 6425c6c1daeSBarry Smith Output Parameter: 6435c6c1daeSBarry Smith . file_id - The file id 6445c6c1daeSBarry Smith 6455c6c1daeSBarry Smith Level: intermediate 6465c6c1daeSBarry Smith 647d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()` 6485c6c1daeSBarry Smith @*/ 649d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 650d71ae5a4SJacob Faibussowitsch { 6515c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6525c6c1daeSBarry Smith 6535c6c1daeSBarry Smith PetscFunctionBegin; 6545c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 6555c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6575c6c1daeSBarry Smith } 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith /*@C 6605c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6615c6c1daeSBarry Smith 66220f4b53cSBarry Smith Not Collective 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith Input Parameters: 665811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 6665c6c1daeSBarry Smith - name - The group name 6675c6c1daeSBarry Smith 6685c6c1daeSBarry Smith Level: intermediate 6695c6c1daeSBarry Smith 6702e361344SVaclav Hapla Notes: 6712e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 672811af0c4SBarry Smith .vb 673811af0c4SBarry 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. 6743f423023SBarry Smith `NULL`, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 675811af0c4SBarry Smith "." means the current group is pushed again. 676811af0c4SBarry Smith .ve 6772e361344SVaclav Hapla 6782e361344SVaclav Hapla Example: 6792e361344SVaclav Hapla Suppose the current group is "/a". 6803f423023SBarry Smith .vb 6813f423023SBarry Smith If name is `NULL`, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6823f423023SBarry Smith If name is ".", then the new group will be "/a". 6833f423023SBarry Smith If name is "b", then the new group will be "/a/b". 6843f423023SBarry Smith If name is "/b", then the new group will be "/b". 6853f423023SBarry Smith .ve 6862e361344SVaclav Hapla 687aec76313SJacob Faibussowitsch Developer Notes: 6883f423023SBarry Smith The root group "/" is internally stored as `NULL`. 689820fc2d1SVaclav Hapla 690d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 6915c6c1daeSBarry Smith @*/ 692d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 693d71ae5a4SJacob Faibussowitsch { 6945c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6957d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6962e361344SVaclav Hapla size_t i, len; 6972e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6982e361344SVaclav Hapla const char *gname; 6995c6c1daeSBarry Smith 7005c6c1daeSBarry Smith PetscFunctionBegin; 7015c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 7024f572ea9SToby Isaac if (name) PetscAssertPointer(name, 2); 7039566063dSJacob Faibussowitsch PetscCall(PetscStrlen(name, &len)); 7042e361344SVaclav Hapla gname = NULL; 7052e361344SVaclav Hapla if (len) { 7062e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 7072e361344SVaclav Hapla /* use current name */ 7082e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 7092e361344SVaclav Hapla } else if (name[0] == '/') { 7102e361344SVaclav Hapla /* absolute */ 7112e361344SVaclav Hapla for (i = 1; i < len; i++) { 7122e361344SVaclav Hapla if (name[i] != '/') { 7132e361344SVaclav Hapla gname = name; 7142e361344SVaclav Hapla break; 7152e361344SVaclav Hapla } 7162e361344SVaclav Hapla } 7172e361344SVaclav Hapla } else { 7182e361344SVaclav Hapla /* relative */ 7192e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 7209566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 7212e361344SVaclav Hapla gname = buf; 7222e361344SVaclav Hapla } 7232e361344SVaclav Hapla } 7249566063dSJacob Faibussowitsch PetscCall(PetscNew(&groupNode)); 7259566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name)); 7265c6c1daeSBarry Smith groupNode->next = hdf5->groups; 7275c6c1daeSBarry Smith hdf5->groups = groupNode; 7283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7295c6c1daeSBarry Smith } 7305c6c1daeSBarry Smith 7313ef9c667SSatish Balay /*@ 7325c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 7335c6c1daeSBarry Smith 73420f4b53cSBarry Smith Not Collective 7355c6c1daeSBarry Smith 7365c6c1daeSBarry Smith Input Parameter: 737811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7385c6c1daeSBarry Smith 7395c6c1daeSBarry Smith Level: intermediate 7405c6c1daeSBarry Smith 741d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 7425c6c1daeSBarry Smith @*/ 743d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 744d71ae5a4SJacob Faibussowitsch { 7455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7467d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7475c6c1daeSBarry Smith 7485c6c1daeSBarry Smith PetscFunctionBegin; 7495c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 7505f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7515c6c1daeSBarry Smith groupNode = hdf5->groups; 7525c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7539566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode->name)); 7549566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode)); 7553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7565c6c1daeSBarry Smith } 7575c6c1daeSBarry Smith 75834e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[]) 759d71ae5a4SJacob Faibussowitsch { 7605c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7615c6c1daeSBarry Smith 7625c6c1daeSBarry Smith PetscFunctionBegin; 7635c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 7644f572ea9SToby Isaac PetscAssertPointer(name, 2); 765a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 7660298fd71SBarry Smith else *name = NULL; 7673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7685c6c1daeSBarry Smith } 7695c6c1daeSBarry Smith 7703014b61aSVaclav Hapla /*@C 771811af0c4SBarry Smith PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`, 772874270d9SVaclav Hapla and return this group's ID and file ID. 773811af0c4SBarry Smith If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID. 774874270d9SVaclav Hapla 77520f4b53cSBarry Smith Not Collective 776874270d9SVaclav Hapla 7773014b61aSVaclav Hapla Input Parameters: 7783014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7793014b61aSVaclav Hapla - path - (Optional) The path relative to the pushed group 780874270d9SVaclav Hapla 781d8d19677SJose E. Roman Output Parameters: 782874270d9SVaclav Hapla + fileId - The HDF5 file ID 783874270d9SVaclav Hapla - groupId - The HDF5 group ID 784874270d9SVaclav Hapla 7853f423023SBarry Smith Level: intermediate 7863f423023SBarry Smith 787811af0c4SBarry Smith Note: 7883014b61aSVaclav 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. 7893f423023SBarry Smith `NULL` or empty path means the current pushed group. 7903014b61aSVaclav Hapla 791e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 792e74616beSVaclav Hapla 793d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()` 794874270d9SVaclav Hapla @*/ 795d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId) 796d71ae5a4SJacob Faibussowitsch { 79790e03537SVaclav Hapla hid_t file_id; 79890e03537SVaclav Hapla H5O_type_t type; 7993014b61aSVaclav Hapla const char *fileName = NULL; 800ce78bad3SBarry Smith const char *groupName = NULL; 80176d59af2SVaclav Hapla PetscBool writable, has; 80254dbf706SVaclav Hapla 80354dbf706SVaclav Hapla PetscFunctionBegin; 8049566063dSJacob Faibussowitsch PetscCall(PetscViewerWritable(viewer, &writable)); 8059566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 8069566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetName(viewer, &fileName)); 8073014b61aSVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName)); 8089566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 80976d59af2SVaclav Hapla if (!has) { 8105f80ce2aSJacob 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); 811f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 81276d59af2SVaclav Hapla } 8135f80ce2aSJacob 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); 8143014b61aSVaclav Hapla PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT)); 8153014b61aSVaclav Hapla PetscCall(PetscFree(groupName)); 81654dbf706SVaclav Hapla *fileId = file_id; 8173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 81854dbf706SVaclav Hapla } 81954dbf706SVaclav Hapla 82063cb69f5SVaclav Hapla /*@C 82163cb69f5SVaclav Hapla PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file 82263cb69f5SVaclav Hapla 82320f4b53cSBarry Smith Not Collective 82463cb69f5SVaclav Hapla 82563cb69f5SVaclav Hapla Input Parameters: 82663cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 82763cb69f5SVaclav Hapla - path - (Optional) The path relative to the pushed group 82863cb69f5SVaclav Hapla 8293f423023SBarry Smith Level: intermediate 8303f423023SBarry Smith 83163cb69f5SVaclav Hapla Note: 83263cb69f5SVaclav 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. 8333f423023SBarry Smith `NULL` or empty path means the current pushed group. 83463cb69f5SVaclav Hapla 83563cb69f5SVaclav Hapla This will fail if the viewer is not writable. 83663cb69f5SVaclav Hapla 837d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 83863cb69f5SVaclav Hapla @*/ 839d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[]) 840d71ae5a4SJacob Faibussowitsch { 84163cb69f5SVaclav Hapla hid_t fileId, groupId; 84263cb69f5SVaclav Hapla 84363cb69f5SVaclav Hapla PetscFunctionBegin; 84463cb69f5SVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created 84563cb69f5SVaclav Hapla PetscCallHDF5(H5Gclose, (groupId)); 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84763cb69f5SVaclav Hapla } 84863cb69f5SVaclav Hapla 8493ef9c667SSatish Balay /*@ 850d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8515c6c1daeSBarry Smith 85220f4b53cSBarry Smith Not Collective 8535c6c1daeSBarry Smith 8545c6c1daeSBarry Smith Input Parameter: 855811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 8565c6c1daeSBarry Smith 8575c6c1daeSBarry Smith Level: intermediate 8585c6c1daeSBarry Smith 859d7dd068bSVaclav Hapla Notes: 860811af0c4SBarry Smith On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0. 861811af0c4SBarry Smith Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`. 862811af0c4SBarry 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()`]. 863811af0c4SBarry Smith Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 864811af0c4SBarry Smith Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`. 865d7dd068bSVaclav Hapla 866d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 867d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 868d7dd068bSVaclav Hapla 869aec76313SJacob Faibussowitsch Developer Notes: 870d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 871d7dd068bSVaclav Hapla 872d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 873d7dd068bSVaclav Hapla @*/ 874d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 875d71ae5a4SJacob Faibussowitsch { 876d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 877d7dd068bSVaclav Hapla 878d7dd068bSVaclav Hapla PetscFunctionBegin; 879d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 8805f80ce2aSJacob Faibussowitsch PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 881d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 882d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 8833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 884d7dd068bSVaclav Hapla } 885d7dd068bSVaclav Hapla 886d7dd068bSVaclav Hapla /*@ 887d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 888d7dd068bSVaclav Hapla 88920f4b53cSBarry Smith Not Collective 890d7dd068bSVaclav Hapla 891d7dd068bSVaclav Hapla Input Parameter: 892811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 893d7dd068bSVaclav Hapla 894d7dd068bSVaclav Hapla Level: intermediate 895d7dd068bSVaclav Hapla 896811af0c4SBarry Smith Note: 897811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 898d7dd068bSVaclav Hapla 899d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 900d7dd068bSVaclav Hapla @*/ 901d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 902d71ae5a4SJacob Faibussowitsch { 903d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 904d7dd068bSVaclav Hapla 905d7dd068bSVaclav Hapla PetscFunctionBegin; 906d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 9075f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 908d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 9093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 910d7dd068bSVaclav Hapla } 911d7dd068bSVaclav Hapla 912d7dd068bSVaclav Hapla /*@ 913d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 914d7dd068bSVaclav Hapla 91520f4b53cSBarry Smith Not Collective 916d7dd068bSVaclav Hapla 917d7dd068bSVaclav Hapla Input Parameter: 918811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 919d7dd068bSVaclav Hapla 920d7dd068bSVaclav Hapla Output Parameter: 921d7dd068bSVaclav Hapla . flg - is timestepping active? 922d7dd068bSVaclav Hapla 923d7dd068bSVaclav Hapla Level: intermediate 924d7dd068bSVaclav Hapla 925811af0c4SBarry Smith Note: 926811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 927d7dd068bSVaclav Hapla 928d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 929d7dd068bSVaclav Hapla @*/ 930d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 931d71ae5a4SJacob Faibussowitsch { 932d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 933d7dd068bSVaclav Hapla 934d7dd068bSVaclav Hapla PetscFunctionBegin; 935d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 936d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 9373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 938d7dd068bSVaclav Hapla } 939d7dd068bSVaclav Hapla 940d7dd068bSVaclav Hapla /*@ 941d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 942d7dd068bSVaclav Hapla 94320f4b53cSBarry Smith Not Collective 944d7dd068bSVaclav Hapla 945d7dd068bSVaclav Hapla Input Parameter: 946811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 947d7dd068bSVaclav Hapla 948d7dd068bSVaclav Hapla Level: intermediate 949d7dd068bSVaclav Hapla 950811af0c4SBarry Smith Note: 951811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 952d7dd068bSVaclav Hapla 953d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 9545c6c1daeSBarry Smith @*/ 955d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 956d71ae5a4SJacob Faibussowitsch { 9575c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9585c6c1daeSBarry Smith 9595c6c1daeSBarry Smith PetscFunctionBegin; 9605c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 9615f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9625c6c1daeSBarry Smith ++hdf5->timestep; 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9645c6c1daeSBarry Smith } 9655c6c1daeSBarry Smith 9663ef9c667SSatish Balay /*@ 967d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9685c6c1daeSBarry Smith 96920f4b53cSBarry Smith Logically Collective 9705c6c1daeSBarry Smith 9715c6c1daeSBarry Smith Input Parameters: 972811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 973d7dd068bSVaclav Hapla - timestep - The timestep 9745c6c1daeSBarry Smith 9755c6c1daeSBarry Smith Level: intermediate 9765c6c1daeSBarry Smith 977811af0c4SBarry Smith Note: 978811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 979d7dd068bSVaclav Hapla 980d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 9815c6c1daeSBarry Smith @*/ 982d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 983d71ae5a4SJacob Faibussowitsch { 9845c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9855c6c1daeSBarry Smith 9865c6c1daeSBarry Smith PetscFunctionBegin; 9875c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 988d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 9895f80ce2aSJacob Faibussowitsch PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 9905f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9915c6c1daeSBarry Smith hdf5->timestep = timestep; 9923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9935c6c1daeSBarry Smith } 9945c6c1daeSBarry Smith 9953ef9c667SSatish Balay /*@ 9965c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 9975c6c1daeSBarry Smith 99820f4b53cSBarry Smith Not Collective 9995c6c1daeSBarry Smith 10005c6c1daeSBarry Smith Input Parameter: 1001811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 10025c6c1daeSBarry Smith 10035c6c1daeSBarry Smith Output Parameter: 1004d7dd068bSVaclav Hapla . timestep - The timestep 10055c6c1daeSBarry Smith 10065c6c1daeSBarry Smith Level: intermediate 10075c6c1daeSBarry Smith 1008811af0c4SBarry Smith Note: 1009811af0c4SBarry Smith This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 1010d7dd068bSVaclav Hapla 1011d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 10125c6c1daeSBarry Smith @*/ 1013d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 1014d71ae5a4SJacob Faibussowitsch { 10155c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 10165c6c1daeSBarry Smith 10175c6c1daeSBarry Smith PetscFunctionBegin; 10185c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 10194f572ea9SToby Isaac PetscAssertPointer(timestep, 2); 10205f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 10215c6c1daeSBarry Smith *timestep = hdf5->timestep; 10223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10235c6c1daeSBarry Smith } 10245c6c1daeSBarry Smith 102536bce6e8SMatthew G. Knepley /*@C 102636bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 102736bce6e8SMatthew G. Knepley 102820f4b53cSBarry Smith Not Collective 102936bce6e8SMatthew G. Knepley 103036bce6e8SMatthew G. Knepley Input Parameter: 1031811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 103236bce6e8SMatthew G. Knepley 103336bce6e8SMatthew G. Knepley Output Parameter: 1034aec76313SJacob Faibussowitsch . htype - the HDF5 datatype 103536bce6e8SMatthew G. Knepley 103636bce6e8SMatthew G. Knepley Level: advanced 103736bce6e8SMatthew G. Knepley 1038d1f92df0SBarry Smith .seealso: [](sec_viewers), `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 103936bce6e8SMatthew G. Knepley @*/ 1040d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 1041d71ae5a4SJacob Faibussowitsch { 104236bce6e8SMatthew G. Knepley PetscFunctionBegin; 104336bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 104436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 104536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 104636bce6e8SMatthew G. Knepley #else 104736bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 104836bce6e8SMatthew G. Knepley #endif 104936bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 105036bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 105136bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1052c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1053de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 105436bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 105536bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 105636bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10577e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 105836bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 10593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 106036bce6e8SMatthew G. Knepley } 106136bce6e8SMatthew G. Knepley 106236bce6e8SMatthew G. Knepley /*@C 106336bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 106436bce6e8SMatthew G. Knepley 106520f4b53cSBarry Smith Not Collective 106636bce6e8SMatthew G. Knepley 106736bce6e8SMatthew G. Knepley Input Parameter: 1068811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...) 106936bce6e8SMatthew G. Knepley 107036bce6e8SMatthew G. Knepley Output Parameter: 1071811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 107236bce6e8SMatthew G. Knepley 107336bce6e8SMatthew G. Knepley Level: advanced 107436bce6e8SMatthew G. Knepley 107542747ad1SJacob Faibussowitsch .seealso: [](sec_viewers), `PetscDataType` 107636bce6e8SMatthew G. Knepley @*/ 1077d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 1078d71ae5a4SJacob Faibussowitsch { 107936bce6e8SMatthew G. Knepley PetscFunctionBegin; 108036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 108136bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 108236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 108336bce6e8SMatthew G. Knepley #else 108436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 108536bce6e8SMatthew G. Knepley #endif 108636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 108736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 108836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 108936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 109036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 109136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 10927e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 109336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 10943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 109536bce6e8SMatthew G. Knepley } 109636bce6e8SMatthew G. Knepley 1097df863907SAlex Fikl /*@C 1098b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 109936bce6e8SMatthew G. Knepley 1100343a20b2SVaclav Hapla Collective 1101343a20b2SVaclav Hapla 110236bce6e8SMatthew G. Knepley Input Parameters: 1103811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 11044302210dSVaclav Hapla . parent - The parent dataset/group name 110536bce6e8SMatthew G. Knepley . name - The attribute name 110636bce6e8SMatthew G. Knepley . datatype - The attribute type 110736bce6e8SMatthew G. Knepley - value - The attribute value 110836bce6e8SMatthew G. Knepley 110936bce6e8SMatthew G. Knepley Level: advanced 111036bce6e8SMatthew G. Knepley 1111811af0c4SBarry Smith Note: 1112343a20b2SVaclav 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. 11134302210dSVaclav Hapla 1114d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, 1115811af0c4SBarry Smith `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 111636bce6e8SMatthew G. Knepley @*/ 1117d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 1118d71ae5a4SJacob Faibussowitsch { 1119ce78bad3SBarry Smith const char *parentAbsPath; 112060568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1121080f144cSVaclav Hapla PetscBool has; 112236bce6e8SMatthew G. Knepley 112336bce6e8SMatthew G. Knepley PetscFunctionBegin; 11245cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 11254f572ea9SToby Isaac if (parent) PetscAssertPointer(parent, 2); 11264f572ea9SToby Isaac PetscAssertPointer(name, 3); 11274302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer, datatype, 4); 11284f572ea9SToby Isaac PetscAssertPointer(value, 5); 112977717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 11309566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 11319566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 11329566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 11337e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 11347e97c476SMatthew G. Knepley size_t len; 11359566063dSJacob Faibussowitsch PetscCall(PetscStrlen((const char *)value, &len)); 1136792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 11377e97c476SMatthew G. Knepley } 11389566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1139792fecdfSBarry Smith PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR)); 1140792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1141080f144cSVaclav Hapla if (has) { 1142792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1143080f144cSVaclav Hapla } else { 1144792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1145080f144cSVaclav Hapla } 1146792fecdfSBarry Smith PetscCallHDF5(H5Awrite, (attribute, dtype, value)); 1147792fecdfSBarry Smith if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype)); 1148792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1149792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 1150792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (dataspace)); 11519566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 11523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 115336bce6e8SMatthew G. Knepley } 115436bce6e8SMatthew G. Knepley 1155df863907SAlex Fikl /*@C 1156811af0c4SBarry Smith PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name 1157577e0e04SVaclav Hapla 1158343a20b2SVaclav Hapla Collective 1159343a20b2SVaclav Hapla 1160577e0e04SVaclav Hapla Input Parameters: 1161811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1162577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1163577e0e04SVaclav Hapla . name - The attribute name 1164577e0e04SVaclav Hapla . datatype - The attribute type 1165577e0e04SVaclav Hapla - value - The attribute value 1166577e0e04SVaclav Hapla 11673f423023SBarry Smith Level: advanced 11683f423023SBarry Smith 1169811af0c4SBarry Smith Note: 1170343a20b2SVaclav 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). 1171811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1172577e0e04SVaclav Hapla 1173d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, 1174811af0c4SBarry Smith `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1175577e0e04SVaclav Hapla @*/ 1176d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1177d71ae5a4SJacob Faibussowitsch { 1178577e0e04SVaclav Hapla PetscFunctionBegin; 1179577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1180577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 11814f572ea9SToby Isaac PetscAssertPointer(name, 3); 11824f572ea9SToby Isaac PetscAssertPointer(value, 5); 11839566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 11849566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 11853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1186577e0e04SVaclav Hapla } 1187577e0e04SVaclav Hapla 1188577e0e04SVaclav Hapla /*@C 1189b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1190ad0c61feSMatthew G. Knepley 1191343a20b2SVaclav Hapla Collective 1192343a20b2SVaclav Hapla 1193ad0c61feSMatthew G. Knepley Input Parameters: 1194811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 11954302210dSVaclav Hapla . parent - The parent dataset/group name 1196ad0c61feSMatthew G. Knepley . name - The attribute name 1197a2d6be1bSVaclav Hapla . datatype - The attribute type 1198a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1199ad0c61feSMatthew G. Knepley 1200ad0c61feSMatthew G. Knepley Output Parameter: 1201a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1202ad0c61feSMatthew G. Knepley 12033f423023SBarry Smith Level: advanced 12043f423023SBarry Smith 1205a2d6be1bSVaclav Hapla Notes: 12063f423023SBarry Smith If defaultValue is `NULL` and the attribute is not found, an error occurs. 12073f423023SBarry Smith 12083f423023SBarry Smith If defaultValue is not `NULL` and the attribute is not found, `defaultValue` is copied to value. 12093f423023SBarry Smith 12103f423023SBarry Smith The pointers `defaultValue` and value can be the same; for instance 12113f423023SBarry Smith .vb 12123f423023SBarry Smith flg = PETSC_FALSE; 12133f423023SBarry Smith PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 12143f423023SBarry Smith .ve 1215a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1216a2d6be1bSVaclav Hapla 1217811af0c4SBarry Smith If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed. 1218708d4cb3SBarry Smith 12193f423023SBarry 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. 12204302210dSVaclav Hapla 1221d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1222ad0c61feSMatthew G. Knepley @*/ 1223d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1224d71ae5a4SJacob Faibussowitsch { 1225ce78bad3SBarry Smith const char *parentAbsPath; 1226a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1227969748fdSVaclav Hapla PetscBool has; 1228ad0c61feSMatthew G. Knepley 1229ad0c61feSMatthew G. Knepley PetscFunctionBegin; 12305cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 12314f572ea9SToby Isaac if (parent) PetscAssertPointer(parent, 2); 12324f572ea9SToby Isaac PetscAssertPointer(name, 3); 12334f572ea9SToby Isaac if (defaultValue) PetscAssertPointer(defaultValue, 5); 12344f572ea9SToby Isaac PetscAssertPointer(value, 6); 12359566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 123677717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 12379566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 12389566063dSJacob Faibussowitsch if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1239a2d6be1bSVaclav Hapla if (!has) { 1240a2d6be1bSVaclav Hapla if (defaultValue) { 1241a2d6be1bSVaclav Hapla if (defaultValue != value) { 1242a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 12439566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value)); 1244a2d6be1bSVaclav Hapla } else { 1245a2d6be1bSVaclav Hapla size_t len; 1246792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype)); 12479566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(value, defaultValue, len)); 1248a2d6be1bSVaclav Hapla } 1249a2d6be1bSVaclav Hapla } 12509566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 12513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125298921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1253a2d6be1bSVaclav Hapla } 12549566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1255792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1256792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1257f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1258f0b58479SMatthew G. Knepley size_t len; 1259a2d6be1bSVaclav Hapla hid_t atype; 1260792fecdfSBarry Smith PetscCallHDF5Return(atype, H5Aget_type, (attribute)); 1261792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype)); 12629566063dSJacob Faibussowitsch PetscCall(PetscMalloc((len + 1) * sizeof(char), value)); 1263792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 1264792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value)); 1265708d4cb3SBarry Smith } else { 1266792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, value)); 1267708d4cb3SBarry Smith } 1268792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1269e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1270792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 12719566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 12723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1273ad0c61feSMatthew G. Knepley } 1274ad0c61feSMatthew G. Knepley 1275577e0e04SVaclav Hapla /*@C 1276811af0c4SBarry Smith PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name 1277577e0e04SVaclav Hapla 1278343a20b2SVaclav Hapla Collective 1279343a20b2SVaclav Hapla 1280577e0e04SVaclav Hapla Input Parameters: 1281811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1282577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1283577e0e04SVaclav Hapla . name - The attribute name 128410450e9eSJacob Faibussowitsch . datatype - The attribute type 128510450e9eSJacob Faibussowitsch - defaultValue - The default attribute value 1286577e0e04SVaclav Hapla 1287577e0e04SVaclav Hapla Output Parameter: 1288577e0e04SVaclav Hapla . value - The attribute value 1289577e0e04SVaclav Hapla 12903f423023SBarry Smith Level: advanced 12913f423023SBarry Smith 1292811af0c4SBarry Smith Note: 1293577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1294811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1295577e0e04SVaclav Hapla 1296d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1297577e0e04SVaclav Hapla @*/ 1298d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1299d71ae5a4SJacob Faibussowitsch { 1300577e0e04SVaclav Hapla PetscFunctionBegin; 1301577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1302577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 13034f572ea9SToby Isaac PetscAssertPointer(name, 3); 13044f572ea9SToby Isaac PetscAssertPointer(value, 6); 13059566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 13069566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 13073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1308577e0e04SVaclav Hapla } 1309577e0e04SVaclav Hapla 131034e79e72SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1311d71ae5a4SJacob Faibussowitsch { 1312a75c3fd4SVaclav Hapla htri_t exists; 1313a75c3fd4SVaclav Hapla hid_t group; 1314a75c3fd4SVaclav Hapla 1315a75c3fd4SVaclav Hapla PetscFunctionBegin; 1316792fecdfSBarry Smith PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT)); 1317792fecdfSBarry Smith if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT)); 1318a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1319792fecdfSBarry Smith PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1320792fecdfSBarry Smith PetscCallHDF5(H5Gclose, (group)); 1321a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1322a75c3fd4SVaclav Hapla } 1323a75c3fd4SVaclav Hapla *exists_ = (PetscBool)exists; 13243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1325a75c3fd4SVaclav Hapla } 1326a75c3fd4SVaclav Hapla 1327d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1328d71ae5a4SJacob Faibussowitsch { 132990e03537SVaclav Hapla const char rootGroupName[] = "/"; 13305cdeb108SMatthew G. Knepley hid_t h5; 1331e5bf9ebcSVaclav Hapla PetscBool exists = PETSC_FALSE; 133215b861d2SVaclav Hapla PetscInt i; 133315b861d2SVaclav Hapla int n; 133485688ae2SVaclav Hapla char **hierarchy; 133585688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 13365cdeb108SMatthew G. Knepley 13375cdeb108SMatthew G. Knepley PetscFunctionBegin; 13385cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 13394f572ea9SToby Isaac if (name) PetscAssertPointer(name, 2); 134090e03537SVaclav Hapla else name = rootGroupName; 1341ccd66a83SVaclav Hapla if (has) { 13424f572ea9SToby Isaac PetscAssertPointer(has, 4); 13435cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1344ccd66a83SVaclav Hapla } 1345ccd66a83SVaclav Hapla if (otype) { 13464f572ea9SToby Isaac PetscAssertPointer(otype, 5); 134756cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1348ccd66a83SVaclav Hapla } 13499566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 135085688ae2SVaclav Hapla 135185688ae2SVaclav Hapla /* 135285688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 135385688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 135485688ae2SVaclav Hapla 1) whether it's a valid link 135585688ae2SVaclav Hapla 2) whether this link resolves to an object 135685688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 135785688ae2SVaclav Hapla */ 13589566063dSJacob Faibussowitsch PetscCall(PetscStrToArray(name, '/', &n, &hierarchy)); 135985688ae2SVaclav Hapla if (!n) { 136085688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1361ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1362ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 13639566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 13643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136585688ae2SVaclav Hapla } 136685688ae2SVaclav Hapla for (i = 0; i < n; i++) { 1367c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, "/", sizeof(buf))); 1368c6a7a370SJeremy L Thompson PetscCall(PetscStrlcat(buf, hierarchy[i], sizeof(buf))); 13699566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1370a75c3fd4SVaclav Hapla if (!exists) break; 137185688ae2SVaclav Hapla } 13729566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 137385688ae2SVaclav Hapla 137485688ae2SVaclav Hapla /* If the object exists, get its type */ 1375ccd66a83SVaclav Hapla if (exists && otype) { 13765cdeb108SMatthew G. Knepley H5O_info_t info; 13775cdeb108SMatthew G. Knepley 137876276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1379792fecdfSBarry Smith PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT)); 138056cc0592SVaclav Hapla *otype = info.type; 13815cdeb108SMatthew G. Knepley } 1382ccd66a83SVaclav Hapla if (has) *has = exists; 13833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13845cdeb108SMatthew G. Knepley } 13855cdeb108SMatthew G. Knepley 13864302210dSVaclav Hapla /*@C 13870a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13880a9f38efSVaclav Hapla 1389343a20b2SVaclav Hapla Collective 1390343a20b2SVaclav Hapla 13910a9f38efSVaclav Hapla Input Parameters: 1392811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1393a0558868SVaclav Hapla - path - (Optional) The path relative to the pushed group 13940a9f38efSVaclav Hapla 13950a9f38efSVaclav Hapla Output Parameter: 13960a9f38efSVaclav Hapla . has - Flag for group existence 13970a9f38efSVaclav Hapla 13980a9f38efSVaclav Hapla Level: advanced 13990a9f38efSVaclav Hapla 14004302210dSVaclav Hapla Notes: 1401785c443eSVaclav 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. 14023f423023SBarry Smith `NULL` or empty path means the current pushed group. 14034302210dSVaclav Hapla 1404811af0c4SBarry Smith If path exists but is not a group, `PETSC_FALSE` is returned. 1405811af0c4SBarry Smith 1406d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 14070a9f38efSVaclav Hapla @*/ 1408d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 1409d71ae5a4SJacob Faibussowitsch { 14100a9f38efSVaclav Hapla H5O_type_t type; 1411ce78bad3SBarry Smith const char *abspath; 14120a9f38efSVaclav Hapla 14130a9f38efSVaclav Hapla PetscFunctionBegin; 14140a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 14154f572ea9SToby Isaac if (path) PetscAssertPointer(path, 2); 14164f572ea9SToby Isaac PetscAssertPointer(has, 3); 141777717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 14189566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 14194302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 14209566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 14213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 14220a9f38efSVaclav Hapla } 14230a9f38efSVaclav Hapla 142489e0ef10SVaclav Hapla /*@C 142589e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 142689e0ef10SVaclav Hapla 1427343a20b2SVaclav Hapla Collective 1428343a20b2SVaclav Hapla 142989e0ef10SVaclav Hapla Input Parameters: 1430811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 143189e0ef10SVaclav Hapla - path - The dataset path 143289e0ef10SVaclav Hapla 143389e0ef10SVaclav Hapla Output Parameter: 143489e0ef10SVaclav Hapla . has - Flag whether dataset exists 143589e0ef10SVaclav Hapla 143689e0ef10SVaclav Hapla Level: advanced 143789e0ef10SVaclav Hapla 143889e0ef10SVaclav Hapla Notes: 1439343a20b2SVaclav 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. 144089e0ef10SVaclav Hapla 14413f423023SBarry Smith If `path` is `NULL` or empty, has is set to `PETSC_FALSE`. 1442811af0c4SBarry Smith 14433f423023SBarry Smith If `path` exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1444811af0c4SBarry Smith 1445d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 144689e0ef10SVaclav Hapla @*/ 1447d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 1448d71ae5a4SJacob Faibussowitsch { 144989e0ef10SVaclav Hapla H5O_type_t type; 1450ce78bad3SBarry Smith const char *abspath; 145189e0ef10SVaclav Hapla 145289e0ef10SVaclav Hapla PetscFunctionBegin; 145389e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 14544f572ea9SToby Isaac if (path) PetscAssertPointer(path, 2); 14554f572ea9SToby Isaac PetscAssertPointer(has, 3); 145677717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 14579566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 145889e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 14599566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 14603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 146189e0ef10SVaclav Hapla } 146289e0ef10SVaclav Hapla 14630a9f38efSVaclav Hapla /*@ 1464e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1465ecce1506SVaclav Hapla 1466343a20b2SVaclav Hapla Collective 1467343a20b2SVaclav Hapla 1468ecce1506SVaclav Hapla Input Parameters: 1469811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1470ecce1506SVaclav Hapla - obj - The named object 1471ecce1506SVaclav Hapla 1472ecce1506SVaclav Hapla Output Parameter: 147389e0ef10SVaclav Hapla . has - Flag for dataset existence 1474ecce1506SVaclav Hapla 14753f423023SBarry Smith Level: advanced 14763f423023SBarry Smith 1477e3f143f7SVaclav Hapla Notes: 147889e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1479811af0c4SBarry Smith 1480811af0c4SBarry Smith If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1481e3f143f7SVaclav Hapla 1482d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1483ecce1506SVaclav Hapla @*/ 1484d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1485d71ae5a4SJacob Faibussowitsch { 148689e0ef10SVaclav Hapla size_t len; 1487ecce1506SVaclav Hapla 1488ecce1506SVaclav Hapla PetscFunctionBegin; 1489c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1490c57a1a9aSVaclav Hapla PetscValidHeader(obj, 2); 14914f572ea9SToby Isaac PetscAssertPointer(has, 3); 14929566063dSJacob Faibussowitsch PetscCall(PetscStrlen(obj->name, &len)); 14935f80ce2aSJacob Faibussowitsch PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 14949566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 14953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1496ecce1506SVaclav Hapla } 1497ecce1506SVaclav Hapla 1498df863907SAlex Fikl /*@C 1499ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1500e7bdbf8eSMatthew G. Knepley 1501343a20b2SVaclav Hapla Collective 1502343a20b2SVaclav Hapla 1503e7bdbf8eSMatthew G. Knepley Input Parameters: 1504811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 15054302210dSVaclav Hapla . parent - The parent dataset/group name 1506e7bdbf8eSMatthew G. Knepley - name - The attribute name 1507e7bdbf8eSMatthew G. Knepley 1508e7bdbf8eSMatthew G. Knepley Output Parameter: 1509e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1510e7bdbf8eSMatthew G. Knepley 1511e7bdbf8eSMatthew G. Knepley Level: advanced 1512e7bdbf8eSMatthew G. Knepley 1513811af0c4SBarry Smith Note: 15143f423023SBarry 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. 15154302210dSVaclav Hapla 1516d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1517e7bdbf8eSMatthew G. Knepley @*/ 1518d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1519d71ae5a4SJacob Faibussowitsch { 1520ce78bad3SBarry Smith const char *parentAbsPath; 1521e7bdbf8eSMatthew G. Knepley 1522e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 15235cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 15244f572ea9SToby Isaac if (parent) PetscAssertPointer(parent, 2); 15254f572ea9SToby Isaac PetscAssertPointer(name, 3); 15264f572ea9SToby Isaac PetscAssertPointer(has, 4); 152777717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 15289566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 15299566063dSJacob Faibussowitsch if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 15309566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 15313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 153206db490cSVaclav Hapla } 153306db490cSVaclav Hapla 1534577e0e04SVaclav Hapla /*@C 1535811af0c4SBarry Smith PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name 1536577e0e04SVaclav Hapla 1537343a20b2SVaclav Hapla Collective 1538343a20b2SVaclav Hapla 1539577e0e04SVaclav Hapla Input Parameters: 1540811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1541577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1542577e0e04SVaclav Hapla - name - The attribute name 1543577e0e04SVaclav Hapla 1544577e0e04SVaclav Hapla Output Parameter: 1545577e0e04SVaclav Hapla . has - Flag for attribute existence 1546577e0e04SVaclav Hapla 15473f423023SBarry Smith Level: advanced 15483f423023SBarry Smith 1549811af0c4SBarry Smith Note: 1550577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1551811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1552577e0e04SVaclav Hapla 1553d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1554577e0e04SVaclav Hapla @*/ 1555d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1556d71ae5a4SJacob Faibussowitsch { 1557577e0e04SVaclav Hapla PetscFunctionBegin; 1558577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1559577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 15604f572ea9SToby Isaac PetscAssertPointer(name, 3); 15614f572ea9SToby Isaac PetscAssertPointer(has, 4); 15629566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 15639566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 15643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1565577e0e04SVaclav Hapla } 1566577e0e04SVaclav Hapla 1567d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1568d71ae5a4SJacob Faibussowitsch { 156906db490cSVaclav Hapla hid_t h5; 157006db490cSVaclav Hapla htri_t hhas; 157106db490cSVaclav Hapla 157206db490cSVaclav Hapla PetscFunctionBegin; 15739566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1574792fecdfSBarry Smith PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT)); 15755cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 15763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1577e7bdbf8eSMatthew G. Knepley } 1578e7bdbf8eSMatthew G. Knepley 1579a75e6a4aSMatthew G. Knepley /* 1580a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1581a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1582a75e6a4aSMatthew G. Knepley */ 1583d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1584a75e6a4aSMatthew G. Knepley 1585a75e6a4aSMatthew G. Knepley /*@C 1586811af0c4SBarry Smith PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator. 1587a75e6a4aSMatthew G. Knepley 1588d083f849SBarry Smith Collective 1589a75e6a4aSMatthew G. Knepley 1590a75e6a4aSMatthew G. Knepley Input Parameter: 15913f423023SBarry Smith . comm - the MPI communicator to share the `PETSCVIEWERHDF5` `PetscViewer` 1592a75e6a4aSMatthew G. Knepley 1593811af0c4SBarry Smith Options Database Key: 159410699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1595a75e6a4aSMatthew G. Knepley 1596811af0c4SBarry Smith Environmental variable: 1597811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file 1598a75e6a4aSMatthew G. Knepley 159920f4b53cSBarry Smith Level: intermediate 160020f4b53cSBarry Smith 1601811af0c4SBarry Smith Note: 1602811af0c4SBarry Smith Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return 1603811af0c4SBarry Smith an error code. The HDF5 `PetscViewer` is usually used in the form 1604*b44f4de4SBarry Smith .vb 1605*b44f4de4SBarry Smith XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1606*b44f4de4SBarry Smith .ve 1607a75e6a4aSMatthew G. Knepley 1608d1f92df0SBarry Smith .seealso: [](sec_viewers), `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1609a75e6a4aSMatthew G. Knepley @*/ 1610d71ae5a4SJacob Faibussowitsch PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1611d71ae5a4SJacob Faibussowitsch { 1612a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 16133ba16761SJacob Faibussowitsch PetscMPIInt mpi_ierr; 1614a75e6a4aSMatthew G. Knepley PetscBool flg; 1615a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1616a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1617a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1618a75e6a4aSMatthew G. Knepley 1619a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 16209371c9d4SSatish Balay ierr = PetscCommDuplicate(comm, &ncomm, NULL); 16219371c9d4SSatish Balay if (ierr) { 16223ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16239371c9d4SSatish Balay PetscFunctionReturn(NULL); 16249371c9d4SSatish Balay } 1625a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 16263ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL); 16273ba16761SJacob Faibussowitsch if (mpi_ierr) { 16283ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16299371c9d4SSatish Balay PetscFunctionReturn(NULL); 16309371c9d4SSatish Balay } 1631a75e6a4aSMatthew G. Knepley } 16323ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg); 16333ba16761SJacob Faibussowitsch if (mpi_ierr) { 16343ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16359371c9d4SSatish Balay PetscFunctionReturn(NULL); 16369371c9d4SSatish Balay } 1637a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1638a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg); 16399371c9d4SSatish Balay if (ierr) { 16403ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16419371c9d4SSatish Balay PetscFunctionReturn(NULL); 16429371c9d4SSatish Balay } 1643a75e6a4aSMatthew G. Knepley if (!flg) { 1644c6a7a370SJeremy L Thompson ierr = PetscStrncpy(fname, "output.h5", sizeof(fname)); 16459371c9d4SSatish Balay if (ierr) { 16463ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16479371c9d4SSatish Balay PetscFunctionReturn(NULL); 16489371c9d4SSatish Balay } 1649a75e6a4aSMatthew G. Knepley } 1650a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer); 16519371c9d4SSatish Balay if (ierr) { 16523ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16539371c9d4SSatish Balay PetscFunctionReturn(NULL); 16549371c9d4SSatish Balay } 1655a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16569371c9d4SSatish Balay if (ierr) { 16573ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16589371c9d4SSatish Balay PetscFunctionReturn(NULL); 16599371c9d4SSatish Balay } 16603ba16761SJacob Faibussowitsch mpi_ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer); 16613ba16761SJacob Faibussowitsch if (mpi_ierr) { 16623ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 16639371c9d4SSatish Balay PetscFunctionReturn(NULL); 16649371c9d4SSatish Balay } 1665a75e6a4aSMatthew G. Knepley } 1666a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 16679371c9d4SSatish Balay if (ierr) { 16683ba16761SJacob Faibussowitsch ierr = PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 16699371c9d4SSatish Balay PetscFunctionReturn(NULL); 16709371c9d4SSatish Balay } 1671a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1672a75e6a4aSMatthew G. Knepley } 1673