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 SETERRQ(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 if (PetscDefined(USE_DEBUG)) { 413 PetscMPIInt rank; 414 PetscBool flg; 415 416 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);CHKERRMPI(ierr); 417 if (rank == 0) { 418 ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 419 PetscCheckFalse(!flg,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"File %s requested for reading does not exist",hdf5->filename); 420 } 421 ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRMPI(ierr); 422 } 423 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 424 break; 425 case FILE_MODE_APPEND: 426 case FILE_MODE_UPDATE: 427 { 428 PetscBool flg; 429 ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 430 if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 431 else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 432 break; 433 } 434 case FILE_MODE_WRITE: 435 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 436 break; 437 case FILE_MODE_UNDEFINED: 438 SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 439 default: 440 SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 441 } 442 PetscCheckFalse(hdf5->file_id < 0,PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 443 PetscStackCallHDF5(H5Pclose,(plist_id)); 444 PetscFunctionReturn(0); 445 } 446 447 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 448 { 449 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 450 451 PetscFunctionBegin; 452 *name = vhdf5->filename; 453 PetscFunctionReturn(0); 454 } 455 456 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 457 { 458 /* 459 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 460 PetscErrorCode ierr; 461 */ 462 463 PetscFunctionBegin; 464 PetscFunctionReturn(0); 465 } 466 467 static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 468 { 469 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 470 471 PetscFunctionBegin; 472 hdf5->defTimestepping = flg; 473 PetscFunctionReturn(0); 474 } 475 476 static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 477 { 478 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 479 480 PetscFunctionBegin; 481 *flg = hdf5->defTimestepping; 482 PetscFunctionReturn(0); 483 } 484 485 /*@ 486 PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 487 488 Logically Collective on PetscViewer 489 490 Input Parameters: 491 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 492 - flg - if PETSC_TRUE we will assume that timestepping is on 493 494 Options Database: 495 . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 496 497 Notes: 498 If the timestepping attribute is not found for an object, then the default timestepping is used 499 500 Level: intermediate 501 502 .seealso: PetscViewerHDF5GetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep() 503 @*/ 504 PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 505 { 506 PetscErrorCode ierr; 507 508 PetscFunctionBegin; 509 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 510 ierr = PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg));CHKERRQ(ierr); 511 PetscFunctionReturn(0); 512 } 513 514 /*@ 515 PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 516 517 Not collective 518 519 Input Parameter: 520 . viewer - the PetscViewer 521 522 Output Parameter: 523 . flg - if PETSC_TRUE we will assume that timestepping is on 524 525 Notes: 526 If the timestepping attribute is not found for an object, then the default timestepping is used 527 528 Level: intermediate 529 530 .seealso: PetscViewerHDF5SetDefaultTimestepping(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5GetTimestep() 531 @*/ 532 PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 533 { 534 PetscErrorCode ierr; 535 536 PetscFunctionBegin; 537 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 538 ierr = PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg));CHKERRQ(ierr); 539 PetscFunctionReturn(0); 540 } 541 542 /*MC 543 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 544 545 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 546 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 547 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 548 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 549 550 Level: beginner 551 M*/ 552 553 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 554 { 555 PetscViewer_HDF5 *hdf5; 556 PetscErrorCode ierr; 557 558 PetscFunctionBegin; 559 #if !defined(H5_HAVE_PARALLEL) 560 { 561 PetscMPIInt size; 562 ierr = MPI_Comm_size(PetscObjectComm((PetscObject)v), &size);CHKERRMPI(ierr); 563 PetscCheckFalse(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)"); 564 } 565 #endif 566 567 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 568 569 v->data = (void*) hdf5; 570 v->ops->destroy = PetscViewerDestroy_HDF5; 571 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 572 v->ops->setup = PetscViewerSetUp_HDF5; 573 v->ops->view = PetscViewerView_HDF5; 574 v->ops->flush = PetscViewerFlush_HDF5; 575 hdf5->btype = FILE_MODE_UNDEFINED; 576 hdf5->filename = NULL; 577 hdf5->timestep = -1; 578 hdf5->groups = NULL; 579 580 PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 581 582 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 583 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 584 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 585 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 586 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 587 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 588 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 589 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 590 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5);CHKERRQ(ierr); 591 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5);CHKERRQ(ierr); 592 PetscFunctionReturn(0); 593 } 594 595 /*@C 596 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 597 598 Collective 599 600 Input Parameters: 601 + comm - MPI communicator 602 . name - name of file 603 - type - type of file 604 605 Output Parameter: 606 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 607 608 Options Database: 609 + -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 610 - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 611 612 Level: beginner 613 614 Notes: 615 Reading is always available, regardless of the mode. Available modes are 616 + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 617 . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 618 . FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL] 619 - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 620 621 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. 622 623 This PetscViewer should be destroyed with PetscViewerDestroy(). 624 625 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 626 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 627 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 628 @*/ 629 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 630 { 631 PetscErrorCode ierr; 632 633 PetscFunctionBegin; 634 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 635 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 636 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 637 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 638 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 639 PetscFunctionReturn(0); 640 } 641 642 /*@C 643 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 644 645 Not collective 646 647 Input Parameter: 648 . viewer - the PetscViewer 649 650 Output Parameter: 651 . file_id - The file id 652 653 Level: intermediate 654 655 .seealso: PetscViewerHDF5Open() 656 @*/ 657 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 658 { 659 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 660 661 PetscFunctionBegin; 662 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 663 if (file_id) *file_id = hdf5->file_id; 664 PetscFunctionReturn(0); 665 } 666 667 /*@C 668 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 669 670 Not collective 671 672 Input Parameters: 673 + viewer - the PetscViewer 674 - name - The group name 675 676 Level: intermediate 677 678 Notes: 679 This is designed to mnemonically resemble the Unix cd command. 680 + 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. 681 . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 682 - "." means the current group is pushed again. 683 684 Example: 685 Suppose the current group is "/a". 686 + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 687 . If name is ".", then the new group will be "/a". 688 . If name is "b", then the new group will be "/a/b". 689 - If name is "/b", then the new group will be "/b". 690 691 Developer Notes: 692 The root group "/" is internally stored as NULL. 693 694 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 695 @*/ 696 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 697 { 698 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 699 PetscViewerHDF5GroupList *groupNode; 700 size_t i,len; 701 char buf[PETSC_MAX_PATH_LEN]; 702 const char *gname; 703 PetscErrorCode ierr; 704 705 PetscFunctionBegin; 706 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 707 if (name) PetscValidCharPointer(name,2); 708 ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 709 gname = NULL; 710 if (len) { 711 if (len == 1 && name[0] == '.') { 712 /* use current name */ 713 gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 714 } else if (name[0] == '/') { 715 /* absolute */ 716 for (i=1; i<len; i++) { 717 if (name[i] != '/') { 718 gname = name; 719 break; 720 } 721 } 722 } else { 723 /* relative */ 724 const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 725 ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 726 gname = buf; 727 } 728 } 729 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 730 ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 731 groupNode->next = hdf5->groups; 732 hdf5->groups = groupNode; 733 PetscFunctionReturn(0); 734 } 735 736 /*@ 737 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 738 739 Not collective 740 741 Input Parameter: 742 . viewer - the PetscViewer 743 744 Level: intermediate 745 746 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 747 @*/ 748 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 749 { 750 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 751 PetscViewerHDF5GroupList *groupNode; 752 PetscErrorCode ierr; 753 754 PetscFunctionBegin; 755 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 756 PetscCheckFalse(!hdf5->groups,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 757 groupNode = hdf5->groups; 758 hdf5->groups = hdf5->groups->next; 759 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 760 ierr = PetscFree(groupNode);CHKERRQ(ierr); 761 PetscFunctionReturn(0); 762 } 763 764 /*@C 765 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 766 If none has been assigned, returns NULL. 767 768 Not collective 769 770 Input Parameter: 771 . viewer - the PetscViewer 772 773 Output Parameter: 774 . name - The group name 775 776 Level: intermediate 777 778 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 779 @*/ 780 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 781 { 782 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 783 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 786 PetscValidPointer(name,2); 787 if (hdf5->groups) *name = hdf5->groups->name; 788 else *name = NULL; 789 PetscFunctionReturn(0); 790 } 791 792 /*@ 793 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 794 and return this group's ID and file ID. 795 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 796 797 Not collective 798 799 Input Parameter: 800 . viewer - the PetscViewer 801 802 Output Parameters: 803 + fileId - The HDF5 file ID 804 - groupId - The HDF5 group ID 805 806 Notes: 807 If the viewer is writable, the group is created if it doesn't exist yet. 808 809 Level: intermediate 810 811 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 812 @*/ 813 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 814 { 815 hid_t file_id; 816 H5O_type_t type; 817 const char *groupName = NULL, *fileName = NULL; 818 PetscBool writable, has; 819 PetscErrorCode ierr; 820 821 PetscFunctionBegin; 822 ierr = PetscViewerWritable(viewer, &writable);CHKERRQ(ierr); 823 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 824 ierr = PetscViewerFileGetName(viewer, &fileName);CHKERRQ(ierr); 825 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 826 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type);CHKERRQ(ierr); 827 if (!has) { 828 PetscCheckFalse(!writable,PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName); 829 else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 830 } 831 PetscCheckFalse(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); 832 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 833 *fileId = file_id; 834 PetscFunctionReturn(0); 835 } 836 837 /*@ 838 PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 839 840 Not collective 841 842 Input Parameter: 843 . viewer - the PetscViewer 844 845 Level: intermediate 846 847 Notes: 848 On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 849 Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 850 Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()]. 851 Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 852 Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 853 854 If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 855 Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 856 857 Developer notes: 858 Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 859 860 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 861 @*/ 862 PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 863 { 864 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 865 866 PetscFunctionBegin; 867 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 868 PetscCheckFalse(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 869 hdf5->timestepping = PETSC_TRUE; 870 if (hdf5->timestep < 0) hdf5->timestep = 0; 871 PetscFunctionReturn(0); 872 } 873 874 /*@ 875 PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 876 877 Not collective 878 879 Input Parameter: 880 . viewer - the PetscViewer 881 882 Level: intermediate 883 884 Notes: 885 See PetscViewerHDF5PushTimestepping() for details. 886 887 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 888 @*/ 889 PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 890 { 891 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 892 893 PetscFunctionBegin; 894 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 895 PetscCheckFalse(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 896 hdf5->timestepping = PETSC_FALSE; 897 PetscFunctionReturn(0); 898 } 899 900 /*@ 901 PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 902 903 Not collective 904 905 Input Parameter: 906 . viewer - the PetscViewer 907 908 Output Parameter: 909 . flg - is timestepping active? 910 911 Level: intermediate 912 913 Notes: 914 See PetscViewerHDF5PushTimestepping() for details. 915 916 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 917 @*/ 918 PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 919 { 920 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 921 922 PetscFunctionBegin; 923 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 924 *flg = hdf5->timestepping; 925 PetscFunctionReturn(0); 926 } 927 928 /*@ 929 PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 930 931 Not collective 932 933 Input Parameter: 934 . viewer - the PetscViewer 935 936 Level: intermediate 937 938 Notes: 939 This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 940 941 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 942 @*/ 943 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 944 { 945 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 946 947 PetscFunctionBegin; 948 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 949 PetscCheckFalse(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 950 ++hdf5->timestep; 951 PetscFunctionReturn(0); 952 } 953 954 /*@ 955 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 956 957 Logically collective 958 959 Input Parameters: 960 + viewer - the PetscViewer 961 - timestep - The timestep 962 963 Level: intermediate 964 965 Notes: 966 This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 967 968 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 969 @*/ 970 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 971 { 972 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 973 974 PetscFunctionBegin; 975 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 976 PetscValidLogicalCollectiveInt(viewer, timestep, 2); 977 PetscCheckFalse(timestep < 0,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 978 PetscCheckFalse(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 979 hdf5->timestep = timestep; 980 PetscFunctionReturn(0); 981 } 982 983 /*@ 984 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 985 986 Not collective 987 988 Input Parameter: 989 . viewer - the PetscViewer 990 991 Output Parameter: 992 . timestep - The timestep 993 994 Level: intermediate 995 996 Notes: 997 This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 998 999 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 1000 @*/ 1001 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 1002 { 1003 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1004 1005 PetscFunctionBegin; 1006 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1007 PetscValidIntPointer(timestep,2); 1008 PetscCheckFalse(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 1009 *timestep = hdf5->timestep; 1010 PetscFunctionReturn(0); 1011 } 1012 1013 /*@C 1014 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 1015 1016 Not collective 1017 1018 Input Parameter: 1019 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 1020 1021 Output Parameter: 1022 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 1023 1024 Level: advanced 1025 1026 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 1027 @*/ 1028 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 1029 { 1030 PetscFunctionBegin; 1031 if (ptype == PETSC_INT) 1032 #if defined(PETSC_USE_64BIT_INDICES) 1033 *htype = H5T_NATIVE_LLONG; 1034 #else 1035 *htype = H5T_NATIVE_INT; 1036 #endif 1037 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 1038 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 1039 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1040 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1041 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 1042 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 1043 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 1044 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 1045 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 1046 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 1047 PetscFunctionReturn(0); 1048 } 1049 1050 /*@C 1051 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 1052 1053 Not collective 1054 1055 Input Parameter: 1056 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 1057 1058 Output Parameter: 1059 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 1060 1061 Level: advanced 1062 1063 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 1064 @*/ 1065 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 1066 { 1067 PetscFunctionBegin; 1068 #if defined(PETSC_USE_64BIT_INDICES) 1069 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 1070 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 1071 #else 1072 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 1073 #endif 1074 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 1075 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 1076 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 1077 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 1078 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 1079 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 1080 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 1081 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 1082 PetscFunctionReturn(0); 1083 } 1084 1085 /*@C 1086 PetscViewerHDF5WriteAttribute - Write an attribute 1087 1088 Collective 1089 1090 Input Parameters: 1091 + viewer - The HDF5 viewer 1092 . parent - The parent dataset/group name 1093 . name - The attribute name 1094 . datatype - The attribute type 1095 - value - The attribute value 1096 1097 Level: advanced 1098 1099 Notes: 1100 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. 1101 1102 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1103 @*/ 1104 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 1105 { 1106 char *parentAbsPath; 1107 hid_t h5, dataspace, obj, attribute, dtype; 1108 PetscBool has; 1109 PetscErrorCode ierr; 1110 1111 PetscFunctionBegin; 1112 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1113 if (parent) PetscValidCharPointer(parent, 2); 1114 PetscValidCharPointer(name, 3); 1115 PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1116 PetscValidPointer(value, 5); 1117 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1118 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 1119 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 1120 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 1121 if (datatype == PETSC_STRING) { 1122 size_t len; 1123 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1124 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1125 } 1126 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1127 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 1128 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1129 if (has) { 1130 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1131 } else { 1132 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1133 } 1134 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1135 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1136 PetscStackCallHDF5(H5Aclose,(attribute)); 1137 PetscStackCallHDF5(H5Oclose,(obj)); 1138 PetscStackCallHDF5(H5Sclose,(dataspace)); 1139 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1140 PetscFunctionReturn(0); 1141 } 1142 1143 /*@C 1144 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1145 1146 Collective 1147 1148 Input Parameters: 1149 + viewer - The HDF5 viewer 1150 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1151 . name - The attribute name 1152 . datatype - The attribute type 1153 - value - The attribute value 1154 1155 Notes: 1156 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). 1157 You might want to check first if it does using PetscViewerHDF5HasObject(). 1158 1159 Level: advanced 1160 1161 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1162 @*/ 1163 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1164 { 1165 PetscErrorCode ierr; 1166 1167 PetscFunctionBegin; 1168 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1169 PetscValidHeader(obj,2); 1170 PetscValidCharPointer(name,3); 1171 PetscValidPointer(value,5); 1172 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1173 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1174 PetscFunctionReturn(0); 1175 } 1176 1177 /*@C 1178 PetscViewerHDF5ReadAttribute - Read an attribute 1179 1180 Collective 1181 1182 Input Parameters: 1183 + viewer - The HDF5 viewer 1184 . parent - The parent dataset/group name 1185 . name - The attribute name 1186 . datatype - The attribute type 1187 - defaultValue - The pointer to the default value 1188 1189 Output Parameter: 1190 . value - The pointer to the read HDF5 attribute value 1191 1192 Notes: 1193 If defaultValue is NULL and the attribute is not found, an error occurs. 1194 If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1195 The pointers defaultValue and value can be the same; for instance 1196 $ flg = PETSC_FALSE; 1197 $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1198 is valid, but make sure the default value is initialized. 1199 1200 If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1201 1202 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. 1203 1204 Level: advanced 1205 1206 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1207 @*/ 1208 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1209 { 1210 char *parentAbsPath; 1211 hid_t h5, obj, attribute, dtype; 1212 PetscBool has; 1213 PetscErrorCode ierr; 1214 1215 PetscFunctionBegin; 1216 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1217 if (parent) PetscValidCharPointer(parent, 2); 1218 PetscValidCharPointer(name, 3); 1219 if (defaultValue) PetscValidPointer(defaultValue, 5); 1220 PetscValidPointer(value, 6); 1221 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 1222 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1223 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 1224 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1225 if (!has) { 1226 if (defaultValue) { 1227 if (defaultValue != value) { 1228 if (datatype == PETSC_STRING) { 1229 ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1230 } else { 1231 size_t len; 1232 PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1233 ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1234 } 1235 } 1236 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1237 PetscFunctionReturn(0); 1238 } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1239 } 1240 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1241 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1242 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1243 if (datatype == PETSC_STRING) { 1244 size_t len; 1245 hid_t atype; 1246 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1247 PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1248 ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1249 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1250 PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1251 } else { 1252 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1253 } 1254 PetscStackCallHDF5(H5Aclose,(attribute)); 1255 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1256 PetscStackCallHDF5(H5Oclose,(obj)); 1257 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1258 PetscFunctionReturn(0); 1259 } 1260 1261 /*@C 1262 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1263 1264 Collective 1265 1266 Input Parameters: 1267 + viewer - The HDF5 viewer 1268 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1269 . name - The attribute name 1270 - datatype - The attribute type 1271 1272 Output Parameter: 1273 . value - The attribute value 1274 1275 Notes: 1276 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1277 You might want to check first if it does using PetscViewerHDF5HasObject(). 1278 1279 Level: advanced 1280 1281 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1282 @*/ 1283 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1284 { 1285 PetscErrorCode ierr; 1286 1287 PetscFunctionBegin; 1288 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1289 PetscValidHeader(obj,2); 1290 PetscValidCharPointer(name,3); 1291 PetscValidPointer(value, 6); 1292 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1293 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1294 PetscFunctionReturn(0); 1295 } 1296 1297 static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1298 { 1299 htri_t exists; 1300 hid_t group; 1301 1302 PetscFunctionBegin; 1303 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1304 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1305 if (!exists && createGroup) { 1306 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1307 PetscStackCallHDF5(H5Gclose,(group)); 1308 exists = PETSC_TRUE; 1309 } 1310 *exists_ = (PetscBool) exists; 1311 PetscFunctionReturn(0); 1312 } 1313 1314 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1315 { 1316 const char rootGroupName[] = "/"; 1317 hid_t h5; 1318 PetscBool exists=PETSC_FALSE; 1319 PetscInt i; 1320 int n; 1321 char **hierarchy; 1322 char buf[PETSC_MAX_PATH_LEN]=""; 1323 PetscErrorCode ierr; 1324 1325 PetscFunctionBegin; 1326 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1327 if (name) PetscValidCharPointer(name, 2); 1328 else name = rootGroupName; 1329 if (has) { 1330 PetscValidBoolPointer(has, 4); 1331 *has = PETSC_FALSE; 1332 } 1333 if (otype) { 1334 PetscValidIntPointer(otype, 5); 1335 *otype = H5O_TYPE_UNKNOWN; 1336 } 1337 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1338 1339 /* 1340 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1341 Hence, each of them needs to be tested separately: 1342 1) whether it's a valid link 1343 2) whether this link resolves to an object 1344 See H5Oexists_by_name() documentation. 1345 */ 1346 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1347 if (!n) { 1348 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1349 if (has) *has = PETSC_TRUE; 1350 if (otype) *otype = H5O_TYPE_GROUP; 1351 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1352 PetscFunctionReturn(0); 1353 } 1354 for (i=0; i<n; i++) { 1355 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1356 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1357 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1358 if (!exists) break; 1359 } 1360 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1361 1362 /* If the object exists, get its type */ 1363 if (exists && otype) { 1364 H5O_info_t info; 1365 1366 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1367 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1368 *otype = info.type; 1369 } 1370 if (has) *has = exists; 1371 PetscFunctionReturn(0); 1372 } 1373 1374 /*@C 1375 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1376 1377 Collective 1378 1379 Input Parameters: 1380 + viewer - The HDF5 viewer 1381 - path - The group path 1382 1383 Output Parameter: 1384 . has - Flag for group existence 1385 1386 Level: advanced 1387 1388 Notes: 1389 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1390 So NULL or empty path means the current pushed group. 1391 If path exists but is not a group, PETSC_FALSE is returned. 1392 1393 .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 1394 @*/ 1395 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 1396 { 1397 H5O_type_t type; 1398 char *abspath; 1399 PetscErrorCode ierr; 1400 1401 PetscFunctionBegin; 1402 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1403 if (path) PetscValidCharPointer(path,2); 1404 PetscValidBoolPointer(has,3); 1405 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 1406 ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 1407 *has = (PetscBool)(type == H5O_TYPE_GROUP); 1408 ierr = PetscFree(abspath);CHKERRQ(ierr); 1409 PetscFunctionReturn(0); 1410 } 1411 1412 /*@C 1413 PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 1414 1415 Collective 1416 1417 Input Parameters: 1418 + viewer - The HDF5 viewer 1419 - path - The dataset path 1420 1421 Output Parameter: 1422 . has - Flag whether dataset exists 1423 1424 Level: advanced 1425 1426 Notes: 1427 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1428 If path is NULL or empty, has is set to PETSC_FALSE. 1429 If path exists but is not a dataset, has is set to PETSC_FALSE as well. 1430 1431 .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1432 @*/ 1433 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 1434 { 1435 H5O_type_t type; 1436 char *abspath; 1437 PetscErrorCode ierr; 1438 1439 PetscFunctionBegin; 1440 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1441 if (path) PetscValidCharPointer(path,2); 1442 PetscValidBoolPointer(has,3); 1443 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 1444 ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 1445 *has = (PetscBool)(type == H5O_TYPE_DATASET); 1446 ierr = PetscFree(abspath);CHKERRQ(ierr); 1447 PetscFunctionReturn(0); 1448 } 1449 1450 /*@ 1451 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1452 1453 Collective 1454 1455 Input Parameters: 1456 + viewer - The HDF5 viewer 1457 - obj - The named object 1458 1459 Output Parameter: 1460 . has - Flag for dataset existence 1461 1462 Notes: 1463 If the object is unnamed, an error occurs. 1464 If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1465 1466 Level: advanced 1467 1468 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1469 @*/ 1470 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1471 { 1472 size_t len; 1473 PetscErrorCode ierr; 1474 1475 PetscFunctionBegin; 1476 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1477 PetscValidHeader(obj,2); 1478 PetscValidBoolPointer(has,3); 1479 ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 1480 PetscCheckFalse(!len,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 1481 ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1482 PetscFunctionReturn(0); 1483 } 1484 1485 /*@C 1486 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1487 1488 Collective 1489 1490 Input Parameters: 1491 + viewer - The HDF5 viewer 1492 . parent - The parent dataset/group name 1493 - name - The attribute name 1494 1495 Output Parameter: 1496 . has - Flag for attribute existence 1497 1498 Level: advanced 1499 1500 Notes: 1501 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. 1502 1503 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1504 @*/ 1505 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1506 { 1507 char *parentAbsPath; 1508 PetscErrorCode ierr; 1509 1510 PetscFunctionBegin; 1511 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1512 if (parent) PetscValidCharPointer(parent,2); 1513 PetscValidCharPointer(name,3); 1514 PetscValidBoolPointer(has,4); 1515 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1516 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1517 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 1518 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1519 PetscFunctionReturn(0); 1520 } 1521 1522 /*@C 1523 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1524 1525 Collective 1526 1527 Input Parameters: 1528 + viewer - The HDF5 viewer 1529 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1530 - name - The attribute name 1531 1532 Output Parameter: 1533 . has - Flag for attribute existence 1534 1535 Notes: 1536 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1537 You might want to check first if it does using PetscViewerHDF5HasObject(). 1538 1539 Level: advanced 1540 1541 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1542 @*/ 1543 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1544 { 1545 PetscErrorCode ierr; 1546 1547 PetscFunctionBegin; 1548 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1549 PetscValidHeader(obj,2); 1550 PetscValidCharPointer(name,3); 1551 PetscValidBoolPointer(has,4); 1552 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1553 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1554 PetscFunctionReturn(0); 1555 } 1556 1557 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1558 { 1559 hid_t h5; 1560 htri_t hhas; 1561 PetscErrorCode ierr; 1562 1563 PetscFunctionBegin; 1564 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1565 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1566 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1567 PetscFunctionReturn(0); 1568 } 1569 1570 /* 1571 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1572 is attached to a communicator, in this case the attribute is a PetscViewer. 1573 */ 1574 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1575 1576 /*@C 1577 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1578 1579 Collective 1580 1581 Input Parameter: 1582 . comm - the MPI communicator to share the HDF5 PetscViewer 1583 1584 Level: intermediate 1585 1586 Options Database Keys: 1587 . -viewer_hdf5_filename <name> - name of the HDF5 file 1588 1589 Environmental variables: 1590 . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1591 1592 Notes: 1593 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1594 an error code. The HDF5 PetscViewer is usually used in the form 1595 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1596 1597 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1598 @*/ 1599 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1600 { 1601 PetscErrorCode ierr; 1602 PetscBool flg; 1603 PetscViewer viewer; 1604 char fname[PETSC_MAX_PATH_LEN]; 1605 MPI_Comm ncomm; 1606 1607 PetscFunctionBegin; 1608 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1609 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1610 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1611 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1612 } 1613 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1614 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1615 if (!flg) { /* PetscViewer not yet created */ 1616 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1617 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1618 if (!flg) { 1619 ierr = PetscStrcpy(fname,"output.h5"); 1620 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1621 } 1622 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1623 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1624 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1625 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1626 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1627 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1628 } 1629 ierr = PetscCommDestroy(&ncomm); 1630 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1631 PetscFunctionReturn(viewer); 1632 } 1633