1 #include <petsc/private/viewerimpl.h> 2 #include <petsc/private/viewerhdf5impl.h> 3 #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 4 5 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 6 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 7 8 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 9 { 10 const char *group; 11 char buf[PETSC_MAX_PATH_LEN]=""; 12 PetscErrorCode ierr; 13 14 PetscFunctionBegin; 15 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 16 ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 17 ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 18 ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 19 ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 20 PetscFunctionReturn(0); 21 } 22 23 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 24 { 25 PetscBool has; 26 const char *group; 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 31 ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 32 if (!has) { 33 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 34 SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 35 } 36 PetscFunctionReturn(0); 37 } 38 39 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 40 { 41 PetscErrorCode ierr; 42 PetscBool flg = PETSC_FALSE, set; 43 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 44 45 PetscFunctionBegin; 46 ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 47 ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 48 ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 49 ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 50 if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 51 ierr = PetscOptionsTail();CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 55 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 56 { 57 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 58 PetscBool flg; 59 PetscErrorCode ierr; 60 61 PetscFunctionBegin; 62 if (hdf5->filename) { 63 ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 64 } 65 ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 66 ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 67 ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 68 ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 73 { 74 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 75 PetscErrorCode ierr; 76 77 PetscFunctionBegin; 78 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 79 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 80 PetscFunctionReturn(0); 81 } 82 83 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 84 { 85 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 86 PetscErrorCode ierr; 87 88 PetscFunctionBegin; 89 PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 90 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 91 while (hdf5->groups) { 92 PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 93 94 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 95 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 96 hdf5->groups = tmp; 97 } 98 ierr = PetscFree(hdf5);CHKERRQ(ierr); 99 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 100 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 101 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 102 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 103 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 104 PetscFunctionReturn(0); 105 } 106 107 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 108 { 109 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 110 111 PetscFunctionBegin; 112 hdf5->btype = type; 113 PetscFunctionReturn(0); 114 } 115 116 static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 117 { 118 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 119 120 PetscFunctionBegin; 121 *type = hdf5->btype; 122 PetscFunctionReturn(0); 123 } 124 125 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 126 { 127 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 128 129 PetscFunctionBegin; 130 hdf5->basedimension2 = flg; 131 PetscFunctionReturn(0); 132 } 133 134 /*@ 135 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 136 dimension of 2. 137 138 Logically Collective on PetscViewer 139 140 Input Parameters: 141 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 142 - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 143 144 Options Database: 145 . -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 146 147 148 Notes: 149 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 150 of one when the dimension is lower. Others think the option is crazy. 151 152 Level: intermediate 153 154 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 155 156 @*/ 157 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 158 { 159 PetscErrorCode ierr; 160 161 PetscFunctionBegin; 162 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 163 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 164 PetscFunctionReturn(0); 165 } 166 167 /*@ 168 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 169 dimension of 2. 170 171 Logically Collective on PetscViewer 172 173 Input Parameter: 174 . viewer - the PetscViewer, must be of type HDF5 175 176 Output Parameter: 177 . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 178 179 Notes: 180 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 181 of one when the dimension is lower. Others think the option is crazy. 182 183 Level: intermediate 184 185 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 186 187 @*/ 188 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 189 { 190 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 191 192 PetscFunctionBegin; 193 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 194 *flg = hdf5->basedimension2; 195 PetscFunctionReturn(0); 196 } 197 198 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 199 { 200 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 201 202 PetscFunctionBegin; 203 hdf5->spoutput = flg; 204 PetscFunctionReturn(0); 205 } 206 207 /*@ 208 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 209 compiled with double precision PetscReal. 210 211 Logically Collective on PetscViewer 212 213 Input Parameters: 214 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 215 - flg - if PETSC_TRUE the data will be written to disk with single precision 216 217 Options Database: 218 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 219 220 221 Notes: 222 Setting this option does not make any difference if PETSc is compiled with single precision 223 in the first place. It does not affect reading datasets (HDF5 handle this internally). 224 225 Level: intermediate 226 227 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 228 PetscReal 229 230 @*/ 231 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 232 { 233 PetscErrorCode ierr; 234 235 PetscFunctionBegin; 236 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 237 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 238 PetscFunctionReturn(0); 239 } 240 241 /*@ 242 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 243 compiled with double precision PetscReal. 244 245 Logically Collective on PetscViewer 246 247 Input Parameter: 248 . viewer - the PetscViewer, must be of type HDF5 249 250 Output Parameter: 251 . flg - if PETSC_TRUE the data will be written to disk with single precision 252 253 Notes: 254 Setting this option does not make any difference if PETSc is compiled with single precision 255 in the first place. It does not affect reading datasets (HDF5 handle this internally). 256 257 Level: intermediate 258 259 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 260 PetscReal 261 262 @*/ 263 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 264 { 265 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 266 267 PetscFunctionBegin; 268 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 269 *flg = hdf5->spoutput; 270 PetscFunctionReturn(0); 271 } 272 273 static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 274 { 275 PetscFunctionBegin; 276 /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 277 - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 278 #if H5_VERSION_GE(1,10,3) 279 { 280 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 281 PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 282 } 283 #else 284 if (flg) { 285 PetscErrorCode ierr; 286 ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3\n");CHKERRQ(ierr); 287 } 288 #endif 289 PetscFunctionReturn(0); 290 } 291 292 /*@ 293 PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 294 295 Logically Collective; flg must contain common value 296 297 Input Parameters: 298 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 299 - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 300 301 Options Database: 302 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 303 304 Notes: 305 Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 306 However, this works correctly only since HDF5 1.10.3; hence, we ignore this setting for older versions. 307 308 Developer notes: 309 In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 310 This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 311 See HDF5 documentation and MPI-IO documentation for details. 312 313 Level: intermediate 314 315 .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 316 317 @*/ 318 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 319 { 320 PetscErrorCode ierr; 321 322 PetscFunctionBegin; 323 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 324 PetscValidLogicalCollectiveBool(viewer,flg,2); 325 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 326 PetscFunctionReturn(0); 327 } 328 329 static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 330 { 331 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 332 H5FD_mpio_xfer_t mode; 333 334 PetscFunctionBegin; 335 PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 336 *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 337 PetscFunctionReturn(0); 338 } 339 340 /*@ 341 PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 342 343 Not Collective 344 345 Input Parameters: 346 . viewer - the HDF5 PetscViewer 347 348 Output Parameters: 349 . flg - the flag 350 351 Level: intermediate 352 353 Notes: 354 This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned. 355 For more details, see PetscViewerHDF5SetCollective(). 356 357 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 358 359 @*/ 360 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 361 { 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 366 PetscValidBoolPointer(flg,2); 367 368 ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 373 { 374 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 375 hid_t plist_id; 376 PetscErrorCode ierr; 377 378 PetscFunctionBegin; 379 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 380 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 381 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 382 /* Set up file access property list with parallel I/O access */ 383 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 384 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 385 /* Create or open the file collectively */ 386 switch (hdf5->btype) { 387 case FILE_MODE_READ: 388 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 389 break; 390 case FILE_MODE_APPEND: 391 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 392 break; 393 case FILE_MODE_WRITE: 394 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 395 break; 396 default: 397 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 398 } 399 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 400 PetscStackCallHDF5(H5Pclose,(plist_id)); 401 PetscFunctionReturn(0); 402 } 403 404 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 405 { 406 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 407 408 PetscFunctionBegin; 409 *name = vhdf5->filename; 410 PetscFunctionReturn(0); 411 } 412 413 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 414 { 415 /* 416 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 417 PetscErrorCode ierr; 418 */ 419 420 PetscFunctionBegin; 421 PetscFunctionReturn(0); 422 } 423 424 /*MC 425 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 426 427 428 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 429 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 430 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 431 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 432 433 Level: beginner 434 M*/ 435 436 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 437 { 438 PetscViewer_HDF5 *hdf5; 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 443 444 v->data = (void*) hdf5; 445 v->ops->destroy = PetscViewerDestroy_HDF5; 446 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 447 v->ops->setup = PetscViewerSetUp_HDF5; 448 v->ops->view = PetscViewerView_HDF5; 449 v->ops->flush = 0; 450 hdf5->btype = (PetscFileMode) -1; 451 hdf5->filename = 0; 452 hdf5->timestep = -1; 453 hdf5->groups = NULL; 454 455 PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 456 457 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 458 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 459 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 460 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 461 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 462 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 463 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 464 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 465 PetscFunctionReturn(0); 466 } 467 468 /*@C 469 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 470 471 Collective 472 473 Input Parameters: 474 + comm - MPI communicator 475 . name - name of file 476 - type - type of file 477 $ FILE_MODE_WRITE - create new file for binary output 478 $ FILE_MODE_READ - open existing file for binary input 479 $ FILE_MODE_APPEND - open existing file for binary output 480 481 Output Parameter: 482 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 483 484 Options Database: 485 + -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 486 - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 487 488 Level: beginner 489 490 Note: 491 This PetscViewer should be destroyed with PetscViewerDestroy(). 492 493 494 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 495 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 496 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 497 @*/ 498 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 499 { 500 PetscErrorCode ierr; 501 502 PetscFunctionBegin; 503 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 504 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 505 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 506 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 507 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 508 PetscFunctionReturn(0); 509 } 510 511 /*@C 512 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 513 514 Not collective 515 516 Input Parameter: 517 . viewer - the PetscViewer 518 519 Output Parameter: 520 . file_id - The file id 521 522 Level: intermediate 523 524 .seealso: PetscViewerHDF5Open() 525 @*/ 526 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 527 { 528 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 529 530 PetscFunctionBegin; 531 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 532 if (file_id) *file_id = hdf5->file_id; 533 PetscFunctionReturn(0); 534 } 535 536 /*@C 537 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 538 539 Not collective 540 541 Input Parameters: 542 + viewer - the PetscViewer 543 - name - The group name 544 545 Level: intermediate 546 547 Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/". 548 549 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 550 @*/ 551 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 552 { 553 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 554 PetscViewerHDF5GroupList *groupNode; 555 PetscErrorCode ierr; 556 557 PetscFunctionBegin; 558 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 559 if (name) PetscValidCharPointer(name,2); 560 if (name && name[0]) { 561 size_t i,len; 562 ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 563 for (i=0; i<len; i++) if (name[i] != '/') break; 564 if (i == len) name = NULL; 565 } else name = NULL; 566 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 567 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 568 groupNode->next = hdf5->groups; 569 hdf5->groups = groupNode; 570 PetscFunctionReturn(0); 571 } 572 573 /*@ 574 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 575 576 Not collective 577 578 Input Parameter: 579 . viewer - the PetscViewer 580 581 Level: intermediate 582 583 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 584 @*/ 585 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 586 { 587 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 588 PetscViewerHDF5GroupList *groupNode; 589 PetscErrorCode ierr; 590 591 PetscFunctionBegin; 592 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 593 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 594 groupNode = hdf5->groups; 595 hdf5->groups = hdf5->groups->next; 596 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 597 ierr = PetscFree(groupNode);CHKERRQ(ierr); 598 PetscFunctionReturn(0); 599 } 600 601 /*@C 602 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 603 If none has been assigned, returns NULL. 604 605 Not collective 606 607 Input Parameter: 608 . viewer - the PetscViewer 609 610 Output Parameter: 611 . name - The group name 612 613 Level: intermediate 614 615 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 616 @*/ 617 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 618 { 619 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 620 621 PetscFunctionBegin; 622 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 623 PetscValidPointer(name,2); 624 if (hdf5->groups) *name = hdf5->groups->name; 625 else *name = NULL; 626 PetscFunctionReturn(0); 627 } 628 629 /*@ 630 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 631 and return this group's ID and file ID. 632 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 633 634 Not collective 635 636 Input Parameter: 637 . viewer - the PetscViewer 638 639 Output Parameter: 640 + fileId - The HDF5 file ID 641 - groupId - The HDF5 group ID 642 643 Notes: 644 If the viewer is writable, the group is created if it doesn't exist yet. 645 646 Level: intermediate 647 648 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 649 @*/ 650 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 651 { 652 hid_t file_id; 653 H5O_type_t type; 654 const char *groupName = NULL; 655 PetscBool create; 656 PetscErrorCode ierr; 657 658 PetscFunctionBegin; 659 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 660 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 661 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 662 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 663 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 664 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 665 *fileId = file_id; 666 PetscFunctionReturn(0); 667 } 668 669 /*@ 670 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 671 672 Not collective 673 674 Input Parameter: 675 . viewer - the PetscViewer 676 677 Level: intermediate 678 679 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 680 @*/ 681 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 682 { 683 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 684 685 PetscFunctionBegin; 686 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 687 ++hdf5->timestep; 688 PetscFunctionReturn(0); 689 } 690 691 /*@ 692 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 693 of -1 disables blocking with timesteps. 694 695 Not collective 696 697 Input Parameters: 698 + viewer - the PetscViewer 699 - timestep - The timestep number 700 701 Level: intermediate 702 703 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 704 @*/ 705 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 706 { 707 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 708 709 PetscFunctionBegin; 710 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 711 hdf5->timestep = timestep; 712 PetscFunctionReturn(0); 713 } 714 715 /*@ 716 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 717 718 Not collective 719 720 Input Parameter: 721 . viewer - the PetscViewer 722 723 Output Parameter: 724 . timestep - The timestep number 725 726 Level: intermediate 727 728 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 729 @*/ 730 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 731 { 732 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 736 PetscValidPointer(timestep,2); 737 *timestep = hdf5->timestep; 738 PetscFunctionReturn(0); 739 } 740 741 /*@C 742 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 743 744 Not collective 745 746 Input Parameter: 747 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 748 749 Output Parameter: 750 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 751 752 Level: advanced 753 754 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 755 @*/ 756 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 757 { 758 PetscFunctionBegin; 759 if (ptype == PETSC_INT) 760 #if defined(PETSC_USE_64BIT_INDICES) 761 *htype = H5T_NATIVE_LLONG; 762 #else 763 *htype = H5T_NATIVE_INT; 764 #endif 765 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 766 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 767 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 768 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 769 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 770 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 771 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 772 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 773 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 774 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 775 PetscFunctionReturn(0); 776 } 777 778 /*@C 779 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 780 781 Not collective 782 783 Input Parameter: 784 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 785 786 Output Parameter: 787 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 788 789 Level: advanced 790 791 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 792 @*/ 793 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 794 { 795 PetscFunctionBegin; 796 #if defined(PETSC_USE_64BIT_INDICES) 797 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 798 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 799 #else 800 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 801 #endif 802 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 803 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 804 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 805 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 806 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 807 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 808 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 809 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 810 PetscFunctionReturn(0); 811 } 812 813 /*@C 814 PetscViewerHDF5WriteAttribute - Write an attribute 815 816 Input Parameters: 817 + viewer - The HDF5 viewer 818 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 819 . name - The attribute name 820 . datatype - The attribute type 821 - value - The attribute value 822 823 Level: advanced 824 825 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 826 @*/ 827 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 828 { 829 char *parent; 830 hid_t h5, dataspace, obj, attribute, dtype; 831 PetscBool has; 832 PetscErrorCode ierr; 833 834 PetscFunctionBegin; 835 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 836 if (dataset) PetscValidCharPointer(dataset, 2); 837 PetscValidCharPointer(name, 3); 838 PetscValidPointer(value, 5); 839 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 840 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 841 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 842 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 843 if (datatype == PETSC_STRING) { 844 size_t len; 845 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 846 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 847 } 848 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 849 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 850 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 851 if (has) { 852 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 853 } else { 854 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 855 } 856 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 857 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 858 PetscStackCallHDF5(H5Aclose,(attribute)); 859 PetscStackCallHDF5(H5Oclose,(obj)); 860 PetscStackCallHDF5(H5Sclose,(dataspace)); 861 ierr = PetscFree(parent);CHKERRQ(ierr); 862 PetscFunctionReturn(0); 863 } 864 865 /*@C 866 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 867 868 Input Parameters: 869 + viewer - The HDF5 viewer 870 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 871 . name - The attribute name 872 . datatype - The attribute type 873 - value - The attribute value 874 875 Notes: 876 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 877 You might want to check first if it does using PetscViewerHDF5HasObject(). 878 879 Level: advanced 880 881 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 882 @*/ 883 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 884 { 885 PetscErrorCode ierr; 886 887 PetscFunctionBegin; 888 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 889 PetscValidHeader(obj,2); 890 PetscValidCharPointer(name,3); 891 PetscValidPointer(value,5); 892 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 893 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 894 PetscFunctionReturn(0); 895 } 896 897 /*@C 898 PetscViewerHDF5ReadAttribute - Read an attribute 899 900 Input Parameters: 901 + viewer - The HDF5 viewer 902 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 903 . name - The attribute name 904 - datatype - The attribute type 905 906 Output Parameter: 907 . value - The attribute value 908 909 Level: advanced 910 911 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 912 @*/ 913 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 914 { 915 char *parent; 916 hid_t h5, obj, attribute, atype, dtype; 917 PetscBool has; 918 PetscErrorCode ierr; 919 920 PetscFunctionBegin; 921 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 922 if (dataset) PetscValidCharPointer(dataset, 2); 923 PetscValidCharPointer(name, 3); 924 PetscValidPointer(value, 5); 925 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 926 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 927 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 928 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 929 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 930 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 931 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 932 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 933 if (datatype == PETSC_STRING) { 934 size_t len; 935 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 936 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 937 PetscStackCallHDF5(H5Tclose,(atype)); 938 ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 939 } 940 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 941 PetscStackCallHDF5(H5Aclose,(attribute)); 942 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 943 PetscStackCallHDF5(H5Oclose,(obj)); 944 ierr = PetscFree(parent);CHKERRQ(ierr); 945 PetscFunctionReturn(0); 946 } 947 948 /*@C 949 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 950 951 Input Parameters: 952 + viewer - The HDF5 viewer 953 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 954 . name - The attribute name 955 - datatype - The attribute type 956 957 Output Parameter: 958 . value - The attribute value 959 960 Notes: 961 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 962 You might want to check first if it does using PetscViewerHDF5HasObject(). 963 964 Level: advanced 965 966 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 967 @*/ 968 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 969 { 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 974 PetscValidHeader(obj,2); 975 PetscValidCharPointer(name,3); 976 PetscValidPointer(value, 5); 977 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 978 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 979 PetscFunctionReturn(0); 980 } 981 982 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 983 { 984 htri_t exists; 985 hid_t group; 986 987 PetscFunctionBegin; 988 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 989 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 990 if (!exists && createGroup) { 991 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 992 PetscStackCallHDF5(H5Gclose,(group)); 993 exists = PETSC_TRUE; 994 } 995 *exists_ = (PetscBool) exists; 996 PetscFunctionReturn(0); 997 } 998 999 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1000 { 1001 const char rootGroupName[] = "/"; 1002 hid_t h5; 1003 PetscBool exists=PETSC_FALSE; 1004 PetscInt i; 1005 int n; 1006 char **hierarchy; 1007 char buf[PETSC_MAX_PATH_LEN]=""; 1008 PetscErrorCode ierr; 1009 1010 PetscFunctionBegin; 1011 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1012 if (name) PetscValidCharPointer(name, 2); 1013 else name = rootGroupName; 1014 if (has) { 1015 PetscValidIntPointer(has, 3); 1016 *has = PETSC_FALSE; 1017 } 1018 if (otype) { 1019 PetscValidIntPointer(otype, 4); 1020 *otype = H5O_TYPE_UNKNOWN; 1021 } 1022 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1023 1024 /* 1025 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1026 Hence, each of them needs to be tested separately: 1027 1) whether it's a valid link 1028 2) whether this link resolves to an object 1029 See H5Oexists_by_name() documentation. 1030 */ 1031 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1032 if (!n) { 1033 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1034 if (has) *has = PETSC_TRUE; 1035 if (otype) *otype = H5O_TYPE_GROUP; 1036 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1037 PetscFunctionReturn(0); 1038 } 1039 for (i=0; i<n; i++) { 1040 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1041 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1042 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1043 if (!exists) break; 1044 } 1045 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1046 1047 /* If the object exists, get its type */ 1048 if (exists && otype) { 1049 H5O_info_t info; 1050 1051 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1052 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1053 *otype = info.type; 1054 } 1055 if (has) *has = exists; 1056 PetscFunctionReturn(0); 1057 } 1058 1059 /*@ 1060 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1061 1062 Input Parameters: 1063 . viewer - The HDF5 viewer 1064 1065 Output Parameter: 1066 . has - Flag for group existence 1067 1068 Notes: 1069 If the path exists but is not a group, this returns PETSC_FALSE as well. 1070 1071 Level: advanced 1072 1073 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 1074 @*/ 1075 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 1076 { 1077 H5O_type_t type; 1078 const char *name; 1079 PetscErrorCode ierr; 1080 1081 PetscFunctionBegin; 1082 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1083 PetscValidIntPointer(has,2); 1084 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 1085 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1086 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 1087 PetscFunctionReturn(0); 1088 } 1089 1090 /*@ 1091 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1092 1093 Input Parameters: 1094 + viewer - The HDF5 viewer 1095 - obj - The named object 1096 1097 Output Parameter: 1098 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1099 1100 Notes: 1101 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1102 1103 Level: advanced 1104 1105 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1106 @*/ 1107 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1108 { 1109 H5O_type_t type; 1110 char *path; 1111 PetscErrorCode ierr; 1112 1113 PetscFunctionBegin; 1114 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1115 PetscValidHeader(obj,2); 1116 PetscValidIntPointer(has,3); 1117 *has = PETSC_FALSE; 1118 if (!obj->name) PetscFunctionReturn(0); 1119 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1120 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1121 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1122 ierr = PetscFree(path);CHKERRQ(ierr); 1123 PetscFunctionReturn(0); 1124 } 1125 1126 /*@C 1127 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1128 1129 Input Parameters: 1130 + viewer - The HDF5 viewer 1131 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1132 - name - The attribute name 1133 1134 Output Parameter: 1135 . has - Flag for attribute existence 1136 1137 Level: advanced 1138 1139 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1140 @*/ 1141 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1142 { 1143 char *parent; 1144 PetscErrorCode ierr; 1145 1146 PetscFunctionBegin; 1147 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1148 if (dataset) PetscValidCharPointer(dataset,2); 1149 PetscValidCharPointer(name,3); 1150 PetscValidIntPointer(has,4); 1151 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1152 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1153 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1154 ierr = PetscFree(parent);CHKERRQ(ierr); 1155 PetscFunctionReturn(0); 1156 } 1157 1158 /*@C 1159 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1160 1161 Input Parameters: 1162 + viewer - The HDF5 viewer 1163 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1164 - name - The attribute name 1165 1166 Output Parameter: 1167 . has - Flag for attribute existence 1168 1169 Notes: 1170 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1171 You might want to check first if it does using PetscViewerHDF5HasObject(). 1172 1173 Level: advanced 1174 1175 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1176 @*/ 1177 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1178 { 1179 PetscErrorCode ierr; 1180 1181 PetscFunctionBegin; 1182 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1183 PetscValidHeader(obj,2); 1184 PetscValidCharPointer(name,3); 1185 PetscValidIntPointer(has,4); 1186 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1187 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1188 PetscFunctionReturn(0); 1189 } 1190 1191 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1192 { 1193 hid_t h5; 1194 htri_t hhas; 1195 PetscErrorCode ierr; 1196 1197 PetscFunctionBegin; 1198 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1199 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1200 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1201 PetscFunctionReturn(0); 1202 } 1203 1204 /* 1205 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1206 is attached to a communicator, in this case the attribute is a PetscViewer. 1207 */ 1208 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1209 1210 /*@C 1211 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1212 1213 Collective 1214 1215 Input Parameter: 1216 . comm - the MPI communicator to share the HDF5 PetscViewer 1217 1218 Level: intermediate 1219 1220 Options Database Keys: 1221 . -viewer_hdf5_filename <name> 1222 1223 Environmental variables: 1224 . PETSC_VIEWER_HDF5_FILENAME 1225 1226 Notes: 1227 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1228 an error code. The HDF5 PetscViewer is usually used in the form 1229 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1230 1231 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1232 @*/ 1233 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1234 { 1235 PetscErrorCode ierr; 1236 PetscBool flg; 1237 PetscViewer viewer; 1238 char fname[PETSC_MAX_PATH_LEN]; 1239 MPI_Comm ncomm; 1240 1241 PetscFunctionBegin; 1242 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1243 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1244 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1245 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1246 } 1247 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1248 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1249 if (!flg) { /* PetscViewer not yet created */ 1250 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1251 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1252 if (!flg) { 1253 ierr = PetscStrcpy(fname,"output.h5"); 1254 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1255 } 1256 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1257 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1258 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1259 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1260 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1261 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1262 } 1263 ierr = PetscCommDestroy(&ncomm); 1264 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1265 PetscFunctionReturn(viewer); 1266 } 1267