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