16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 25c6c1daeSBarry Smith 3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *); 406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *); 506db490cSVaclav Hapla 677717648SVaclav Hapla /*@C 777717648SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`. 877717648SVaclav Hapla 977717648SVaclav Hapla Not collective 1077717648SVaclav Hapla 1177717648SVaclav Hapla Input Parameters: 1277717648SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 1377717648SVaclav Hapla - path - (Optional) The path relative to the pushed group 1477717648SVaclav Hapla 1577717648SVaclav Hapla Output Parameter: 1677717648SVaclav Hapla . abspath - The absolute HDF5 path (group) 1777717648SVaclav Hapla 1877717648SVaclav Hapla Level: intermediate 1977717648SVaclav Hapla 2077717648SVaclav Hapla Notes: 2177717648SVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 2277717648SVaclav Hapla So NULL or empty path means the current pushed group. 2377717648SVaclav Hapla 2477717648SVaclav Hapla The output abspath is newly allocated so needs to be freed. 2577717648SVaclav Hapla 26*63cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 2777717648SVaclav Hapla @*/ 2877717648SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[]) { 2977717648SVaclav Hapla size_t len; 304302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 316c132bc1SVaclav Hapla const char *group; 326c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 336c132bc1SVaclav Hapla 346c132bc1SVaclav Hapla PetscFunctionBegin; 3577717648SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 3677717648SVaclav Hapla if (path) PetscValidCharPointer(path, 2); 3777717648SVaclav Hapla PetscValidPointer(abspath, 3); 3877717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group)); 3977717648SVaclav Hapla PetscCall(PetscStrlen(path, &len)); 4077717648SVaclav Hapla relative = (PetscBool)(!len || path[0] != '/'); 414302210dSVaclav Hapla if (relative) { 429566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(buf, group)); 4377717648SVaclav Hapla if (!group || len) PetscCall(PetscStrcat(buf, "/")); 449566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf, path)); 459566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(buf, abspath)); 464302210dSVaclav Hapla } else { 479566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(path, abspath)); 484302210dSVaclav Hapla } 496c132bc1SVaclav Hapla PetscFunctionReturn(0); 506c132bc1SVaclav Hapla } 516c132bc1SVaclav Hapla 529371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) { 53577e0e04SVaclav Hapla PetscBool has; 54577e0e04SVaclav Hapla 55577e0e04SVaclav Hapla PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has)); 57577e0e04SVaclav Hapla if (!has) { 5877717648SVaclav Hapla char *group; 5977717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 6077717648SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group); 61577e0e04SVaclav Hapla } 62577e0e04SVaclav Hapla PetscFunctionReturn(0); 63577e0e04SVaclav Hapla } 64577e0e04SVaclav Hapla 659371c9d4SSatish Balay static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject) { 66ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 6782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 6882ea9b62SBarry Smith 6982ea9b62SBarry Smith PetscFunctionBegin; 70d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options"); 719566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL)); 729566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL)); 739566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set)); 749566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg)); 7519a20e4cSMatthew G. Knepley flg = PETSC_FALSE; 769566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set)); 779566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg)); 78d0609cedSBarry Smith PetscOptionsHeadEnd(); 7982ea9b62SBarry Smith PetscFunctionReturn(0); 8082ea9b62SBarry Smith } 8182ea9b62SBarry Smith 829371c9d4SSatish Balay static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer) { 831b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 8403fa8834SVaclav Hapla PetscBool flg; 851b793a25SVaclav Hapla 861b793a25SVaclav Hapla PetscFunctionBegin; 8748a46eb9SPierre Jolivet if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename)); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2])); 899566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput])); 909566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetCollective(v, &flg)); 919566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent")); 929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping])); 931b793a25SVaclav Hapla PetscFunctionReturn(0); 941b793a25SVaclav Hapla } 951b793a25SVaclav Hapla 969371c9d4SSatish Balay static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) { 975c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 985c6c1daeSBarry Smith 995c6c1daeSBarry Smith PetscFunctionBegin; 1009566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->filename)); 101792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 1025c6c1daeSBarry Smith PetscFunctionReturn(0); 1035c6c1daeSBarry Smith } 1045c6c1daeSBarry Smith 1059371c9d4SSatish Balay static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) { 1066226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1076226335aSVaclav Hapla 1086226335aSVaclav Hapla PetscFunctionBegin; 109792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL)); 1106226335aSVaclav Hapla PetscFunctionReturn(0); 1116226335aSVaclav Hapla } 1126226335aSVaclav Hapla 1139371c9d4SSatish Balay static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) { 1145c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1155c6c1daeSBarry Smith 1165c6c1daeSBarry Smith PetscFunctionBegin; 117792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (hdf5->dxpl_id)); 1189566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_HDF5(viewer)); 1195c6c1daeSBarry Smith while (hdf5->groups) { 1207d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1215c6c1daeSBarry Smith 1229566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups->name)); 1239566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups)); 1245c6c1daeSBarry Smith hdf5->groups = tmp; 1255c6c1daeSBarry Smith } 1269566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5)); 1279566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 1289566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 1299566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 1302e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 1319566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL)); 1329566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL)); 1332e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL)); 1342e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL)); 1359566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL)); 1369566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL)); 1375c6c1daeSBarry Smith PetscFunctionReturn(0); 1385c6c1daeSBarry Smith } 1395c6c1daeSBarry Smith 1409371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) { 1415c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1425c6c1daeSBarry Smith 1435c6c1daeSBarry Smith PetscFunctionBegin; 1445c6c1daeSBarry Smith hdf5->btype = type; 1455c6c1daeSBarry Smith PetscFunctionReturn(0); 1465c6c1daeSBarry Smith } 1475c6c1daeSBarry Smith 1489371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) { 1490b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 1500b62783dSVaclav Hapla 1510b62783dSVaclav Hapla PetscFunctionBegin; 1520b62783dSVaclav Hapla *type = hdf5->btype; 1530b62783dSVaclav Hapla PetscFunctionReturn(0); 1540b62783dSVaclav Hapla } 1550b62783dSVaclav Hapla 1569371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) { 15782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith PetscFunctionBegin; 16082ea9b62SBarry Smith hdf5->basedimension2 = flg; 16182ea9b62SBarry Smith PetscFunctionReturn(0); 16282ea9b62SBarry Smith } 16382ea9b62SBarry Smith 164df863907SAlex Fikl /*@ 16582ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 16682ea9b62SBarry Smith dimension of 2. 16782ea9b62SBarry Smith 168811af0c4SBarry Smith Logically Collective on viewer 16982ea9b62SBarry Smith 17082ea9b62SBarry Smith Input Parameters: 171811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored 172811af0c4SBarry 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 17382ea9b62SBarry Smith 174811af0c4SBarry Smith Options Database Key: 17582ea9b62SBarry 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 17682ea9b62SBarry Smith 177811af0c4SBarry Smith Note: 17895452b02SPatrick 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 17982ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18082ea9b62SBarry Smith 18182ea9b62SBarry Smith Level: intermediate 18282ea9b62SBarry Smith 183811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 18482ea9b62SBarry Smith @*/ 1859371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg) { 18682ea9b62SBarry Smith PetscFunctionBegin; 18782ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 188cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg)); 18982ea9b62SBarry Smith PetscFunctionReturn(0); 19082ea9b62SBarry Smith } 19182ea9b62SBarry Smith 192df863907SAlex Fikl /*@ 19382ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 19482ea9b62SBarry Smith dimension of 2. 19582ea9b62SBarry Smith 196811af0c4SBarry Smith Logically Collective on viewer 19782ea9b62SBarry Smith 19882ea9b62SBarry Smith Input Parameter: 199811af0c4SBarry Smith . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5` 20082ea9b62SBarry Smith 20182ea9b62SBarry Smith Output Parameter: 202811af0c4SBarry 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 20382ea9b62SBarry Smith 204811af0c4SBarry Smith Note: 20595452b02SPatrick 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 20682ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 20782ea9b62SBarry Smith 20882ea9b62SBarry Smith Level: intermediate 20982ea9b62SBarry Smith 210811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 21182ea9b62SBarry Smith @*/ 2129371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg) { 21382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 21482ea9b62SBarry Smith 21582ea9b62SBarry Smith PetscFunctionBegin; 21682ea9b62SBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 21782ea9b62SBarry Smith *flg = hdf5->basedimension2; 21882ea9b62SBarry Smith PetscFunctionReturn(0); 21982ea9b62SBarry Smith } 22082ea9b62SBarry Smith 2219371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) { 2229a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes PetscFunctionBegin; 2259a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2269a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2279a0502c6SHåkon Strandenes } 2289a0502c6SHåkon Strandenes 229df863907SAlex Fikl /*@ 2309a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 231811af0c4SBarry Smith compiled with double precision `PetscReal`. 2329a0502c6SHåkon Strandenes 233811af0c4SBarry Smith Logically Collective on viewer 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes Input Parameters: 236811af0c4SBarry Smith + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored 237811af0c4SBarry Smith - flg - if `PETSC_TRUE` the data will be written to disk with single precision 2389a0502c6SHåkon Strandenes 239811af0c4SBarry Smith Options Database Key: 2409a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2419a0502c6SHåkon Strandenes 242811af0c4SBarry Smith Note: 24395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2449a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2459a0502c6SHåkon Strandenes 2469a0502c6SHåkon Strandenes Level: intermediate 2479a0502c6SHåkon Strandenes 248811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 249811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5GetSPOutput()` 2509a0502c6SHåkon Strandenes @*/ 2519371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg) { 2529a0502c6SHåkon Strandenes PetscFunctionBegin; 2539a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 254cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg)); 2559a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2569a0502c6SHåkon Strandenes } 2579a0502c6SHåkon Strandenes 258df863907SAlex Fikl /*@ 2599a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 260811af0c4SBarry Smith compiled with double precision `PetscReal`. 2619a0502c6SHåkon Strandenes 262811af0c4SBarry Smith Logically Collective on viewer 2639a0502c6SHåkon Strandenes 2649a0502c6SHåkon Strandenes Input Parameter: 265811af0c4SBarry Smith . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5` 2669a0502c6SHåkon Strandenes 2679a0502c6SHåkon Strandenes Output Parameter: 268811af0c4SBarry Smith . flg - if `PETSC_TRUE` the data will be written to disk with single precision 2699a0502c6SHåkon Strandenes 2709a0502c6SHåkon Strandenes Level: intermediate 2719a0502c6SHåkon Strandenes 272db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 273811af0c4SBarry Smith `PetscReal`, `PetscViewerHDF5SetSPOutput()` 2749a0502c6SHåkon Strandenes @*/ 2759371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg) { 2769a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 2779a0502c6SHåkon Strandenes 2789a0502c6SHåkon Strandenes PetscFunctionBegin; 2799a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 2809a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2819a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2829a0502c6SHåkon Strandenes } 2839a0502c6SHåkon Strandenes 2849371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) { 285ee8b9732SVaclav Hapla PetscFunctionBegin; 286ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 287ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 288c796909eSBarry Smith #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL) 289ee8b9732SVaclav Hapla { 290ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 291792fecdfSBarry Smith PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 292ee8b9732SVaclav Hapla } 293ee8b9732SVaclav Hapla #else 2949566063dSJacob 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")); 295ee8b9732SVaclav Hapla #endif 296ee8b9732SVaclav Hapla PetscFunctionReturn(0); 297ee8b9732SVaclav Hapla } 298ee8b9732SVaclav Hapla 299ee8b9732SVaclav Hapla /*@ 300ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 301ee8b9732SVaclav Hapla 302ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 303ee8b9732SVaclav Hapla 304ee8b9732SVaclav Hapla Input Parameters: 305811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 306811af0c4SBarry Smith - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default) 307ee8b9732SVaclav Hapla 308811af0c4SBarry Smith Options Database Key: 309ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 310ee8b9732SVaclav Hapla 311811af0c4SBarry Smith Note: 312ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 31353effed7SBarry 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. 314ee8b9732SVaclav Hapla 315811af0c4SBarry Smith Developer Note: 316811af0c4SBarry Smith In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively. 317ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 318ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 319ee8b9732SVaclav Hapla 320ee8b9732SVaclav Hapla Level: intermediate 321ee8b9732SVaclav Hapla 322811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 323ee8b9732SVaclav Hapla @*/ 3249371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg) { 325ee8b9732SVaclav Hapla PetscFunctionBegin; 326ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 327ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer, flg, 2); 328cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg)); 329ee8b9732SVaclav Hapla PetscFunctionReturn(0); 330ee8b9732SVaclav Hapla } 331ee8b9732SVaclav Hapla 3329371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) { 333c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 334ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 335ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 336c796909eSBarry Smith #endif 337ee8b9732SVaclav Hapla 338ee8b9732SVaclav Hapla PetscFunctionBegin; 339c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 340c796909eSBarry Smith *flg = PETSC_FALSE; 341c796909eSBarry Smith #else 342792fecdfSBarry Smith PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode)); 343ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 344c796909eSBarry Smith #endif 345ee8b9732SVaclav Hapla PetscFunctionReturn(0); 346ee8b9732SVaclav Hapla } 347ee8b9732SVaclav Hapla 348ee8b9732SVaclav Hapla /*@ 349ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 350ee8b9732SVaclav Hapla 351ee8b9732SVaclav Hapla Not Collective 352ee8b9732SVaclav Hapla 353ee8b9732SVaclav Hapla Input Parameters: 354811af0c4SBarry Smith . viewer - the `PETSCVIEWERHDF5` `PetscViewer` 355ee8b9732SVaclav Hapla 356ee8b9732SVaclav Hapla Output Parameters: 357ee8b9732SVaclav Hapla . flg - the flag 358ee8b9732SVaclav Hapla 359ee8b9732SVaclav Hapla Level: intermediate 360ee8b9732SVaclav Hapla 361811af0c4SBarry Smith Note: 362811af0c4SBarry 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. 363811af0c4SBarry Smith For more details, see `PetscViewerHDF5SetCollective()`. 364ee8b9732SVaclav Hapla 365811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 366ee8b9732SVaclav Hapla @*/ 3679371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg) { 368ee8b9732SVaclav Hapla PetscFunctionBegin; 369ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 370534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 371ee8b9732SVaclav Hapla 372cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg)); 373ee8b9732SVaclav Hapla PetscFunctionReturn(0); 374ee8b9732SVaclav Hapla } 375ee8b9732SVaclav Hapla 3769371c9d4SSatish Balay static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) { 3775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 3785c6c1daeSBarry Smith hid_t plist_id; 3795c6c1daeSBarry Smith 3805c6c1daeSBarry Smith PetscFunctionBegin; 381792fecdfSBarry Smith if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 3829566063dSJacob Faibussowitsch if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 3839566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &hdf5->filename)); 3845c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 385792fecdfSBarry Smith PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS)); 386c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 387792fecdfSBarry Smith PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 388c796909eSBarry Smith #endif 3895c6c1daeSBarry Smith /* Create or open the file collectively */ 3905c6c1daeSBarry Smith switch (hdf5->btype) { 3915c6c1daeSBarry Smith case FILE_MODE_READ: 3928a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 3938a2871f6SBarry Smith PetscMPIInt rank; 3948a2871f6SBarry Smith PetscBool flg; 3958a2871f6SBarry Smith 3969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 3978a2871f6SBarry Smith if (rank == 0) { 3989566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 3995f80ce2aSJacob Faibussowitsch PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename); 4008a2871f6SBarry Smith } 4019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 4028a2871f6SBarry Smith } 403792fecdfSBarry Smith PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id)); 4045c6c1daeSBarry Smith break; 4055c6c1daeSBarry Smith case FILE_MODE_APPEND: 4069371c9d4SSatish Balay case FILE_MODE_UPDATE: { 4077e4fd573SVaclav Hapla PetscBool flg; 4089566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 409792fecdfSBarry Smith if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id)); 410792fecdfSBarry Smith else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4115c6c1daeSBarry Smith break; 4127e4fd573SVaclav Hapla } 4139371c9d4SSatish Balay case FILE_MODE_WRITE: PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); break; 4149371c9d4SSatish Balay case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4159371c9d4SSatish Balay default: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]); 4165c6c1daeSBarry Smith } 4175f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 418792fecdfSBarry Smith PetscCallHDF5(H5Pclose, (plist_id)); 4195c6c1daeSBarry Smith PetscFunctionReturn(0); 4205c6c1daeSBarry Smith } 4215c6c1daeSBarry Smith 4229371c9d4SSatish Balay static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name) { 423d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data; 424d1232d7fSToby Isaac 425d1232d7fSToby Isaac PetscFunctionBegin; 426d1232d7fSToby Isaac *name = vhdf5->filename; 427d1232d7fSToby Isaac PetscFunctionReturn(0); 428d1232d7fSToby Isaac } 429d1232d7fSToby Isaac 4309371c9d4SSatish Balay static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) { 4315dc64a97SVaclav Hapla /* 432b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4335dc64a97SVaclav Hapla */ 434b723ab35SVaclav Hapla 435b723ab35SVaclav Hapla PetscFunctionBegin; 436b723ab35SVaclav Hapla PetscFunctionReturn(0); 437b723ab35SVaclav Hapla } 438b723ab35SVaclav Hapla 4399371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) { 44019a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 44119a20e4cSMatthew G. Knepley 44219a20e4cSMatthew G. Knepley PetscFunctionBegin; 44319a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 44419a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 44519a20e4cSMatthew G. Knepley } 44619a20e4cSMatthew G. Knepley 4479371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) { 44819a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 44919a20e4cSMatthew G. Knepley 45019a20e4cSMatthew G. Knepley PetscFunctionBegin; 45119a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 45219a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 45319a20e4cSMatthew G. Knepley } 45419a20e4cSMatthew G. Knepley 45519a20e4cSMatthew G. Knepley /*@ 45619a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 45719a20e4cSMatthew G. Knepley 458811af0c4SBarry Smith Logically Collective on viewer 45919a20e4cSMatthew G. Knepley 46019a20e4cSMatthew G. Knepley Input Parameters: 461811af0c4SBarry Smith + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 462811af0c4SBarry Smith - flg - if `PETSC_TRUE` we will assume that timestepping is on 46319a20e4cSMatthew G. Knepley 464811af0c4SBarry Smith Options Database Key: 46519a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 46619a20e4cSMatthew G. Knepley 467811af0c4SBarry Smith Note: 46819a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 46919a20e4cSMatthew G. Knepley 47019a20e4cSMatthew G. Knepley Level: intermediate 47119a20e4cSMatthew G. Knepley 472811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 47319a20e4cSMatthew G. Knepley @*/ 4749371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) { 47519a20e4cSMatthew G. Knepley PetscFunctionBegin; 47619a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 477cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg)); 47819a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 47919a20e4cSMatthew G. Knepley } 48019a20e4cSMatthew G. Knepley 48119a20e4cSMatthew G. Knepley /*@ 48219a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 48319a20e4cSMatthew G. Knepley 48419a20e4cSMatthew G. Knepley Not collective 48519a20e4cSMatthew G. Knepley 48619a20e4cSMatthew G. Knepley Input Parameter: 487811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 48819a20e4cSMatthew G. Knepley 48919a20e4cSMatthew G. Knepley Output Parameter: 490811af0c4SBarry Smith . flg - if `PETSC_TRUE` we will assume that timestepping is on 49119a20e4cSMatthew G. Knepley 49219a20e4cSMatthew G. Knepley Level: intermediate 49319a20e4cSMatthew G. Knepley 494811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 49519a20e4cSMatthew G. Knepley @*/ 4969371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) { 49719a20e4cSMatthew G. Knepley PetscFunctionBegin; 49819a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 499cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg)); 50019a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 50119a20e4cSMatthew G. Knepley } 50219a20e4cSMatthew G. Knepley 5038556b5ebSBarry Smith /*MC 5048556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5058556b5ebSBarry Smith 506811af0c4SBarry Smith Level: beginner 507811af0c4SBarry Smith 508db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 509db781477SPatrick Sanan `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 510db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 511db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 5128556b5ebSBarry Smith M*/ 513d1232d7fSToby Isaac 5149371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) { 5155c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5165c6c1daeSBarry Smith 5175c6c1daeSBarry Smith PetscFunctionBegin; 51899335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 51999335c34SVaclav Hapla { 52099335c34SVaclav Hapla PetscMPIInt size; 5219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 5225f80ce2aSJacob 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)"); 52399335c34SVaclav Hapla } 52499335c34SVaclav Hapla #endif 52599335c34SVaclav Hapla 5269566063dSJacob Faibussowitsch PetscCall(PetscNewLog(v, &hdf5)); 5275c6c1daeSBarry Smith 5285c6c1daeSBarry Smith v->data = (void *)hdf5; 5295c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 53082ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 531b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5321b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5336226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5347e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 535908793a3SLisandro Dalcin hdf5->filename = NULL; 5365c6c1daeSBarry Smith hdf5->timestep = -1; 5370298fd71SBarry Smith hdf5->groups = NULL; 5385c6c1daeSBarry Smith 539792fecdfSBarry Smith PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER)); 5409c5072fbSVaclav Hapla 5419566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5)); 5429566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5)); 5439566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5)); 5449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5)); 5459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5)); 5469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5)); 5479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5)); 5489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5)); 5499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5)); 5509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5)); 5515c6c1daeSBarry Smith PetscFunctionReturn(0); 5525c6c1daeSBarry Smith } 5535c6c1daeSBarry Smith 5545c6c1daeSBarry Smith /*@C 555811af0c4SBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer` 5565c6c1daeSBarry Smith 557d083f849SBarry Smith Collective 5585c6c1daeSBarry Smith 5595c6c1daeSBarry Smith Input Parameters: 5605c6c1daeSBarry Smith + comm - MPI communicator 5615c6c1daeSBarry Smith . name - name of file 5625c6c1daeSBarry Smith - type - type of file 5635c6c1daeSBarry Smith 5645c6c1daeSBarry Smith Output Parameter: 565811af0c4SBarry Smith . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file 5665c6c1daeSBarry Smith 567811af0c4SBarry Smith Options Database Keys: 568a2b725a8SWilliam 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 569a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 57082ea9b62SBarry Smith 5715c6c1daeSBarry Smith Level: beginner 5725c6c1daeSBarry Smith 5737e4fd573SVaclav Hapla Notes: 5747e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 575811af0c4SBarry Smith .vb 576811af0c4SBarry Smith FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 577811af0c4SBarry Smith FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 578811af0c4SBarry 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] 579811af0c4SBarry Smith FILE_MODE_UPDATE - same as FILE_MODE_APPEND 580811af0c4SBarry Smith .ve 5817e4fd573SVaclav Hapla 582811af0c4SBarry Smith In case of `FILE_MODE_APPEND` / `FILE_MODE_UPDATE`, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified. 5837e4fd573SVaclav Hapla 5845c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5855c6c1daeSBarry Smith 586811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 587db781477SPatrick Sanan `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 588db781477SPatrick Sanan `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 5895c6c1daeSBarry Smith @*/ 5909371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) { 5915c6c1daeSBarry Smith PetscFunctionBegin; 5929566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, hdf5v)); 5939566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 5949566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 5959566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*hdf5v, name)); 5969566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*hdf5v)); 5975c6c1daeSBarry Smith PetscFunctionReturn(0); 5985c6c1daeSBarry Smith } 5995c6c1daeSBarry Smith 6005c6c1daeSBarry Smith /*@C 6015c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6025c6c1daeSBarry Smith 6035c6c1daeSBarry Smith Not collective 6045c6c1daeSBarry Smith 6055c6c1daeSBarry Smith Input Parameter: 606811af0c4SBarry Smith . viewer - the `PetscViewer` 6075c6c1daeSBarry Smith 6085c6c1daeSBarry Smith Output Parameter: 6095c6c1daeSBarry Smith . file_id - The file id 6105c6c1daeSBarry Smith 6115c6c1daeSBarry Smith Level: intermediate 6125c6c1daeSBarry Smith 613811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()` 6145c6c1daeSBarry Smith @*/ 6159371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) { 6165c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6175c6c1daeSBarry Smith 6185c6c1daeSBarry Smith PetscFunctionBegin; 6195c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 6205c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6215c6c1daeSBarry Smith PetscFunctionReturn(0); 6225c6c1daeSBarry Smith } 6235c6c1daeSBarry Smith 6245c6c1daeSBarry Smith /*@C 6255c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6265c6c1daeSBarry Smith 6275c6c1daeSBarry Smith Not collective 6285c6c1daeSBarry Smith 6295c6c1daeSBarry Smith Input Parameters: 630811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 6315c6c1daeSBarry Smith - name - The group name 6325c6c1daeSBarry Smith 6335c6c1daeSBarry Smith Level: intermediate 6345c6c1daeSBarry Smith 6352e361344SVaclav Hapla Notes: 6362e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 637811af0c4SBarry Smith .vb 638811af0c4SBarry 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. 639811af0c4SBarry Smith NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 640811af0c4SBarry Smith "." means the current group is pushed again. 641811af0c4SBarry Smith .ve 6422e361344SVaclav Hapla 6432e361344SVaclav Hapla Example: 6442e361344SVaclav Hapla Suppose the current group is "/a". 6452e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6462e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 6472e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 6482e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 6492e361344SVaclav Hapla 650811af0c4SBarry Smith Developer Note: 6512e361344SVaclav Hapla The root group "/" is internally stored as NULL. 652820fc2d1SVaclav Hapla 653*63cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 6545c6c1daeSBarry Smith @*/ 6559371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) { 6565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 6577d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6582e361344SVaclav Hapla size_t i, len; 6592e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6602e361344SVaclav Hapla const char *gname; 6615c6c1daeSBarry Smith 6625c6c1daeSBarry Smith PetscFunctionBegin; 6635c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 664820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 6659566063dSJacob Faibussowitsch PetscCall(PetscStrlen(name, &len)); 6662e361344SVaclav Hapla gname = NULL; 6672e361344SVaclav Hapla if (len) { 6682e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6692e361344SVaclav Hapla /* use current name */ 6702e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6712e361344SVaclav Hapla } else if (name[0] == '/') { 6722e361344SVaclav Hapla /* absolute */ 6732e361344SVaclav Hapla for (i = 1; i < len; i++) { 6742e361344SVaclav Hapla if (name[i] != '/') { 6752e361344SVaclav Hapla gname = name; 6762e361344SVaclav Hapla break; 6772e361344SVaclav Hapla } 6782e361344SVaclav Hapla } 6792e361344SVaclav Hapla } else { 6802e361344SVaclav Hapla /* relative */ 6812e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6829566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 6832e361344SVaclav Hapla gname = buf; 6842e361344SVaclav Hapla } 6852e361344SVaclav Hapla } 6869566063dSJacob Faibussowitsch PetscCall(PetscNew(&groupNode)); 6879566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name)); 6885c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6895c6c1daeSBarry Smith hdf5->groups = groupNode; 6905c6c1daeSBarry Smith PetscFunctionReturn(0); 6915c6c1daeSBarry Smith } 6925c6c1daeSBarry Smith 6933ef9c667SSatish Balay /*@ 6945c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6955c6c1daeSBarry Smith 6965c6c1daeSBarry Smith Not collective 6975c6c1daeSBarry Smith 6985c6c1daeSBarry Smith Input Parameter: 699811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7005c6c1daeSBarry Smith 7015c6c1daeSBarry Smith Level: intermediate 7025c6c1daeSBarry Smith 703*63cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 7045c6c1daeSBarry Smith @*/ 7059371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) { 7065c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7077d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7085c6c1daeSBarry Smith 7095c6c1daeSBarry Smith PetscFunctionBegin; 7105c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 7115f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7125c6c1daeSBarry Smith groupNode = hdf5->groups; 7135c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7149566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode->name)); 7159566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode)); 7165c6c1daeSBarry Smith PetscFunctionReturn(0); 7175c6c1daeSBarry Smith } 7185c6c1daeSBarry Smith 71977717648SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[]) { 7205c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 7215c6c1daeSBarry Smith 7225c6c1daeSBarry Smith PetscFunctionBegin; 7235c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 724c959eef4SJed Brown PetscValidPointer(name, 2); 725a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 7260298fd71SBarry Smith else *name = NULL; 7275c6c1daeSBarry Smith PetscFunctionReturn(0); 7285c6c1daeSBarry Smith } 7295c6c1daeSBarry Smith 7303014b61aSVaclav Hapla /*@C 731811af0c4SBarry Smith PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`, 732874270d9SVaclav Hapla and return this group's ID and file ID. 733811af0c4SBarry Smith If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID. 734874270d9SVaclav Hapla 735874270d9SVaclav Hapla Not collective 736874270d9SVaclav Hapla 7373014b61aSVaclav Hapla Input Parameters: 7383014b61aSVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 7393014b61aSVaclav Hapla - path - (Optional) The path relative to the pushed group 740874270d9SVaclav Hapla 741d8d19677SJose E. Roman Output Parameters: 742874270d9SVaclav Hapla + fileId - The HDF5 file ID 743874270d9SVaclav Hapla - groupId - The HDF5 group ID 744874270d9SVaclav Hapla 745811af0c4SBarry Smith Note: 7463014b61aSVaclav 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. 7473014b61aSVaclav Hapla So NULL or empty path means the current pushed group. 7483014b61aSVaclav Hapla 749e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 750e74616beSVaclav Hapla 751874270d9SVaclav Hapla Level: intermediate 752874270d9SVaclav Hapla 753*63cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()` 754874270d9SVaclav Hapla @*/ 7553014b61aSVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId) { 75690e03537SVaclav Hapla hid_t file_id; 75790e03537SVaclav Hapla H5O_type_t type; 7583014b61aSVaclav Hapla const char *fileName = NULL; 7593014b61aSVaclav Hapla char *groupName = NULL; 76076d59af2SVaclav Hapla PetscBool writable, has; 76154dbf706SVaclav Hapla 76254dbf706SVaclav Hapla PetscFunctionBegin; 7639566063dSJacob Faibussowitsch PetscCall(PetscViewerWritable(viewer, &writable)); 7649566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 7659566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetName(viewer, &fileName)); 7663014b61aSVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName)); 7679566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 76876d59af2SVaclav Hapla if (!has) { 7695f80ce2aSJacob 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); 770f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 77176d59af2SVaclav Hapla } 7725f80ce2aSJacob 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); 7733014b61aSVaclav Hapla PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT)); 7743014b61aSVaclav Hapla PetscCall(PetscFree(groupName)); 77554dbf706SVaclav Hapla *fileId = file_id; 77654dbf706SVaclav Hapla PetscFunctionReturn(0); 77754dbf706SVaclav Hapla } 77854dbf706SVaclav Hapla 779*63cb69f5SVaclav Hapla /*@C 780*63cb69f5SVaclav Hapla PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file 781*63cb69f5SVaclav Hapla 782*63cb69f5SVaclav Hapla Not collective 783*63cb69f5SVaclav Hapla 784*63cb69f5SVaclav Hapla Input Parameters: 785*63cb69f5SVaclav Hapla + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 786*63cb69f5SVaclav Hapla - path - (Optional) The path relative to the pushed group 787*63cb69f5SVaclav Hapla 788*63cb69f5SVaclav Hapla Note: 789*63cb69f5SVaclav 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. 790*63cb69f5SVaclav Hapla So NULL or empty path means the current pushed group. 791*63cb69f5SVaclav Hapla 792*63cb69f5SVaclav Hapla This will fail if the viewer is not writable. 793*63cb69f5SVaclav Hapla 794*63cb69f5SVaclav Hapla Level: intermediate 795*63cb69f5SVaclav Hapla 796*63cb69f5SVaclav Hapla .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 797*63cb69f5SVaclav Hapla @*/ 798*63cb69f5SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[]) { 799*63cb69f5SVaclav Hapla hid_t fileId, groupId; 800*63cb69f5SVaclav Hapla 801*63cb69f5SVaclav Hapla PetscFunctionBegin; 802*63cb69f5SVaclav Hapla PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created 803*63cb69f5SVaclav Hapla PetscCallHDF5(H5Gclose, (groupId)); 804*63cb69f5SVaclav Hapla PetscFunctionReturn(0); 805*63cb69f5SVaclav Hapla } 806*63cb69f5SVaclav Hapla 8073ef9c667SSatish Balay /*@ 808d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8095c6c1daeSBarry Smith 8105c6c1daeSBarry Smith Not collective 8115c6c1daeSBarry Smith 8125c6c1daeSBarry Smith Input Parameter: 813811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 8145c6c1daeSBarry Smith 8155c6c1daeSBarry Smith Level: intermediate 8165c6c1daeSBarry Smith 817d7dd068bSVaclav Hapla Notes: 818811af0c4SBarry Smith On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0. 819811af0c4SBarry Smith Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`. 820811af0c4SBarry 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()`]. 821811af0c4SBarry Smith Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 822811af0c4SBarry Smith Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`. 823d7dd068bSVaclav Hapla 824d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 825d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 826d7dd068bSVaclav Hapla 827811af0c4SBarry Smith Developer note: 828d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 829d7dd068bSVaclav Hapla 830811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 831d7dd068bSVaclav Hapla @*/ 8329371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) { 833d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 834d7dd068bSVaclav Hapla 835d7dd068bSVaclav Hapla PetscFunctionBegin; 836d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 8375f80ce2aSJacob Faibussowitsch PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 838d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 839d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 840d7dd068bSVaclav Hapla PetscFunctionReturn(0); 841d7dd068bSVaclav Hapla } 842d7dd068bSVaclav Hapla 843d7dd068bSVaclav Hapla /*@ 844d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 845d7dd068bSVaclav Hapla 846d7dd068bSVaclav Hapla Not collective 847d7dd068bSVaclav Hapla 848d7dd068bSVaclav Hapla Input Parameter: 849811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 850d7dd068bSVaclav Hapla 851d7dd068bSVaclav Hapla Level: intermediate 852d7dd068bSVaclav Hapla 853811af0c4SBarry Smith Note: 854811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 855d7dd068bSVaclav Hapla 856811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 857d7dd068bSVaclav Hapla @*/ 8589371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) { 859d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 860d7dd068bSVaclav Hapla 861d7dd068bSVaclav Hapla PetscFunctionBegin; 862d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 8635f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 864d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 865d7dd068bSVaclav Hapla PetscFunctionReturn(0); 866d7dd068bSVaclav Hapla } 867d7dd068bSVaclav Hapla 868d7dd068bSVaclav Hapla /*@ 869d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 870d7dd068bSVaclav Hapla 871d7dd068bSVaclav Hapla Not collective 872d7dd068bSVaclav Hapla 873d7dd068bSVaclav Hapla Input Parameter: 874811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 875d7dd068bSVaclav Hapla 876d7dd068bSVaclav Hapla Output Parameter: 877d7dd068bSVaclav Hapla . flg - is timestepping active? 878d7dd068bSVaclav Hapla 879d7dd068bSVaclav Hapla Level: intermediate 880d7dd068bSVaclav Hapla 881811af0c4SBarry Smith Note: 882811af0c4SBarry Smith See `PetscViewerHDF5PushTimestepping()` for details. 883d7dd068bSVaclav Hapla 884811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 885d7dd068bSVaclav Hapla @*/ 8869371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) { 887d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 888d7dd068bSVaclav Hapla 889d7dd068bSVaclav Hapla PetscFunctionBegin; 890d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 891d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 892d7dd068bSVaclav Hapla PetscFunctionReturn(0); 893d7dd068bSVaclav Hapla } 894d7dd068bSVaclav Hapla 895d7dd068bSVaclav Hapla /*@ 896d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 897d7dd068bSVaclav Hapla 898d7dd068bSVaclav Hapla Not collective 899d7dd068bSVaclav Hapla 900d7dd068bSVaclav Hapla Input Parameter: 901811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 902d7dd068bSVaclav Hapla 903d7dd068bSVaclav Hapla Level: intermediate 904d7dd068bSVaclav Hapla 905811af0c4SBarry Smith Note: 906811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 907d7dd068bSVaclav Hapla 908811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 9095c6c1daeSBarry Smith @*/ 9109371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) { 9115c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9125c6c1daeSBarry Smith 9135c6c1daeSBarry Smith PetscFunctionBegin; 9145c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 9155f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9165c6c1daeSBarry Smith ++hdf5->timestep; 9175c6c1daeSBarry Smith PetscFunctionReturn(0); 9185c6c1daeSBarry Smith } 9195c6c1daeSBarry Smith 9203ef9c667SSatish Balay /*@ 921d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9225c6c1daeSBarry Smith 923d7dd068bSVaclav Hapla Logically collective 9245c6c1daeSBarry Smith 9255c6c1daeSBarry Smith Input Parameters: 926811af0c4SBarry Smith + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 927d7dd068bSVaclav Hapla - timestep - The timestep 9285c6c1daeSBarry Smith 9295c6c1daeSBarry Smith Level: intermediate 9305c6c1daeSBarry Smith 931811af0c4SBarry Smith Note: 932811af0c4SBarry Smith This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 933d7dd068bSVaclav Hapla 934811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 9355c6c1daeSBarry Smith @*/ 9369371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) { 9375c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9385c6c1daeSBarry Smith 9395c6c1daeSBarry Smith PetscFunctionBegin; 9405c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 941d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 9425f80ce2aSJacob Faibussowitsch PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 9435f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9445c6c1daeSBarry Smith hdf5->timestep = timestep; 9455c6c1daeSBarry Smith PetscFunctionReturn(0); 9465c6c1daeSBarry Smith } 9475c6c1daeSBarry Smith 9483ef9c667SSatish Balay /*@ 9495c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 9505c6c1daeSBarry Smith 9515c6c1daeSBarry Smith Not collective 9525c6c1daeSBarry Smith 9535c6c1daeSBarry Smith Input Parameter: 954811af0c4SBarry Smith . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 9555c6c1daeSBarry Smith 9565c6c1daeSBarry Smith Output Parameter: 957d7dd068bSVaclav Hapla . timestep - The timestep 9585c6c1daeSBarry Smith 9595c6c1daeSBarry Smith Level: intermediate 9605c6c1daeSBarry Smith 961811af0c4SBarry Smith Note: 962811af0c4SBarry Smith This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 963d7dd068bSVaclav Hapla 964811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 9655c6c1daeSBarry Smith @*/ 9669371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) { 9675c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 9685c6c1daeSBarry Smith 9695c6c1daeSBarry Smith PetscFunctionBegin; 9705c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 971d7dd068bSVaclav Hapla PetscValidIntPointer(timestep, 2); 9725f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9735c6c1daeSBarry Smith *timestep = hdf5->timestep; 9745c6c1daeSBarry Smith PetscFunctionReturn(0); 9755c6c1daeSBarry Smith } 9765c6c1daeSBarry Smith 97736bce6e8SMatthew G. Knepley /*@C 97836bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 97936bce6e8SMatthew G. Knepley 98036bce6e8SMatthew G. Knepley Not collective 98136bce6e8SMatthew G. Knepley 98236bce6e8SMatthew G. Knepley Input Parameter: 983811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 98436bce6e8SMatthew G. Knepley 98536bce6e8SMatthew G. Knepley Output Parameter: 986811af0c4SBarry Smith . mtype - the HDF5 datatype 98736bce6e8SMatthew G. Knepley 98836bce6e8SMatthew G. Knepley Level: advanced 98936bce6e8SMatthew G. Knepley 990db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 99136bce6e8SMatthew G. Knepley @*/ 9929371c9d4SSatish Balay PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) { 99336bce6e8SMatthew G. Knepley PetscFunctionBegin; 99436bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 99536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 99636bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 99736bce6e8SMatthew G. Knepley #else 99836bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 99936bce6e8SMatthew G. Knepley #endif 100036bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 100136bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 100236bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1003c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1004de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 100536bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 100636bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 100736bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10087e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 100936bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 101036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 101136bce6e8SMatthew G. Knepley } 101236bce6e8SMatthew G. Knepley 101336bce6e8SMatthew G. Knepley /*@C 101436bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 101536bce6e8SMatthew G. Knepley 101636bce6e8SMatthew G. Knepley Not collective 101736bce6e8SMatthew G. Knepley 101836bce6e8SMatthew G. Knepley Input Parameter: 1019811af0c4SBarry Smith . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...) 102036bce6e8SMatthew G. Knepley 102136bce6e8SMatthew G. Knepley Output Parameter: 1022811af0c4SBarry Smith . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 102336bce6e8SMatthew G. Knepley 102436bce6e8SMatthew G. Knepley Level: advanced 102536bce6e8SMatthew G. Knepley 1026db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 102736bce6e8SMatthew G. Knepley @*/ 10289371c9d4SSatish Balay PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) { 102936bce6e8SMatthew G. Knepley PetscFunctionBegin; 103036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 103136bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 103236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 103336bce6e8SMatthew G. Knepley #else 103436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 103536bce6e8SMatthew G. Knepley #endif 103636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 103736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 103836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 103936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 104036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 104136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 10427e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 104336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 104436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 104536bce6e8SMatthew G. Knepley } 104636bce6e8SMatthew G. Knepley 1047df863907SAlex Fikl /*@C 1048b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 104936bce6e8SMatthew G. Knepley 1050343a20b2SVaclav Hapla Collective 1051343a20b2SVaclav Hapla 105236bce6e8SMatthew G. Knepley Input Parameters: 1053811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 10544302210dSVaclav Hapla . parent - The parent dataset/group name 105536bce6e8SMatthew G. Knepley . name - The attribute name 105636bce6e8SMatthew G. Knepley . datatype - The attribute type 105736bce6e8SMatthew G. Knepley - value - The attribute value 105836bce6e8SMatthew G. Knepley 105936bce6e8SMatthew G. Knepley Level: advanced 106036bce6e8SMatthew G. Knepley 1061811af0c4SBarry Smith Note: 1062343a20b2SVaclav 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. 10634302210dSVaclav Hapla 1064811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, 1065811af0c4SBarry Smith `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 106636bce6e8SMatthew G. Knepley @*/ 10679371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) { 10684302210dSVaclav Hapla char *parentAbsPath; 106960568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1070080f144cSVaclav Hapla PetscBool has; 107136bce6e8SMatthew G. Knepley 107236bce6e8SMatthew G. Knepley PetscFunctionBegin; 10735cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 10744302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1075c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 10764302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer, datatype, 4); 1077b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 107877717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 10799566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 10809566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 10819566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 10827e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 10837e97c476SMatthew G. Knepley size_t len; 10849566063dSJacob Faibussowitsch PetscCall(PetscStrlen((const char *)value, &len)); 1085792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 10867e97c476SMatthew G. Knepley } 10879566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1088792fecdfSBarry Smith PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR)); 1089792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1090080f144cSVaclav Hapla if (has) { 1091792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1092080f144cSVaclav Hapla } else { 1093792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1094080f144cSVaclav Hapla } 1095792fecdfSBarry Smith PetscCallHDF5(H5Awrite, (attribute, dtype, value)); 1096792fecdfSBarry Smith if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype)); 1097792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1098792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 1099792fecdfSBarry Smith PetscCallHDF5(H5Sclose, (dataspace)); 11009566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 110136bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 110236bce6e8SMatthew G. Knepley } 110336bce6e8SMatthew G. Knepley 1104df863907SAlex Fikl /*@C 1105811af0c4SBarry Smith PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name 1106577e0e04SVaclav Hapla 1107343a20b2SVaclav Hapla Collective 1108343a20b2SVaclav Hapla 1109577e0e04SVaclav Hapla Input Parameters: 1110811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1111577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1112577e0e04SVaclav Hapla . name - The attribute name 1113577e0e04SVaclav Hapla . datatype - The attribute type 1114577e0e04SVaclav Hapla - value - The attribute value 1115577e0e04SVaclav Hapla 1116811af0c4SBarry Smith Note: 1117343a20b2SVaclav 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). 1118811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1119577e0e04SVaclav Hapla 1120577e0e04SVaclav Hapla Level: advanced 1121577e0e04SVaclav Hapla 1122811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, 1123811af0c4SBarry Smith `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1124577e0e04SVaclav Hapla @*/ 11259371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) { 1126577e0e04SVaclav Hapla PetscFunctionBegin; 1127577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1128577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1129577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 1130b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 11319566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 11329566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 1133577e0e04SVaclav Hapla PetscFunctionReturn(0); 1134577e0e04SVaclav Hapla } 1135577e0e04SVaclav Hapla 1136577e0e04SVaclav Hapla /*@C 1137b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1138ad0c61feSMatthew G. Knepley 1139343a20b2SVaclav Hapla Collective 1140343a20b2SVaclav Hapla 1141ad0c61feSMatthew G. Knepley Input Parameters: 1142811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 11434302210dSVaclav Hapla . parent - The parent dataset/group name 1144ad0c61feSMatthew G. Knepley . name - The attribute name 1145a2d6be1bSVaclav Hapla . datatype - The attribute type 1146a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1147ad0c61feSMatthew G. Knepley 1148ad0c61feSMatthew G. Knepley Output Parameter: 1149a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1150ad0c61feSMatthew G. Knepley 1151a2d6be1bSVaclav Hapla Notes: 1152a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1153a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1154a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1155811af0c4SBarry Smith $ flg = `PETSC_FALSE`; 1156811af0c4SBarry Smith $ PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 1157a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1158a2d6be1bSVaclav Hapla 1159811af0c4SBarry Smith If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed. 1160708d4cb3SBarry Smith 1161343a20b2SVaclav 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. 1162ad0c61feSMatthew G. Knepley 1163343a20b2SVaclav Hapla Level: advanced 11644302210dSVaclav Hapla 1165811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1166ad0c61feSMatthew G. Knepley @*/ 11679371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) { 11684302210dSVaclav Hapla char *parentAbsPath; 1169a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1170969748fdSVaclav Hapla PetscBool has; 1171ad0c61feSMatthew G. Knepley 1172ad0c61feSMatthew G. Knepley PetscFunctionBegin; 11735cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 11744302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1175c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1176a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1177a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 11789566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 117977717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 11809566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 11819566063dSJacob Faibussowitsch if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1182a2d6be1bSVaclav Hapla if (!has) { 1183a2d6be1bSVaclav Hapla if (defaultValue) { 1184a2d6be1bSVaclav Hapla if (defaultValue != value) { 1185a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 11869566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value)); 1187a2d6be1bSVaclav Hapla } else { 1188a2d6be1bSVaclav Hapla size_t len; 1189792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype)); 11909566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(value, defaultValue, len)); 1191a2d6be1bSVaclav Hapla } 1192a2d6be1bSVaclav Hapla } 11939566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1194a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 119598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1196a2d6be1bSVaclav Hapla } 11979566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1198792fecdfSBarry Smith PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1199792fecdfSBarry Smith PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1200f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1201f0b58479SMatthew G. Knepley size_t len; 1202a2d6be1bSVaclav Hapla hid_t atype; 1203792fecdfSBarry Smith PetscCallHDF5Return(atype, H5Aget_type, (attribute)); 1204792fecdfSBarry Smith PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype)); 12059566063dSJacob Faibussowitsch PetscCall(PetscMalloc((len + 1) * sizeof(char), value)); 1206792fecdfSBarry Smith PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 1207792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value)); 1208708d4cb3SBarry Smith } else { 1209792fecdfSBarry Smith PetscCallHDF5(H5Aread, (attribute, dtype, value)); 1210708d4cb3SBarry Smith } 1211792fecdfSBarry Smith PetscCallHDF5(H5Aclose, (attribute)); 1212e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1213792fecdfSBarry Smith PetscCallHDF5(H5Oclose, (obj)); 12149566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1215ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1216ad0c61feSMatthew G. Knepley } 1217ad0c61feSMatthew G. Knepley 1218577e0e04SVaclav Hapla /*@C 1219811af0c4SBarry Smith PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name 1220577e0e04SVaclav Hapla 1221343a20b2SVaclav Hapla Collective 1222343a20b2SVaclav Hapla 1223577e0e04SVaclav Hapla Input Parameters: 1224811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1225577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1226577e0e04SVaclav Hapla . name - The attribute name 1227577e0e04SVaclav Hapla - datatype - The attribute type 1228577e0e04SVaclav Hapla 1229577e0e04SVaclav Hapla Output Parameter: 1230577e0e04SVaclav Hapla . value - The attribute value 1231577e0e04SVaclav Hapla 1232811af0c4SBarry Smith Note: 1233577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1234811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1235577e0e04SVaclav Hapla 1236577e0e04SVaclav Hapla Level: advanced 1237577e0e04SVaclav Hapla 1238811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1239577e0e04SVaclav Hapla @*/ 12409371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) { 1241577e0e04SVaclav Hapla PetscFunctionBegin; 1242577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1243577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1244577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 1245064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 12469566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 12479566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 1248577e0e04SVaclav Hapla PetscFunctionReturn(0); 1249577e0e04SVaclav Hapla } 1250577e0e04SVaclav Hapla 12519371c9d4SSatish Balay static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) { 1252a75c3fd4SVaclav Hapla htri_t exists; 1253a75c3fd4SVaclav Hapla hid_t group; 1254a75c3fd4SVaclav Hapla 1255a75c3fd4SVaclav Hapla PetscFunctionBegin; 1256792fecdfSBarry Smith PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT)); 1257792fecdfSBarry Smith if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT)); 1258a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1259792fecdfSBarry Smith PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1260792fecdfSBarry Smith PetscCallHDF5(H5Gclose, (group)); 1261a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1262a75c3fd4SVaclav Hapla } 1263a75c3fd4SVaclav Hapla *exists_ = (PetscBool)exists; 1264a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1265a75c3fd4SVaclav Hapla } 1266a75c3fd4SVaclav Hapla 12679371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) { 126890e03537SVaclav Hapla const char rootGroupName[] = "/"; 12695cdeb108SMatthew G. Knepley hid_t h5; 1270e5bf9ebcSVaclav Hapla PetscBool exists = PETSC_FALSE; 127115b861d2SVaclav Hapla PetscInt i; 127215b861d2SVaclav Hapla int n; 127385688ae2SVaclav Hapla char **hierarchy; 127485688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 12755cdeb108SMatthew G. Knepley 12765cdeb108SMatthew G. Knepley PetscFunctionBegin; 12775cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 127890e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 127990e03537SVaclav Hapla else name = rootGroupName; 1280ccd66a83SVaclav Hapla if (has) { 1281064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 12825cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1283ccd66a83SVaclav Hapla } 1284ccd66a83SVaclav Hapla if (otype) { 1285064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 128656cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1287ccd66a83SVaclav Hapla } 12889566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 128985688ae2SVaclav Hapla 129085688ae2SVaclav Hapla /* 129185688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 129285688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 129385688ae2SVaclav Hapla 1) whether it's a valid link 129485688ae2SVaclav Hapla 2) whether this link resolves to an object 129585688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 129685688ae2SVaclav Hapla */ 12979566063dSJacob Faibussowitsch PetscCall(PetscStrToArray(name, '/', &n, &hierarchy)); 129885688ae2SVaclav Hapla if (!n) { 129985688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1300ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1301ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 13029566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 130385688ae2SVaclav Hapla PetscFunctionReturn(0); 130485688ae2SVaclav Hapla } 130585688ae2SVaclav Hapla for (i = 0; i < n; i++) { 13069566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf, "/")); 13079566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf, hierarchy[i])); 13089566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1309a75c3fd4SVaclav Hapla if (!exists) break; 131085688ae2SVaclav Hapla } 13119566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 131285688ae2SVaclav Hapla 131385688ae2SVaclav Hapla /* If the object exists, get its type */ 1314ccd66a83SVaclav Hapla if (exists && otype) { 13155cdeb108SMatthew G. Knepley H5O_info_t info; 13165cdeb108SMatthew G. Knepley 131776276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1318792fecdfSBarry Smith PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT)); 131956cc0592SVaclav Hapla *otype = info.type; 13205cdeb108SMatthew G. Knepley } 1321ccd66a83SVaclav Hapla if (has) *has = exists; 13225cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 13235cdeb108SMatthew G. Knepley } 13245cdeb108SMatthew G. Knepley 13254302210dSVaclav Hapla /*@C 13260a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13270a9f38efSVaclav Hapla 1328343a20b2SVaclav Hapla Collective 1329343a20b2SVaclav Hapla 13300a9f38efSVaclav Hapla Input Parameters: 1331811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 133277717648SVaclav Hapla - path - The path relative to the pushed group, or NULL 13330a9f38efSVaclav Hapla 13340a9f38efSVaclav Hapla Output Parameter: 13350a9f38efSVaclav Hapla . has - Flag for group existence 13360a9f38efSVaclav Hapla 13370a9f38efSVaclav Hapla Level: advanced 13380a9f38efSVaclav Hapla 13394302210dSVaclav Hapla Notes: 1340785c443eSVaclav 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. 1341785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 13424302210dSVaclav Hapla 1343811af0c4SBarry Smith If path exists but is not a group, `PETSC_FALSE` is returned. 1344811af0c4SBarry Smith 1345811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 13460a9f38efSVaclav Hapla @*/ 13479371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) { 13480a9f38efSVaclav Hapla H5O_type_t type; 13494302210dSVaclav Hapla char *abspath; 13500a9f38efSVaclav Hapla 13510a9f38efSVaclav Hapla PetscFunctionBegin; 13520a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 13534302210dSVaclav Hapla if (path) PetscValidCharPointer(path, 2); 13544302210dSVaclav Hapla PetscValidBoolPointer(has, 3); 135577717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 13569566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 13574302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 13589566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 13590a9f38efSVaclav Hapla PetscFunctionReturn(0); 13600a9f38efSVaclav Hapla } 13610a9f38efSVaclav Hapla 136289e0ef10SVaclav Hapla /*@C 136389e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 136489e0ef10SVaclav Hapla 1365343a20b2SVaclav Hapla Collective 1366343a20b2SVaclav Hapla 136789e0ef10SVaclav Hapla Input Parameters: 1368811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 136989e0ef10SVaclav Hapla - path - The dataset path 137089e0ef10SVaclav Hapla 137189e0ef10SVaclav Hapla Output Parameter: 137289e0ef10SVaclav Hapla . has - Flag whether dataset exists 137389e0ef10SVaclav Hapla 137489e0ef10SVaclav Hapla Level: advanced 137589e0ef10SVaclav Hapla 137689e0ef10SVaclav Hapla Notes: 1377343a20b2SVaclav 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. 137889e0ef10SVaclav Hapla 1379811af0c4SBarry Smith If path is NULL or empty, has is set to `PETSC_FALSE`. 1380811af0c4SBarry Smith 1381811af0c4SBarry Smith If path exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1382811af0c4SBarry Smith 1383811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 138489e0ef10SVaclav Hapla @*/ 13859371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) { 138689e0ef10SVaclav Hapla H5O_type_t type; 138789e0ef10SVaclav Hapla char *abspath; 138889e0ef10SVaclav Hapla 138989e0ef10SVaclav Hapla PetscFunctionBegin; 139089e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 139189e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path, 2); 139289e0ef10SVaclav Hapla PetscValidBoolPointer(has, 3); 139377717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 13949566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 139589e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 13969566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 139789e0ef10SVaclav Hapla PetscFunctionReturn(0); 139889e0ef10SVaclav Hapla } 139989e0ef10SVaclav Hapla 14000a9f38efSVaclav Hapla /*@ 1401e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1402ecce1506SVaclav Hapla 1403343a20b2SVaclav Hapla Collective 1404343a20b2SVaclav Hapla 1405ecce1506SVaclav Hapla Input Parameters: 1406811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1407ecce1506SVaclav Hapla - obj - The named object 1408ecce1506SVaclav Hapla 1409ecce1506SVaclav Hapla Output Parameter: 141089e0ef10SVaclav Hapla . has - Flag for dataset existence 1411ecce1506SVaclav Hapla 1412e3f143f7SVaclav Hapla Notes: 141389e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1414811af0c4SBarry Smith 1415811af0c4SBarry Smith If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1416e3f143f7SVaclav Hapla 1417ecce1506SVaclav Hapla Level: advanced 1418ecce1506SVaclav Hapla 1419811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1420ecce1506SVaclav Hapla @*/ 14219371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) { 142289e0ef10SVaclav Hapla size_t len; 1423ecce1506SVaclav Hapla 1424ecce1506SVaclav Hapla PetscFunctionBegin; 1425c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1426c57a1a9aSVaclav Hapla PetscValidHeader(obj, 2); 14274302210dSVaclav Hapla PetscValidBoolPointer(has, 3); 14289566063dSJacob Faibussowitsch PetscCall(PetscStrlen(obj->name, &len)); 14295f80ce2aSJacob Faibussowitsch PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 14309566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 1431ecce1506SVaclav Hapla PetscFunctionReturn(0); 1432ecce1506SVaclav Hapla } 1433ecce1506SVaclav Hapla 1434df863907SAlex Fikl /*@C 1435ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1436e7bdbf8eSMatthew G. Knepley 1437343a20b2SVaclav Hapla Collective 1438343a20b2SVaclav Hapla 1439e7bdbf8eSMatthew G. Knepley Input Parameters: 1440811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 14414302210dSVaclav Hapla . parent - The parent dataset/group name 1442e7bdbf8eSMatthew G. Knepley - name - The attribute name 1443e7bdbf8eSMatthew G. Knepley 1444e7bdbf8eSMatthew G. Knepley Output Parameter: 1445e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1446e7bdbf8eSMatthew G. Knepley 1447e7bdbf8eSMatthew G. Knepley Level: advanced 1448e7bdbf8eSMatthew G. Knepley 1449811af0c4SBarry Smith Note: 1450343a20b2SVaclav 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. 14514302210dSVaclav Hapla 1452811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1453e7bdbf8eSMatthew G. Knepley @*/ 14549371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) { 14554302210dSVaclav Hapla char *parentAbsPath; 1456e7bdbf8eSMatthew G. Knepley 1457e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 14585cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 14594302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1460c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 14614302210dSVaclav Hapla PetscValidBoolPointer(has, 4); 146277717648SVaclav Hapla PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 14639566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 14649566063dSJacob Faibussowitsch if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 14659566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 146606db490cSVaclav Hapla PetscFunctionReturn(0); 146706db490cSVaclav Hapla } 146806db490cSVaclav Hapla 1469577e0e04SVaclav Hapla /*@C 1470811af0c4SBarry Smith PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name 1471577e0e04SVaclav Hapla 1472343a20b2SVaclav Hapla Collective 1473343a20b2SVaclav Hapla 1474577e0e04SVaclav Hapla Input Parameters: 1475811af0c4SBarry Smith + viewer - The `PETSCVIEWERHDF5` viewer 1476577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1477577e0e04SVaclav Hapla - name - The attribute name 1478577e0e04SVaclav Hapla 1479577e0e04SVaclav Hapla Output Parameter: 1480577e0e04SVaclav Hapla . has - Flag for attribute existence 1481577e0e04SVaclav Hapla 1482811af0c4SBarry Smith Note: 1483577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1484811af0c4SBarry Smith You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1485577e0e04SVaclav Hapla 1486577e0e04SVaclav Hapla Level: advanced 1487577e0e04SVaclav Hapla 1488811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1489577e0e04SVaclav Hapla @*/ 14909371c9d4SSatish Balay PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) { 1491577e0e04SVaclav Hapla PetscFunctionBegin; 1492577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1493577e0e04SVaclav Hapla PetscValidHeader(obj, 2); 1494577e0e04SVaclav Hapla PetscValidCharPointer(name, 3); 14954302210dSVaclav Hapla PetscValidBoolPointer(has, 4); 14969566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 14979566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 1498577e0e04SVaclav Hapla PetscFunctionReturn(0); 1499577e0e04SVaclav Hapla } 1500577e0e04SVaclav Hapla 15019371c9d4SSatish Balay static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) { 150206db490cSVaclav Hapla hid_t h5; 150306db490cSVaclav Hapla htri_t hhas; 150406db490cSVaclav Hapla 150506db490cSVaclav Hapla PetscFunctionBegin; 15069566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1507792fecdfSBarry Smith PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT)); 15085cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1509e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1510e7bdbf8eSMatthew G. Knepley } 1511e7bdbf8eSMatthew G. Knepley 1512a75e6a4aSMatthew G. Knepley /* 1513a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1514a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1515a75e6a4aSMatthew G. Knepley */ 1516d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1517a75e6a4aSMatthew G. Knepley 1518a75e6a4aSMatthew G. Knepley /*@C 1519811af0c4SBarry Smith PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator. 1520a75e6a4aSMatthew G. Knepley 1521d083f849SBarry Smith Collective 1522a75e6a4aSMatthew G. Knepley 1523a75e6a4aSMatthew G. Knepley Input Parameter: 1524811af0c4SBarry Smith . comm - the MPI communicator to share the HDF5 `PetscViewer` 1525a75e6a4aSMatthew G. Knepley 1526a75e6a4aSMatthew G. Knepley Level: intermediate 1527a75e6a4aSMatthew G. Knepley 1528811af0c4SBarry Smith Options Database Key: 152910699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1530a75e6a4aSMatthew G. Knepley 1531811af0c4SBarry Smith Environmental variable: 1532811af0c4SBarry Smith . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file 1533a75e6a4aSMatthew G. Knepley 1534811af0c4SBarry Smith Note: 1535811af0c4SBarry Smith Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return 1536811af0c4SBarry Smith an error code. The HDF5 `PetscViewer` is usually used in the form 1537a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1538a75e6a4aSMatthew G. Knepley 1539811af0c4SBarry Smith .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1540a75e6a4aSMatthew G. Knepley @*/ 15419371c9d4SSatish Balay PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) { 1542a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1543a75e6a4aSMatthew G. Knepley PetscBool flg; 1544a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1545a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1546a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1547a75e6a4aSMatthew G. Knepley 1548a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 15499371c9d4SSatish Balay ierr = PetscCommDuplicate(comm, &ncomm, NULL); 15509371c9d4SSatish Balay if (ierr) { 15519371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 15529371c9d4SSatish Balay PetscFunctionReturn(NULL); 15539371c9d4SSatish Balay } 1554a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1555908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL); 15569371c9d4SSatish Balay if (ierr) { 15579371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 15589371c9d4SSatish Balay PetscFunctionReturn(NULL); 15599371c9d4SSatish Balay } 1560a75e6a4aSMatthew G. Knepley } 156147435625SJed Brown ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg); 15629371c9d4SSatish Balay if (ierr) { 15639371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 15649371c9d4SSatish Balay PetscFunctionReturn(NULL); 15659371c9d4SSatish Balay } 1566a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1567a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg); 15689371c9d4SSatish Balay if (ierr) { 15699371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 15709371c9d4SSatish Balay PetscFunctionReturn(NULL); 15719371c9d4SSatish Balay } 1572a75e6a4aSMatthew G. Knepley if (!flg) { 1573a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname, "output.h5"); 15749371c9d4SSatish Balay if (ierr) { 15759371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 15769371c9d4SSatish Balay PetscFunctionReturn(NULL); 15779371c9d4SSatish Balay } 1578a75e6a4aSMatthew G. Knepley } 1579a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer); 15809371c9d4SSatish Balay if (ierr) { 15819371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 15829371c9d4SSatish Balay PetscFunctionReturn(NULL); 15839371c9d4SSatish Balay } 1584a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 15859371c9d4SSatish Balay if (ierr) { 15869371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 15879371c9d4SSatish Balay PetscFunctionReturn(NULL); 15889371c9d4SSatish Balay } 158947435625SJed Brown ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer); 15909371c9d4SSatish Balay if (ierr) { 15919371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 15929371c9d4SSatish Balay PetscFunctionReturn(NULL); 15939371c9d4SSatish Balay } 1594a75e6a4aSMatthew G. Knepley } 1595a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 15969371c9d4SSatish Balay if (ierr) { 15979371c9d4SSatish Balay PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 15989371c9d4SSatish Balay PetscFunctionReturn(NULL); 15999371c9d4SSatish Balay } 1600a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1601a75e6a4aSMatthew G. Knepley } 1602