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