1 #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 2 3 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool *, H5O_type_t *); 4 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool *); 5 6 /*@C 7 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with `PetscViewerHDF5PushGroup()`/`PetscViewerHDF5PopGroup()`. 8 9 Not collective 10 11 Input Parameters: 12 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 13 - path - (Optional) The path relative to the pushed group 14 15 Output Parameter: 16 . abspath - The absolute HDF5 path (group) 17 18 Level: intermediate 19 20 Notes: 21 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 22 So NULL or empty path means the current pushed group. 23 24 The output abspath is newly allocated so needs to be freed. 25 26 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 27 @*/ 28 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char path[], char *abspath[]) { 29 size_t len; 30 PetscBool relative = PETSC_FALSE; 31 const char *group; 32 char buf[PETSC_MAX_PATH_LEN] = ""; 33 34 PetscFunctionBegin; 35 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 36 if (path) PetscValidCharPointer(path, 2); 37 PetscValidPointer(abspath, 3); 38 PetscCall(PetscViewerHDF5GetGroup_Internal(viewer, &group)); 39 PetscCall(PetscStrlen(path, &len)); 40 relative = (PetscBool)(!len || path[0] != '/'); 41 if (relative) { 42 PetscCall(PetscStrcpy(buf, group)); 43 if (!group || len) PetscCall(PetscStrcat(buf, "/")); 44 PetscCall(PetscStrcat(buf, path)); 45 PetscCall(PetscStrallocpy(buf, abspath)); 46 } else { 47 PetscCall(PetscStrallocpy(path, abspath)); 48 } 49 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) { 53 PetscBool has; 54 55 PetscFunctionBegin; 56 PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has)); 57 if (!has) { 58 char *group; 59 PetscCall(PetscViewerHDF5GetGroup(viewer, NULL, &group)); 60 SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group); 61 } 62 PetscFunctionReturn(0); 63 } 64 65 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscViewer v, PetscOptionItems *PetscOptionsObject) { 66 PetscBool flg = PETSC_FALSE, set; 67 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 68 69 PetscFunctionBegin; 70 PetscOptionsHeadBegin(PetscOptionsObject, "HDF5 PetscViewer Options"); 71 PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2", "1d Vectors get 2 dimensions in HDF5", "PetscViewerHDF5SetBaseDimension2", hdf5->basedimension2, &hdf5->basedimension2, NULL)); 72 PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output", "Force data to be written in single precision", "PetscViewerHDF5SetSPOutput", hdf5->spoutput, &hdf5->spoutput, NULL)); 73 PetscCall(PetscOptionsBool("-viewer_hdf5_collective", "Enable collective transfer mode", "PetscViewerHDF5SetCollective", flg, &flg, &set)); 74 if (set) PetscCall(PetscViewerHDF5SetCollective(v, flg)); 75 flg = PETSC_FALSE; 76 PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping", "Set default timestepping state", "PetscViewerHDF5SetDefaultTimestepping", flg, &flg, &set)); 77 if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v, flg)); 78 PetscOptionsHeadEnd(); 79 PetscFunctionReturn(0); 80 } 81 82 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v, PetscViewer viewer) { 83 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)v->data; 84 PetscBool flg; 85 86 PetscFunctionBegin; 87 if (hdf5->filename) PetscCall(PetscViewerASCIIPrintf(viewer, "Filename: %s\n", hdf5->filename)); 88 PetscCall(PetscViewerASCIIPrintf(viewer, "Vectors with blocksize 1 saved as 2D datasets: %s\n", PetscBools[hdf5->basedimension2])); 89 PetscCall(PetscViewerASCIIPrintf(viewer, "Enforce single precision storage: %s\n", PetscBools[hdf5->spoutput])); 90 PetscCall(PetscViewerHDF5GetCollective(v, &flg)); 91 PetscCall(PetscViewerASCIIPrintf(viewer, "MPI-IO transfer mode: %s\n", flg ? "collective" : "independent")); 92 PetscCall(PetscViewerASCIIPrintf(viewer, "Default timestepping: %s\n", PetscBools[hdf5->defTimestepping])); 93 PetscFunctionReturn(0); 94 } 95 96 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) { 97 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 98 99 PetscFunctionBegin; 100 PetscCall(PetscFree(hdf5->filename)); 101 if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 102 PetscFunctionReturn(0); 103 } 104 105 static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) { 106 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 107 108 PetscFunctionBegin; 109 if (hdf5->file_id) PetscCallHDF5(H5Fflush, (hdf5->file_id, H5F_SCOPE_LOCAL)); 110 PetscFunctionReturn(0); 111 } 112 113 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) { 114 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 115 116 PetscFunctionBegin; 117 PetscCallHDF5(H5Pclose, (hdf5->dxpl_id)); 118 PetscCall(PetscViewerFileClose_HDF5(viewer)); 119 while (hdf5->groups) { 120 PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 121 122 PetscCall(PetscFree(hdf5->groups->name)); 123 PetscCall(PetscFree(hdf5->groups)); 124 hdf5->groups = tmp; 125 } 126 PetscCall(PetscFree(hdf5)); 127 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL)); 128 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL)); 129 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL)); 130 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL)); 131 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetBaseDimension2_C", NULL)); 132 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetSPOutput_C", NULL)); 133 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetCollective_C", NULL)); 134 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetCollective_C", NULL)); 135 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5GetDefaultTimestepping_C", NULL)); 136 PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerHDF5SetDefaultTimestepping_C", NULL)); 137 PetscFunctionReturn(0); 138 } 139 140 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) { 141 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 142 143 PetscFunctionBegin; 144 hdf5->btype = type; 145 PetscFunctionReturn(0); 146 } 147 148 static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) { 149 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 150 151 PetscFunctionBegin; 152 *type = hdf5->btype; 153 PetscFunctionReturn(0); 154 } 155 156 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) { 157 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 158 159 PetscFunctionBegin; 160 hdf5->basedimension2 = flg; 161 PetscFunctionReturn(0); 162 } 163 164 /*@ 165 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 166 dimension of 2. 167 168 Logically Collective on viewer 169 170 Input Parameters: 171 + viewer - the `PetscViewer`; if it is a `PETSCVIEWERHDF5` then this command is ignored 172 - flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1 173 174 Options Database Key: 175 . -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 176 177 Note: 178 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 179 of one when the dimension is lower. Others think the option is crazy. 180 181 Level: intermediate 182 183 .seealso: `PETSCVIEWERHDF5`, PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 184 @*/ 185 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer, PetscBool flg) { 186 PetscFunctionBegin; 187 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 188 PetscTryMethod(viewer, "PetscViewerHDF5SetBaseDimension2_C", (PetscViewer, PetscBool), (viewer, flg)); 189 PetscFunctionReturn(0); 190 } 191 192 /*@ 193 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 194 dimension of 2. 195 196 Logically Collective on viewer 197 198 Input Parameter: 199 . viewer - the `PetscViewer`, must be `PETSCVIEWERHDF5` 200 201 Output Parameter: 202 . flg - if `PETSC_TRUE` the vector will always have at least a dimension of 2 even if that first dimension is of size 1 203 204 Note: 205 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 206 of one when the dimension is lower. Others think the option is crazy. 207 208 Level: intermediate 209 210 .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 211 @*/ 212 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer, PetscBool *flg) { 213 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 214 215 PetscFunctionBegin; 216 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 217 *flg = hdf5->basedimension2; 218 PetscFunctionReturn(0); 219 } 220 221 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) { 222 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 223 224 PetscFunctionBegin; 225 hdf5->spoutput = flg; 226 PetscFunctionReturn(0); 227 } 228 229 /*@ 230 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 231 compiled with double precision `PetscReal`. 232 233 Logically Collective on viewer 234 235 Input Parameters: 236 + viewer - the PetscViewer; if it is a `PETSCVIEWERHDF5` then this command is ignored 237 - flg - if `PETSC_TRUE` the data will be written to disk with single precision 238 239 Options Database Key: 240 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 241 242 Note: 243 Setting this option does not make any difference if PETSc is compiled with single precision 244 in the first place. It does not affect reading datasets (HDF5 handle this internally). 245 246 Level: intermediate 247 248 .seealso: `PETSCVIEWERHDF5`, `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 249 `PetscReal`, `PetscViewerHDF5GetSPOutput()` 250 @*/ 251 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer, PetscBool flg) { 252 PetscFunctionBegin; 253 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 254 PetscTryMethod(viewer, "PetscViewerHDF5SetSPOutput_C", (PetscViewer, PetscBool), (viewer, flg)); 255 PetscFunctionReturn(0); 256 } 257 258 /*@ 259 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 260 compiled with double precision `PetscReal`. 261 262 Logically Collective on viewer 263 264 Input Parameter: 265 . viewer - the PetscViewer, must be of type `PETSCVIEWERHDF5` 266 267 Output Parameter: 268 . flg - if `PETSC_TRUE` the data will be written to disk with single precision 269 270 Level: intermediate 271 272 .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 273 `PetscReal`, `PetscViewerHDF5SetSPOutput()` 274 @*/ 275 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer, PetscBool *flg) { 276 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 277 278 PetscFunctionBegin; 279 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 280 *flg = hdf5->spoutput; 281 PetscFunctionReturn(0); 282 } 283 284 static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) { 285 PetscFunctionBegin; 286 /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 287 - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 288 #if H5_VERSION_GE(1, 10, 3) && defined(H5_HAVE_PARALLEL) 289 { 290 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 291 PetscCallHDF5(H5Pset_dxpl_mpio, (hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 292 } 293 #else 294 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")); 295 #endif 296 PetscFunctionReturn(0); 297 } 298 299 /*@ 300 PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 301 302 Logically Collective; flg must contain common value 303 304 Input Parameters: 305 + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 306 - flg - `PETSC_TRUE` for collective mode; `PETSC_FALSE` for independent mode (default) 307 308 Options Database Key: 309 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 310 311 Note: 312 Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 313 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. 314 315 Developer Note: 316 In the HDF5 layer, `PETSC_TRUE` / `PETSC_FALSE` means `H5Pset_dxpl_mpio()` is called with `H5FD_MPIO_COLLECTIVE` / `H5FD_MPIO_INDEPENDENT`, respectively. 317 This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 318 See HDF5 documentation and MPI-IO documentation for details. 319 320 Level: intermediate 321 322 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 323 @*/ 324 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer, PetscBool flg) { 325 PetscFunctionBegin; 326 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 327 PetscValidLogicalCollectiveBool(viewer, flg, 2); 328 PetscTryMethod(viewer, "PetscViewerHDF5SetCollective_C", (PetscViewer, PetscBool), (viewer, flg)); 329 PetscFunctionReturn(0); 330 } 331 332 static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) { 333 #if defined(H5_HAVE_PARALLEL) 334 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 335 H5FD_mpio_xfer_t mode; 336 #endif 337 338 PetscFunctionBegin; 339 #if !defined(H5_HAVE_PARALLEL) 340 *flg = PETSC_FALSE; 341 #else 342 PetscCallHDF5(H5Pget_dxpl_mpio, (hdf5->dxpl_id, &mode)); 343 *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 344 #endif 345 PetscFunctionReturn(0); 346 } 347 348 /*@ 349 PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 350 351 Not Collective 352 353 Input Parameters: 354 . viewer - the `PETSCVIEWERHDF5` `PetscViewer` 355 356 Output Parameters: 357 . flg - the flag 358 359 Level: intermediate 360 361 Note: 362 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. 363 For more details, see `PetscViewerHDF5SetCollective()`. 364 365 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 366 @*/ 367 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer, PetscBool *flg) { 368 PetscFunctionBegin; 369 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 370 PetscValidBoolPointer(flg, 2); 371 372 PetscUseMethod(viewer, "PetscViewerHDF5GetCollective_C", (PetscViewer, PetscBool *), (viewer, flg)); 373 PetscFunctionReturn(0); 374 } 375 376 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) { 377 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 378 hid_t plist_id; 379 380 PetscFunctionBegin; 381 if (hdf5->file_id) PetscCallHDF5(H5Fclose, (hdf5->file_id)); 382 if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 383 PetscCall(PetscStrallocpy(name, &hdf5->filename)); 384 /* Set up file access property list with parallel I/O access */ 385 PetscCallHDF5Return(plist_id, H5Pcreate, (H5P_FILE_ACCESS)); 386 #if defined(H5_HAVE_PARALLEL) 387 PetscCallHDF5(H5Pset_fapl_mpio, (plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 388 #endif 389 /* Create or open the file collectively */ 390 switch (hdf5->btype) { 391 case FILE_MODE_READ: 392 if (PetscDefined(USE_DEBUG)) { 393 PetscMPIInt rank; 394 PetscBool flg; 395 396 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 397 if (rank == 0) { 398 PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 399 PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_FILE_OPEN, "File %s requested for reading does not exist", hdf5->filename); 400 } 401 PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 402 } 403 PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDONLY, plist_id)); 404 break; 405 case FILE_MODE_APPEND: 406 case FILE_MODE_UPDATE: { 407 PetscBool flg; 408 PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 409 if (flg) PetscCallHDF5Return(hdf5->file_id, H5Fopen, (name, H5F_ACC_RDWR, plist_id)); 410 else PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 411 break; 412 } 413 case FILE_MODE_WRITE: PetscCallHDF5Return(hdf5->file_id, H5Fcreate, (name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); break; 414 case FILE_MODE_UNDEFINED: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 415 default: SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unsupported file mode %s", PetscFileModes[hdf5->btype]); 416 } 417 PetscCheck(hdf5->file_id >= 0, PETSC_COMM_SELF, PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 418 PetscCallHDF5(H5Pclose, (plist_id)); 419 PetscFunctionReturn(0); 420 } 421 422 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer, const char **name) { 423 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5 *)viewer->data; 424 425 PetscFunctionBegin; 426 *name = vhdf5->filename; 427 PetscFunctionReturn(0); 428 } 429 430 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) { 431 /* 432 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 433 */ 434 435 PetscFunctionBegin; 436 PetscFunctionReturn(0); 437 } 438 439 static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) { 440 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 441 442 PetscFunctionBegin; 443 hdf5->defTimestepping = flg; 444 PetscFunctionReturn(0); 445 } 446 447 static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) { 448 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 449 450 PetscFunctionBegin; 451 *flg = hdf5->defTimestepping; 452 PetscFunctionReturn(0); 453 } 454 455 /*@ 456 PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 457 458 Logically Collective on viewer 459 460 Input Parameters: 461 + viewer - the `PetscViewer`; if it is not `PETSCVIEWERHDF5` then this command is ignored 462 - flg - if `PETSC_TRUE` we will assume that timestepping is on 463 464 Options Database Key: 465 . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 466 467 Note: 468 If the timestepping attribute is not found for an object, then the default timestepping is used 469 470 Level: intermediate 471 472 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 473 @*/ 474 PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) { 475 PetscFunctionBegin; 476 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 477 PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer, PetscBool), (viewer, flg)); 478 PetscFunctionReturn(0); 479 } 480 481 /*@ 482 PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 483 484 Not collective 485 486 Input Parameter: 487 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 488 489 Output Parameter: 490 . flg - if `PETSC_TRUE` we will assume that timestepping is on 491 492 Level: intermediate 493 494 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 495 @*/ 496 PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) { 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 499 PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer, PetscBool *), (viewer, flg)); 500 PetscFunctionReturn(0); 501 } 502 503 /*MC 504 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 505 506 Level: beginner 507 508 .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 509 `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 510 `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 511 `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 512 M*/ 513 514 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) { 515 PetscViewer_HDF5 *hdf5; 516 517 PetscFunctionBegin; 518 #if !defined(H5_HAVE_PARALLEL) 519 { 520 PetscMPIInt size; 521 PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 522 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)"); 523 } 524 #endif 525 526 PetscCall(PetscNew(&hdf5)); 527 528 v->data = (void *)hdf5; 529 v->ops->destroy = PetscViewerDestroy_HDF5; 530 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 531 v->ops->setup = PetscViewerSetUp_HDF5; 532 v->ops->view = PetscViewerView_HDF5; 533 v->ops->flush = PetscViewerFlush_HDF5; 534 hdf5->btype = FILE_MODE_UNDEFINED; 535 hdf5->filename = NULL; 536 hdf5->timestep = -1; 537 hdf5->groups = NULL; 538 539 PetscCallHDF5Return(hdf5->dxpl_id, H5Pcreate, (H5P_DATASET_XFER)); 540 541 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_HDF5)); 542 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_HDF5)); 543 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_HDF5)); 544 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_HDF5)); 545 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetBaseDimension2_C", PetscViewerHDF5SetBaseDimension2_HDF5)); 546 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetSPOutput_C", PetscViewerHDF5SetSPOutput_HDF5)); 547 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetCollective_C", PetscViewerHDF5SetCollective_HDF5)); 548 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetCollective_C", PetscViewerHDF5GetCollective_HDF5)); 549 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5GetDefaultTimestepping_C", PetscViewerHDF5GetDefaultTimestepping_HDF5)); 550 PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerHDF5SetDefaultTimestepping_C", PetscViewerHDF5SetDefaultTimestepping_HDF5)); 551 PetscFunctionReturn(0); 552 } 553 554 /*@C 555 PetscViewerHDF5Open - Opens a file for HDF5 input/output as a `PETSCVIEWERHDF5` `PetscViewer` 556 557 Collective 558 559 Input Parameters: 560 + comm - MPI communicator 561 . name - name of file 562 - type - type of file 563 564 Output Parameter: 565 . hdf5v - `PetscViewer` for HDF5 input/output to use with the specified file 566 567 Options Database Keys: 568 + -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 569 - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 570 571 Level: beginner 572 573 Notes: 574 Reading is always available, regardless of the mode. Available modes are 575 .vb 576 FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 577 FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 578 FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL] 579 FILE_MODE_UPDATE - same as FILE_MODE_APPEND 580 .ve 581 582 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. 583 584 This PetscViewer should be destroyed with PetscViewerDestroy(). 585 586 .seealso: `PETSCVIEWERHDF5`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 587 `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 588 `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 589 @*/ 590 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) { 591 PetscFunctionBegin; 592 PetscCall(PetscViewerCreate(comm, hdf5v)); 593 PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 594 PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 595 PetscCall(PetscViewerFileSetName(*hdf5v, name)); 596 PetscCall(PetscViewerSetFromOptions(*hdf5v)); 597 PetscFunctionReturn(0); 598 } 599 600 /*@C 601 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 602 603 Not collective 604 605 Input Parameter: 606 . viewer - the `PetscViewer` 607 608 Output Parameter: 609 . file_id - The file id 610 611 Level: intermediate 612 613 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()` 614 @*/ 615 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) { 616 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 617 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 620 if (file_id) *file_id = hdf5->file_id; 621 PetscFunctionReturn(0); 622 } 623 624 /*@C 625 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 626 627 Not collective 628 629 Input Parameters: 630 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 631 - name - The group name 632 633 Level: intermediate 634 635 Notes: 636 This is designed to mnemonically resemble the Unix cd command. 637 .vb 638 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. 639 NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 640 "." means the current group is pushed again. 641 .ve 642 643 Example: 644 Suppose the current group is "/a". 645 + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 646 . If name is ".", then the new group will be "/a". 647 . If name is "b", then the new group will be "/a/b". 648 - If name is "/b", then the new group will be "/b". 649 650 Developer Note: 651 The root group "/" is internally stored as NULL. 652 653 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 654 @*/ 655 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) { 656 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 657 PetscViewerHDF5GroupList *groupNode; 658 size_t i, len; 659 char buf[PETSC_MAX_PATH_LEN]; 660 const char *gname; 661 662 PetscFunctionBegin; 663 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 664 if (name) PetscValidCharPointer(name, 2); 665 PetscCall(PetscStrlen(name, &len)); 666 gname = NULL; 667 if (len) { 668 if (len == 1 && name[0] == '.') { 669 /* use current name */ 670 gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 671 } else if (name[0] == '/') { 672 /* absolute */ 673 for (i = 1; i < len; i++) { 674 if (name[i] != '/') { 675 gname = name; 676 break; 677 } 678 } 679 } else { 680 /* relative */ 681 const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 682 PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 683 gname = buf; 684 } 685 } 686 PetscCall(PetscNew(&groupNode)); 687 PetscCall(PetscStrallocpy(gname, (char **)&groupNode->name)); 688 groupNode->next = hdf5->groups; 689 hdf5->groups = groupNode; 690 PetscFunctionReturn(0); 691 } 692 693 /*@ 694 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 695 696 Not collective 697 698 Input Parameter: 699 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 700 701 Level: intermediate 702 703 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5WriteGroup()` 704 @*/ 705 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) { 706 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 707 PetscViewerHDF5GroupList *groupNode; 708 709 PetscFunctionBegin; 710 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 711 PetscCheck(hdf5->groups, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 712 groupNode = hdf5->groups; 713 hdf5->groups = hdf5->groups->next; 714 PetscCall(PetscFree(groupNode->name)); 715 PetscCall(PetscFree(groupNode)); 716 PetscFunctionReturn(0); 717 } 718 719 PetscErrorCode PetscViewerHDF5GetGroup_Internal(PetscViewer viewer, const char *name[]) { 720 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 721 722 PetscFunctionBegin; 723 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 724 PetscValidPointer(name, 2); 725 if (hdf5->groups) *name = hdf5->groups->name; 726 else *name = NULL; 727 PetscFunctionReturn(0); 728 } 729 730 /*@C 731 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by `PetscViewerHDF5GetGroup()`, 732 and return this group's ID and file ID. 733 If `PetscViewerHDF5GetGroup()` yields NULL, then group ID is file ID. 734 735 Not collective 736 737 Input Parameters: 738 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 739 - path - (Optional) The path relative to the pushed group 740 741 Output Parameters: 742 + fileId - The HDF5 file ID 743 - groupId - The HDF5 group ID 744 745 Note: 746 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 747 So NULL or empty path means the current pushed group. 748 749 If the viewer is writable, the group is created if it doesn't exist yet. 750 751 Level: intermediate 752 753 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5WriteGroup()` 754 @*/ 755 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, const char path[], hid_t *fileId, hid_t *groupId) { 756 hid_t file_id; 757 H5O_type_t type; 758 const char *fileName = NULL; 759 char *groupName = NULL; 760 PetscBool writable, has; 761 762 PetscFunctionBegin; 763 PetscCall(PetscViewerWritable(viewer, &writable)); 764 PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 765 PetscCall(PetscViewerFileGetName(viewer, &fileName)); 766 PetscCall(PetscViewerHDF5GetGroup(viewer, path, &groupName)); 767 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 768 if (!has) { 769 PetscCheck(writable, PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName); 770 SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 771 } 772 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); 773 PetscCallHDF5Return(*groupId, H5Gopen2, (file_id, groupName, H5P_DEFAULT)); 774 PetscCall(PetscFree(groupName)); 775 *fileId = file_id; 776 PetscFunctionReturn(0); 777 } 778 779 /*@C 780 PetscViewerHDF5WriteGroup - Ensure the HDF5 group exists in the HDF5 file 781 782 Not collective 783 784 Input Parameters: 785 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 786 - path - (Optional) The path relative to the pushed group 787 788 Note: 789 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 So NULL or empty path means the current pushed group. 791 792 This will fail if the viewer is not writable. 793 794 Level: intermediate 795 796 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 797 @*/ 798 PetscErrorCode PetscViewerHDF5WriteGroup(PetscViewer viewer, const char path[]) { 799 hid_t fileId, groupId; 800 801 PetscFunctionBegin; 802 PetscCall(PetscViewerHDF5OpenGroup(viewer, path, &fileId, &groupId)); // make sure group is actually created 803 PetscCallHDF5(H5Gclose, (groupId)); 804 PetscFunctionReturn(0); 805 } 806 807 /*@ 808 PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 809 810 Not collective 811 812 Input Parameter: 813 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 814 815 Level: intermediate 816 817 Notes: 818 On first `PetscViewerHDF5PushTimestepping()`, the initial time step is set to 0. 819 Next timesteps can then be set using `PetscViewerHDF5IncrementTimestep()` or `PetscViewerHDF5SetTimestep()`. 820 Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. `VecView()`]. 821 Use `PetscViewerHDF5PopTimestepping()` to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 822 Current timestep is remembered between `PetscViewerHDF5PopTimestepping()` and the next `PetscViewerHDF5PushTimestepping()`. 823 824 If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 825 Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 826 827 Developer note: 828 Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 829 830 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 831 @*/ 832 PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) { 833 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 834 835 PetscFunctionBegin; 836 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 837 PetscCheck(!hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 838 hdf5->timestepping = PETSC_TRUE; 839 if (hdf5->timestep < 0) hdf5->timestep = 0; 840 PetscFunctionReturn(0); 841 } 842 843 /*@ 844 PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 845 846 Not collective 847 848 Input Parameter: 849 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 850 851 Level: intermediate 852 853 Note: 854 See `PetscViewerHDF5PushTimestepping()` for details. 855 856 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 857 @*/ 858 PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) { 859 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 860 861 PetscFunctionBegin; 862 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 863 PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 864 hdf5->timestepping = PETSC_FALSE; 865 PetscFunctionReturn(0); 866 } 867 868 /*@ 869 PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 870 871 Not collective 872 873 Input Parameter: 874 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 875 876 Output Parameter: 877 . flg - is timestepping active? 878 879 Level: intermediate 880 881 Note: 882 See `PetscViewerHDF5PushTimestepping()` for details. 883 884 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 885 @*/ 886 PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) { 887 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 888 889 PetscFunctionBegin; 890 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 891 *flg = hdf5->timestepping; 892 PetscFunctionReturn(0); 893 } 894 895 /*@ 896 PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 897 898 Not collective 899 900 Input Parameter: 901 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 902 903 Level: intermediate 904 905 Note: 906 This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 907 908 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 909 @*/ 910 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) { 911 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 912 913 PetscFunctionBegin; 914 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 915 PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 916 ++hdf5->timestep; 917 PetscFunctionReturn(0); 918 } 919 920 /*@ 921 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 922 923 Logically collective 924 925 Input Parameters: 926 + viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 927 - timestep - The timestep 928 929 Level: intermediate 930 931 Note: 932 This can be called only if the viewer is in timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 933 934 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 935 @*/ 936 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) { 937 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 938 939 PetscFunctionBegin; 940 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 941 PetscValidLogicalCollectiveInt(viewer, timestep, 2); 942 PetscCheck(timestep >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 943 PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 944 hdf5->timestep = timestep; 945 PetscFunctionReturn(0); 946 } 947 948 /*@ 949 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 950 951 Not collective 952 953 Input Parameter: 954 . viewer - the `PetscViewer` of type `PETSCVIEWERHDF5` 955 956 Output Parameter: 957 . timestep - The timestep 958 959 Level: intermediate 960 961 Note: 962 This can be called only if the viewer is in the timestepping mode. See `PetscViewerHDF5PushTimestepping()` for details. 963 964 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 965 @*/ 966 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) { 967 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *)viewer->data; 968 969 PetscFunctionBegin; 970 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 971 PetscValidIntPointer(timestep, 2); 972 PetscCheck(hdf5->timestepping, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 973 *timestep = hdf5->timestep; 974 PetscFunctionReturn(0); 975 } 976 977 /*@C 978 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 979 980 Not collective 981 982 Input Parameter: 983 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 984 985 Output Parameter: 986 . mtype - the HDF5 datatype 987 988 Level: advanced 989 990 .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 991 @*/ 992 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) { 993 PetscFunctionBegin; 994 if (ptype == PETSC_INT) 995 #if defined(PETSC_USE_64BIT_INDICES) 996 *htype = H5T_NATIVE_LLONG; 997 #else 998 *htype = H5T_NATIVE_INT; 999 #endif 1000 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 1001 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 1002 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1003 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1004 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 1005 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 1006 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 1007 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 1008 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 1009 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 1010 PetscFunctionReturn(0); 1011 } 1012 1013 /*@C 1014 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 1015 1016 Not collective 1017 1018 Input Parameter: 1019 . htype - the HDF5 datatype (for example `H5T_NATIVE_DOUBLE`, ...) 1020 1021 Output Parameter: 1022 . ptype - the PETSc datatype name (for example `PETSC_DOUBLE`) 1023 1024 Level: advanced 1025 1026 .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 1027 @*/ 1028 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) { 1029 PetscFunctionBegin; 1030 #if defined(PETSC_USE_64BIT_INDICES) 1031 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 1032 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 1033 #else 1034 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 1035 #endif 1036 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 1037 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 1038 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 1039 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 1040 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 1041 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 1042 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 1043 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 1044 PetscFunctionReturn(0); 1045 } 1046 1047 /*@C 1048 PetscViewerHDF5WriteAttribute - Write an attribute 1049 1050 Collective 1051 1052 Input Parameters: 1053 + viewer - The `PETSCVIEWERHDF5` viewer 1054 . parent - The parent dataset/group name 1055 . name - The attribute name 1056 . datatype - The attribute type 1057 - value - The attribute value 1058 1059 Level: advanced 1060 1061 Note: 1062 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. 1063 1064 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, 1065 `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1066 @*/ 1067 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) { 1068 char *parentAbsPath; 1069 hid_t h5, dataspace, obj, attribute, dtype; 1070 PetscBool has; 1071 1072 PetscFunctionBegin; 1073 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1074 if (parent) PetscValidCharPointer(parent, 2); 1075 PetscValidCharPointer(name, 3); 1076 PetscValidLogicalCollectiveEnum(viewer, datatype, 4); 1077 PetscValidPointer(value, 5); 1078 PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 1079 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 1080 PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1081 PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 1082 if (datatype == PETSC_STRING) { 1083 size_t len; 1084 PetscCall(PetscStrlen((const char *)value, &len)); 1085 PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 1086 } 1087 PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1088 PetscCallHDF5Return(dataspace, H5Screate, (H5S_SCALAR)); 1089 PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1090 if (has) { 1091 PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1092 } else { 1093 PetscCallHDF5Return(attribute, H5Acreate2, (obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1094 } 1095 PetscCallHDF5(H5Awrite, (attribute, dtype, value)); 1096 if (datatype == PETSC_STRING) PetscCallHDF5(H5Tclose, (dtype)); 1097 PetscCallHDF5(H5Aclose, (attribute)); 1098 PetscCallHDF5(H5Oclose, (obj)); 1099 PetscCallHDF5(H5Sclose, (dataspace)); 1100 PetscCall(PetscFree(parentAbsPath)); 1101 PetscFunctionReturn(0); 1102 } 1103 1104 /*@C 1105 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given `PetscObject` by name 1106 1107 Collective 1108 1109 Input Parameters: 1110 + viewer - The `PETSCVIEWERHDF5` viewer 1111 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1112 . name - The attribute name 1113 . datatype - The attribute type 1114 - value - The attribute value 1115 1116 Note: 1117 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). 1118 You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1119 1120 Level: advanced 1121 1122 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, 1123 `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1124 @*/ 1125 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) { 1126 PetscFunctionBegin; 1127 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1128 PetscValidHeader(obj, 2); 1129 PetscValidCharPointer(name, 3); 1130 PetscValidPointer(value, 5); 1131 PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 1132 PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 1133 PetscFunctionReturn(0); 1134 } 1135 1136 /*@C 1137 PetscViewerHDF5ReadAttribute - Read an attribute 1138 1139 Collective 1140 1141 Input Parameters: 1142 + viewer - The `PETSCVIEWERHDF5` viewer 1143 . parent - The parent dataset/group name 1144 . name - The attribute name 1145 . datatype - The attribute type 1146 - defaultValue - The pointer to the default value 1147 1148 Output Parameter: 1149 . value - The pointer to the read HDF5 attribute value 1150 1151 Notes: 1152 If defaultValue is NULL and the attribute is not found, an error occurs. 1153 If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1154 The pointers defaultValue and value can be the same; for instance 1155 $ flg = `PETSC_FALSE`; 1156 $ PetscCall(`PetscViewerHDF5ReadAttribute`(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 1157 is valid, but make sure the default value is initialized. 1158 1159 If the datatype is `PETSC_STRING`, the output string is newly allocated so one must `PetscFree()` it when no longer needed. 1160 1161 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. 1162 1163 Level: advanced 1164 1165 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1166 @*/ 1167 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) { 1168 char *parentAbsPath; 1169 hid_t h5, obj, attribute, dtype; 1170 PetscBool has; 1171 1172 PetscFunctionBegin; 1173 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1174 if (parent) PetscValidCharPointer(parent, 2); 1175 PetscValidCharPointer(name, 3); 1176 if (defaultValue) PetscValidPointer(defaultValue, 5); 1177 PetscValidPointer(value, 6); 1178 PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 1179 PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 1180 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 1181 if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1182 if (!has) { 1183 if (defaultValue) { 1184 if (defaultValue != value) { 1185 if (datatype == PETSC_STRING) { 1186 PetscCall(PetscStrallocpy(*(char **)defaultValue, (char **)value)); 1187 } else { 1188 size_t len; 1189 PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (dtype)); 1190 PetscCall(PetscMemcpy(value, defaultValue, len)); 1191 } 1192 } 1193 PetscCall(PetscFree(parentAbsPath)); 1194 PetscFunctionReturn(0); 1195 } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1196 } 1197 PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1198 PetscCallHDF5Return(obj, H5Oopen, (h5, parentAbsPath, H5P_DEFAULT)); 1199 PetscCallHDF5Return(attribute, H5Aopen_name, (obj, name)); 1200 if (datatype == PETSC_STRING) { 1201 size_t len; 1202 hid_t atype; 1203 PetscCallHDF5Return(atype, H5Aget_type, (attribute)); 1204 PetscCallHDF5ReturnNoCheck(len, H5Tget_size, (atype)); 1205 PetscCall(PetscMalloc((len + 1) * sizeof(char), value)); 1206 PetscCallHDF5(H5Tset_size, (dtype, len + 1)); 1207 PetscCallHDF5(H5Aread, (attribute, dtype, *(char **)value)); 1208 } else { 1209 PetscCallHDF5(H5Aread, (attribute, dtype, value)); 1210 } 1211 PetscCallHDF5(H5Aclose, (attribute)); 1212 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1213 PetscCallHDF5(H5Oclose, (obj)); 1214 PetscCall(PetscFree(parentAbsPath)); 1215 PetscFunctionReturn(0); 1216 } 1217 1218 /*@C 1219 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given `PetscObject` by name 1220 1221 Collective 1222 1223 Input Parameters: 1224 + viewer - The `PETSCVIEWERHDF5` viewer 1225 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1226 . name - The attribute name 1227 - datatype - The attribute type 1228 1229 Output Parameter: 1230 . value - The attribute value 1231 1232 Note: 1233 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1234 You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1235 1236 Level: advanced 1237 1238 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1239 @*/ 1240 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) { 1241 PetscFunctionBegin; 1242 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1243 PetscValidHeader(obj, 2); 1244 PetscValidCharPointer(name, 3); 1245 PetscValidPointer(value, 6); 1246 PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 1247 PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 1248 PetscFunctionReturn(0); 1249 } 1250 1251 static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) { 1252 htri_t exists; 1253 hid_t group; 1254 1255 PetscFunctionBegin; 1256 PetscCallHDF5Return(exists, H5Lexists, (h5, name, H5P_DEFAULT)); 1257 if (exists) PetscCallHDF5Return(exists, H5Oexists_by_name, (h5, name, H5P_DEFAULT)); 1258 if (!exists && createGroup) { 1259 PetscCallHDF5Return(group, H5Gcreate2, (h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1260 PetscCallHDF5(H5Gclose, (group)); 1261 exists = PETSC_TRUE; 1262 } 1263 *exists_ = (PetscBool)exists; 1264 PetscFunctionReturn(0); 1265 } 1266 1267 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) { 1268 const char rootGroupName[] = "/"; 1269 hid_t h5; 1270 PetscBool exists = PETSC_FALSE; 1271 PetscInt i; 1272 int n; 1273 char **hierarchy; 1274 char buf[PETSC_MAX_PATH_LEN] = ""; 1275 1276 PetscFunctionBegin; 1277 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1278 if (name) PetscValidCharPointer(name, 2); 1279 else name = rootGroupName; 1280 if (has) { 1281 PetscValidBoolPointer(has, 4); 1282 *has = PETSC_FALSE; 1283 } 1284 if (otype) { 1285 PetscValidIntPointer(otype, 5); 1286 *otype = H5O_TYPE_UNKNOWN; 1287 } 1288 PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1289 1290 /* 1291 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1292 Hence, each of them needs to be tested separately: 1293 1) whether it's a valid link 1294 2) whether this link resolves to an object 1295 See H5Oexists_by_name() documentation. 1296 */ 1297 PetscCall(PetscStrToArray(name, '/', &n, &hierarchy)); 1298 if (!n) { 1299 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1300 if (has) *has = PETSC_TRUE; 1301 if (otype) *otype = H5O_TYPE_GROUP; 1302 PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 1303 PetscFunctionReturn(0); 1304 } 1305 for (i = 0; i < n; i++) { 1306 PetscCall(PetscStrcat(buf, "/")); 1307 PetscCall(PetscStrcat(buf, hierarchy[i])); 1308 PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1309 if (!exists) break; 1310 } 1311 PetscCall(PetscStrToArrayDestroy(n, hierarchy)); 1312 1313 /* If the object exists, get its type */ 1314 if (exists && otype) { 1315 H5O_info_t info; 1316 1317 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1318 PetscCallHDF5(H5Oget_info_by_name, (h5, name, &info, H5P_DEFAULT)); 1319 *otype = info.type; 1320 } 1321 if (has) *has = exists; 1322 PetscFunctionReturn(0); 1323 } 1324 1325 /*@C 1326 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1327 1328 Collective 1329 1330 Input Parameters: 1331 + viewer - The `PETSCVIEWERHDF5` viewer 1332 - path - (Optional) The path relative to the pushed group 1333 1334 Output Parameter: 1335 . has - Flag for group existence 1336 1337 Level: advanced 1338 1339 Notes: 1340 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1341 So NULL or empty path means the current pushed group. 1342 1343 If path exists but is not a group, `PETSC_FALSE` is returned. 1344 1345 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 1346 @*/ 1347 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) { 1348 H5O_type_t type; 1349 char *abspath; 1350 1351 PetscFunctionBegin; 1352 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1353 if (path) PetscValidCharPointer(path, 2); 1354 PetscValidBoolPointer(has, 3); 1355 PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 1356 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 1357 *has = (PetscBool)(type == H5O_TYPE_GROUP); 1358 PetscCall(PetscFree(abspath)); 1359 PetscFunctionReturn(0); 1360 } 1361 1362 /*@C 1363 PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 1364 1365 Collective 1366 1367 Input Parameters: 1368 + viewer - The `PETSCVIEWERHDF5` viewer 1369 - path - The dataset path 1370 1371 Output Parameter: 1372 . has - Flag whether dataset exists 1373 1374 Level: advanced 1375 1376 Notes: 1377 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1378 1379 If path is NULL or empty, has is set to `PETSC_FALSE`. 1380 1381 If path exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1382 1383 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1384 @*/ 1385 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) { 1386 H5O_type_t type; 1387 char *abspath; 1388 1389 PetscFunctionBegin; 1390 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1391 if (path) PetscValidCharPointer(path, 2); 1392 PetscValidBoolPointer(has, 3); 1393 PetscCall(PetscViewerHDF5GetGroup(viewer, path, &abspath)); 1394 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 1395 *has = (PetscBool)(type == H5O_TYPE_DATASET); 1396 PetscCall(PetscFree(abspath)); 1397 PetscFunctionReturn(0); 1398 } 1399 1400 /*@ 1401 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1402 1403 Collective 1404 1405 Input Parameters: 1406 + viewer - The `PETSCVIEWERHDF5` viewer 1407 - obj - The named object 1408 1409 Output Parameter: 1410 . has - Flag for dataset existence 1411 1412 Notes: 1413 If the object is unnamed, an error occurs. 1414 1415 If the path current_group/object_name exists but is not a dataset, has is set to `PETSC_FALSE` as well. 1416 1417 Level: advanced 1418 1419 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1420 @*/ 1421 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) { 1422 size_t len; 1423 1424 PetscFunctionBegin; 1425 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1426 PetscValidHeader(obj, 2); 1427 PetscValidBoolPointer(has, 3); 1428 PetscCall(PetscStrlen(obj->name, &len)); 1429 PetscCheck(len, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 1430 PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 1431 PetscFunctionReturn(0); 1432 } 1433 1434 /*@C 1435 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1436 1437 Collective 1438 1439 Input Parameters: 1440 + viewer - The `PETSCVIEWERHDF5` viewer 1441 . parent - The parent dataset/group name 1442 - name - The attribute name 1443 1444 Output Parameter: 1445 . has - Flag for attribute existence 1446 1447 Level: advanced 1448 1449 Note: 1450 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. 1451 1452 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1453 @*/ 1454 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) { 1455 char *parentAbsPath; 1456 1457 PetscFunctionBegin; 1458 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1459 if (parent) PetscValidCharPointer(parent, 2); 1460 PetscValidCharPointer(name, 3); 1461 PetscValidBoolPointer(has, 4); 1462 PetscCall(PetscViewerHDF5GetGroup(viewer, parent, &parentAbsPath)); 1463 PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 1464 if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 1465 PetscCall(PetscFree(parentAbsPath)); 1466 PetscFunctionReturn(0); 1467 } 1468 1469 /*@C 1470 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given `PetscObject` by name 1471 1472 Collective 1473 1474 Input Parameters: 1475 + viewer - The `PETSCVIEWERHDF5` viewer 1476 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1477 - name - The attribute name 1478 1479 Output Parameter: 1480 . has - Flag for attribute existence 1481 1482 Note: 1483 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1484 You might want to check first if it does using `PetscViewerHDF5HasObject()`. 1485 1486 Level: advanced 1487 1488 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1489 @*/ 1490 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) { 1491 PetscFunctionBegin; 1492 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1493 PetscValidHeader(obj, 2); 1494 PetscValidCharPointer(name, 3); 1495 PetscValidBoolPointer(has, 4); 1496 PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 1497 PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 1498 PetscFunctionReturn(0); 1499 } 1500 1501 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) { 1502 hid_t h5; 1503 htri_t hhas; 1504 1505 PetscFunctionBegin; 1506 PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1507 PetscCallHDF5Return(hhas, H5Aexists_by_name, (h5, parent, name, H5P_DEFAULT)); 1508 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1509 PetscFunctionReturn(0); 1510 } 1511 1512 /* 1513 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1514 is attached to a communicator, in this case the attribute is a PetscViewer. 1515 */ 1516 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1517 1518 /*@C 1519 PETSC_VIEWER_HDF5_ - Creates an `PETSCVIEWERHDF5` `PetscViewer` shared by all processors in a communicator. 1520 1521 Collective 1522 1523 Input Parameter: 1524 . comm - the MPI communicator to share the HDF5 `PetscViewer` 1525 1526 Level: intermediate 1527 1528 Options Database Key: 1529 . -viewer_hdf5_filename <name> - name of the HDF5 file 1530 1531 Environmental variable: 1532 . `PETSC_VIEWER_HDF5_FILENAME` - name of the HDF5 file 1533 1534 Note: 1535 Unlike almost all other PETSc routines, `PETSC_VIEWER_HDF5_()` does not return 1536 an error code. The HDF5 `PetscViewer` is usually used in the form 1537 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1538 1539 .seealso: `PETSCVIEWERHDF5`, `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1540 @*/ 1541 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) { 1542 PetscErrorCode ierr; 1543 PetscBool flg; 1544 PetscViewer viewer; 1545 char fname[PETSC_MAX_PATH_LEN]; 1546 MPI_Comm ncomm; 1547 1548 PetscFunctionBegin; 1549 ierr = PetscCommDuplicate(comm, &ncomm, NULL); 1550 if (ierr) { 1551 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 1552 PetscFunctionReturn(NULL); 1553 } 1554 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1555 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN, MPI_COMM_NULL_DELETE_FN, &Petsc_Viewer_HDF5_keyval, NULL); 1556 if (ierr) { 1557 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 1558 PetscFunctionReturn(NULL); 1559 } 1560 } 1561 ierr = MPI_Comm_get_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void **)&viewer, (int *)&flg); 1562 if (ierr) { 1563 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 1564 PetscFunctionReturn(NULL); 1565 } 1566 if (!flg) { /* PetscViewer not yet created */ 1567 ierr = PetscOptionsGetenv(ncomm, "PETSC_VIEWER_HDF5_FILENAME", fname, PETSC_MAX_PATH_LEN, &flg); 1568 if (ierr) { 1569 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 1570 PetscFunctionReturn(NULL); 1571 } 1572 if (!flg) { 1573 ierr = PetscStrcpy(fname, "output.h5"); 1574 if (ierr) { 1575 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 1576 PetscFunctionReturn(NULL); 1577 } 1578 } 1579 ierr = PetscViewerHDF5Open(ncomm, fname, FILE_MODE_WRITE, &viewer); 1580 if (ierr) { 1581 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 1582 PetscFunctionReturn(NULL); 1583 } 1584 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1585 if (ierr) { 1586 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 1587 PetscFunctionReturn(NULL); 1588 } 1589 ierr = MPI_Comm_set_attr(ncomm, Petsc_Viewer_HDF5_keyval, (void *)viewer); 1590 if (ierr) { 1591 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_INITIAL, " "); 1592 PetscFunctionReturn(NULL); 1593 } 1594 } 1595 ierr = PetscCommDestroy(&ncomm); 1596 if (ierr) { 1597 PetscError(PETSC_COMM_SELF, __LINE__, "PETSC_VIEWER_HDF5_", __FILE__, PETSC_ERR_PLIB, PETSC_ERROR_REPEAT, " "); 1598 PetscFunctionReturn(NULL); 1599 } 1600 PetscFunctionReturn(viewer); 1601 } 1602