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 PetscValidPointer(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 on MPI_Comm 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 Concepts: HDF5 files 494 Concepts: PetscViewerHDF5^creating 495 496 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 497 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 498 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 499 @*/ 500 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 501 { 502 PetscErrorCode ierr; 503 504 PetscFunctionBegin; 505 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 506 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 507 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 508 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 509 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 510 PetscFunctionReturn(0); 511 } 512 513 /*@C 514 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 515 516 Not collective 517 518 Input Parameter: 519 . viewer - the PetscViewer 520 521 Output Parameter: 522 . file_id - The file id 523 524 Level: intermediate 525 526 .seealso: PetscViewerHDF5Open() 527 @*/ 528 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 529 { 530 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 531 532 PetscFunctionBegin; 533 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 534 if (file_id) *file_id = hdf5->file_id; 535 PetscFunctionReturn(0); 536 } 537 538 /*@C 539 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 540 541 Not collective 542 543 Input Parameters: 544 + viewer - the PetscViewer 545 - name - The group name 546 547 Level: intermediate 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 PetscValidCharPointer(name,2); 560 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 561 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 562 563 groupNode->next = hdf5->groups; 564 hdf5->groups = groupNode; 565 PetscFunctionReturn(0); 566 } 567 568 /*@ 569 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 570 571 Not collective 572 573 Input Parameter: 574 . viewer - the PetscViewer 575 576 Level: intermediate 577 578 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 579 @*/ 580 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 581 { 582 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 583 PetscViewerHDF5GroupList *groupNode; 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 588 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 589 groupNode = hdf5->groups; 590 hdf5->groups = hdf5->groups->next; 591 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 592 ierr = PetscFree(groupNode);CHKERRQ(ierr); 593 PetscFunctionReturn(0); 594 } 595 596 /*@C 597 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 598 If none has been assigned, returns NULL. 599 600 Not collective 601 602 Input Parameter: 603 . viewer - the PetscViewer 604 605 Output Parameter: 606 . name - The group name 607 608 Level: intermediate 609 610 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 611 @*/ 612 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 613 { 614 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 615 616 PetscFunctionBegin; 617 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 618 PetscValidPointer(name,2); 619 if (hdf5->groups) *name = hdf5->groups->name; 620 else *name = NULL; 621 PetscFunctionReturn(0); 622 } 623 624 /*@ 625 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 626 and return this group's ID and file ID. 627 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 628 629 Not collective 630 631 Input Parameter: 632 . viewer - the PetscViewer 633 634 Output Parameter: 635 + fileId - The HDF5 file ID 636 - groupId - The HDF5 group ID 637 638 Notes: 639 If the viewer is writable, the group is created if it doesn't exist yet. 640 641 Level: intermediate 642 643 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 644 @*/ 645 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 646 { 647 hid_t file_id; 648 H5O_type_t type; 649 const char *groupName = NULL; 650 PetscBool create; 651 PetscErrorCode ierr; 652 653 PetscFunctionBegin; 654 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 655 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 656 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 657 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 658 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 659 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 660 *fileId = file_id; 661 PetscFunctionReturn(0); 662 } 663 664 /*@ 665 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 666 667 Not collective 668 669 Input Parameter: 670 . viewer - the PetscViewer 671 672 Level: intermediate 673 674 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 675 @*/ 676 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 677 { 678 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 679 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 682 ++hdf5->timestep; 683 PetscFunctionReturn(0); 684 } 685 686 /*@ 687 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 688 of -1 disables blocking with timesteps. 689 690 Not collective 691 692 Input Parameters: 693 + viewer - the PetscViewer 694 - timestep - The timestep number 695 696 Level: intermediate 697 698 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 699 @*/ 700 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 701 { 702 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 703 704 PetscFunctionBegin; 705 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 706 hdf5->timestep = timestep; 707 PetscFunctionReturn(0); 708 } 709 710 /*@ 711 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 712 713 Not collective 714 715 Input Parameter: 716 . viewer - the PetscViewer 717 718 Output Parameter: 719 . timestep - The timestep number 720 721 Level: intermediate 722 723 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 724 @*/ 725 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 726 { 727 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 728 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 731 PetscValidPointer(timestep,2); 732 *timestep = hdf5->timestep; 733 PetscFunctionReturn(0); 734 } 735 736 /*@C 737 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 738 739 Not collective 740 741 Input Parameter: 742 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 743 744 Output Parameter: 745 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 746 747 Level: advanced 748 749 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 750 @*/ 751 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 752 { 753 PetscFunctionBegin; 754 if (ptype == PETSC_INT) 755 #if defined(PETSC_USE_64BIT_INDICES) 756 *htype = H5T_NATIVE_LLONG; 757 #else 758 *htype = H5T_NATIVE_INT; 759 #endif 760 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 761 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 762 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 763 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 764 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 765 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 766 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 767 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 768 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 769 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 770 PetscFunctionReturn(0); 771 } 772 773 /*@C 774 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 775 776 Not collective 777 778 Input Parameter: 779 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 780 781 Output Parameter: 782 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 783 784 Level: advanced 785 786 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 787 @*/ 788 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 789 { 790 PetscFunctionBegin; 791 #if defined(PETSC_USE_64BIT_INDICES) 792 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 793 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 794 #else 795 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 796 #endif 797 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 798 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 799 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 800 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 801 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 802 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 803 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 804 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 805 PetscFunctionReturn(0); 806 } 807 808 /*@C 809 PetscViewerHDF5WriteAttribute - Write an attribute 810 811 Input Parameters: 812 + viewer - The HDF5 viewer 813 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 814 . name - The attribute name 815 . datatype - The attribute type 816 - value - The attribute value 817 818 Level: advanced 819 820 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 821 @*/ 822 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 823 { 824 char *parent; 825 hid_t h5, dataspace, obj, attribute, dtype; 826 PetscBool has; 827 PetscErrorCode ierr; 828 829 PetscFunctionBegin; 830 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 831 if (dataset) PetscValidCharPointer(dataset, 2); 832 PetscValidCharPointer(name, 3); 833 PetscValidPointer(value, 5); 834 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 835 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 836 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 837 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 838 if (datatype == PETSC_STRING) { 839 size_t len; 840 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 841 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 842 } 843 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 844 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 845 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 846 if (has) { 847 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 848 } else { 849 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 850 } 851 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 852 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 853 PetscStackCallHDF5(H5Aclose,(attribute)); 854 PetscStackCallHDF5(H5Oclose,(obj)); 855 PetscStackCallHDF5(H5Sclose,(dataspace)); 856 ierr = PetscFree(parent);CHKERRQ(ierr); 857 PetscFunctionReturn(0); 858 } 859 860 /*@C 861 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 862 863 Input Parameters: 864 + viewer - The HDF5 viewer 865 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 866 . name - The attribute name 867 . datatype - The attribute type 868 - value - The attribute value 869 870 Notes: 871 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 872 You might want to check first if it does using PetscViewerHDF5HasObject(). 873 874 Level: advanced 875 876 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 877 @*/ 878 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 879 { 880 PetscErrorCode ierr; 881 882 PetscFunctionBegin; 883 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 884 PetscValidHeader(obj,2); 885 PetscValidCharPointer(name,3); 886 PetscValidPointer(value,5); 887 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 888 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 889 PetscFunctionReturn(0); 890 } 891 892 /*@C 893 PetscViewerHDF5ReadAttribute - Read an attribute 894 895 Input Parameters: 896 + viewer - The HDF5 viewer 897 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 898 . name - The attribute name 899 - datatype - The attribute type 900 901 Output Parameter: 902 . value - The attribute value 903 904 Level: advanced 905 906 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 907 @*/ 908 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 909 { 910 char *parent; 911 hid_t h5, obj, attribute, atype, dtype; 912 PetscBool has; 913 PetscErrorCode ierr; 914 915 PetscFunctionBegin; 916 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 917 if (dataset) PetscValidCharPointer(dataset, 2); 918 PetscValidCharPointer(name, 3); 919 PetscValidPointer(value, 5); 920 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 921 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 922 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 923 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 924 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 925 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 926 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 927 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 928 if (datatype == PETSC_STRING) { 929 size_t len; 930 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 931 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 932 PetscStackCallHDF5(H5Tclose,(atype)); 933 ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 934 } 935 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 936 PetscStackCallHDF5(H5Aclose,(attribute)); 937 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 938 PetscStackCallHDF5(H5Oclose,(obj)); 939 ierr = PetscFree(parent);CHKERRQ(ierr); 940 PetscFunctionReturn(0); 941 } 942 943 /*@C 944 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 945 946 Input Parameters: 947 + viewer - The HDF5 viewer 948 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 949 . name - The attribute name 950 - datatype - The attribute type 951 952 Output Parameter: 953 . value - The attribute value 954 955 Notes: 956 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 957 You might want to check first if it does using PetscViewerHDF5HasObject(). 958 959 Level: advanced 960 961 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 962 @*/ 963 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 964 { 965 PetscErrorCode ierr; 966 967 PetscFunctionBegin; 968 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 969 PetscValidHeader(obj,2); 970 PetscValidCharPointer(name,3); 971 PetscValidPointer(value, 5); 972 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 973 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 974 PetscFunctionReturn(0); 975 } 976 977 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 978 { 979 htri_t exists; 980 hid_t group; 981 982 PetscFunctionBegin; 983 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 984 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 985 if (!exists && createGroup) { 986 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 987 PetscStackCallHDF5(H5Gclose,(group)); 988 exists = PETSC_TRUE; 989 } 990 *exists_ = (PetscBool) exists; 991 PetscFunctionReturn(0); 992 } 993 994 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 995 { 996 const char rootGroupName[] = "/"; 997 hid_t h5; 998 PetscBool exists=PETSC_FALSE; 999 PetscInt i; 1000 int n; 1001 char **hierarchy; 1002 char buf[PETSC_MAX_PATH_LEN]=""; 1003 PetscErrorCode ierr; 1004 1005 PetscFunctionBegin; 1006 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1007 if (name) PetscValidCharPointer(name, 2); 1008 else name = rootGroupName; 1009 if (has) { 1010 PetscValidIntPointer(has, 3); 1011 *has = PETSC_FALSE; 1012 } 1013 if (otype) { 1014 PetscValidIntPointer(otype, 4); 1015 *otype = H5O_TYPE_UNKNOWN; 1016 } 1017 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1018 1019 /* 1020 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1021 Hence, each of them needs to be tested separately: 1022 1) whether it's a valid link 1023 2) whether this link resolves to an object 1024 See H5Oexists_by_name() documentation. 1025 */ 1026 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1027 if (!n) { 1028 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1029 if (has) *has = PETSC_TRUE; 1030 if (otype) *otype = H5O_TYPE_GROUP; 1031 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1032 PetscFunctionReturn(0); 1033 } 1034 for (i=0; i<n; i++) { 1035 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1036 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1037 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1038 if (!exists) break; 1039 } 1040 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1041 1042 /* If the object exists, get its type */ 1043 if (exists && otype) { 1044 H5O_info_t info; 1045 1046 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1047 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1048 *otype = info.type; 1049 } 1050 if (has) *has = exists; 1051 PetscFunctionReturn(0); 1052 } 1053 1054 /*@ 1055 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1056 1057 Input Parameters: 1058 . viewer - The HDF5 viewer 1059 1060 Output Parameter: 1061 . has - Flag for group existence 1062 1063 Notes: 1064 If the path exists but is not a group, this returns PETSC_FALSE as well. 1065 1066 Level: advanced 1067 1068 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 1069 @*/ 1070 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 1071 { 1072 H5O_type_t type; 1073 const char *name; 1074 PetscErrorCode ierr; 1075 1076 PetscFunctionBegin; 1077 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1078 PetscValidIntPointer(has,2); 1079 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 1080 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1081 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 1082 PetscFunctionReturn(0); 1083 } 1084 1085 /*@ 1086 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1087 1088 Input Parameters: 1089 + viewer - The HDF5 viewer 1090 - obj - The named object 1091 1092 Output Parameter: 1093 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1094 1095 Notes: 1096 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1097 1098 Level: advanced 1099 1100 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1101 @*/ 1102 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1103 { 1104 H5O_type_t type; 1105 char *path; 1106 PetscErrorCode ierr; 1107 1108 PetscFunctionBegin; 1109 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1110 PetscValidHeader(obj,2); 1111 PetscValidIntPointer(has,3); 1112 *has = PETSC_FALSE; 1113 if (!obj->name) PetscFunctionReturn(0); 1114 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1115 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1116 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1117 ierr = PetscFree(path);CHKERRQ(ierr); 1118 PetscFunctionReturn(0); 1119 } 1120 1121 /*@C 1122 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1123 1124 Input Parameters: 1125 + viewer - The HDF5 viewer 1126 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1127 - name - The attribute name 1128 1129 Output Parameter: 1130 . has - Flag for attribute existence 1131 1132 Level: advanced 1133 1134 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1135 @*/ 1136 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1137 { 1138 char *parent; 1139 PetscErrorCode ierr; 1140 1141 PetscFunctionBegin; 1142 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1143 if (dataset) PetscValidCharPointer(dataset,2); 1144 PetscValidCharPointer(name,3); 1145 PetscValidIntPointer(has,4); 1146 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1147 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1148 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1149 ierr = PetscFree(parent);CHKERRQ(ierr); 1150 PetscFunctionReturn(0); 1151 } 1152 1153 /*@C 1154 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1155 1156 Input Parameters: 1157 + viewer - The HDF5 viewer 1158 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1159 - name - The attribute name 1160 1161 Output Parameter: 1162 . has - Flag for attribute existence 1163 1164 Notes: 1165 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1166 You might want to check first if it does using PetscViewerHDF5HasObject(). 1167 1168 Level: advanced 1169 1170 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1171 @*/ 1172 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1173 { 1174 PetscErrorCode ierr; 1175 1176 PetscFunctionBegin; 1177 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1178 PetscValidHeader(obj,2); 1179 PetscValidCharPointer(name,3); 1180 PetscValidIntPointer(has,4); 1181 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1182 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1183 PetscFunctionReturn(0); 1184 } 1185 1186 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1187 { 1188 hid_t h5; 1189 htri_t hhas; 1190 PetscErrorCode ierr; 1191 1192 PetscFunctionBegin; 1193 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1194 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1195 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1196 PetscFunctionReturn(0); 1197 } 1198 1199 /* 1200 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1201 is attached to a communicator, in this case the attribute is a PetscViewer. 1202 */ 1203 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1204 1205 /*@C 1206 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1207 1208 Collective on MPI_Comm 1209 1210 Input Parameter: 1211 . comm - the MPI communicator to share the HDF5 PetscViewer 1212 1213 Level: intermediate 1214 1215 Options Database Keys: 1216 . -viewer_hdf5_filename <name> 1217 1218 Environmental variables: 1219 . PETSC_VIEWER_HDF5_FILENAME 1220 1221 Notes: 1222 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1223 an error code. The HDF5 PetscViewer is usually used in the form 1224 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1225 1226 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1227 @*/ 1228 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1229 { 1230 PetscErrorCode ierr; 1231 PetscBool flg; 1232 PetscViewer viewer; 1233 char fname[PETSC_MAX_PATH_LEN]; 1234 MPI_Comm ncomm; 1235 1236 PetscFunctionBegin; 1237 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1238 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1239 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1240 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1241 } 1242 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1243 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1244 if (!flg) { /* PetscViewer not yet created */ 1245 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1246 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1247 if (!flg) { 1248 ierr = PetscStrcpy(fname,"output.h5"); 1249 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1250 } 1251 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1252 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1253 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1254 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1255 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1256 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1257 } 1258 ierr = PetscCommDestroy(&ncomm); 1259 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1260 PetscFunctionReturn(viewer); 1261 } 1262