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