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