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