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