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) && defined(H5_HAVE_PARALLEL) 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 or if built without MPI support\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 and if HDF5 is installed for MPI; 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 #if defined(H5_HAVE_PARALLEL) 332 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 333 H5FD_mpio_xfer_t mode; 334 #endif 335 336 PetscFunctionBegin; 337 #if !defined(H5_HAVE_PARALLEL) 338 *flg = PETSC_FALSE; 339 #else 340 PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 341 *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 342 #endif 343 PetscFunctionReturn(0); 344 } 345 346 /*@ 347 PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 348 349 Not Collective 350 351 Input Parameters: 352 . viewer - the HDF5 PetscViewer 353 354 Output Parameters: 355 . flg - the flag 356 357 Level: intermediate 358 359 Notes: 360 This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned. 361 For more details, see PetscViewerHDF5SetCollective(). 362 363 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 364 365 @*/ 366 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 367 { 368 PetscErrorCode ierr; 369 370 PetscFunctionBegin; 371 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 372 PetscValidBoolPointer(flg,2); 373 374 ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 375 PetscFunctionReturn(0); 376 } 377 378 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 379 { 380 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 381 hid_t plist_id; 382 PetscErrorCode ierr; 383 384 PetscFunctionBegin; 385 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 386 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 387 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 388 /* Set up file access property list with parallel I/O access */ 389 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 390 #if defined(H5_HAVE_PARALLEL) 391 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 392 #endif 393 /* Create or open the file collectively */ 394 switch (hdf5->btype) { 395 case FILE_MODE_READ: 396 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 397 break; 398 case FILE_MODE_APPEND: 399 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 400 break; 401 case FILE_MODE_WRITE: 402 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 403 break; 404 default: 405 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 406 } 407 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 408 PetscStackCallHDF5(H5Pclose,(plist_id)); 409 PetscFunctionReturn(0); 410 } 411 412 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 413 { 414 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 415 416 PetscFunctionBegin; 417 *name = vhdf5->filename; 418 PetscFunctionReturn(0); 419 } 420 421 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 422 { 423 /* 424 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 425 PetscErrorCode ierr; 426 */ 427 428 PetscFunctionBegin; 429 PetscFunctionReturn(0); 430 } 431 432 /*MC 433 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 434 435 436 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 437 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 438 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 439 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 440 441 Level: beginner 442 M*/ 443 444 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 445 { 446 PetscViewer_HDF5 *hdf5; 447 PetscErrorCode ierr; 448 449 PetscFunctionBegin; 450 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 451 452 v->data = (void*) hdf5; 453 v->ops->destroy = PetscViewerDestroy_HDF5; 454 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 455 v->ops->setup = PetscViewerSetUp_HDF5; 456 v->ops->view = PetscViewerView_HDF5; 457 v->ops->flush = 0; 458 hdf5->btype = (PetscFileMode) -1; 459 hdf5->filename = 0; 460 hdf5->timestep = -1; 461 hdf5->groups = NULL; 462 463 PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 464 465 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 466 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 467 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 468 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 469 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 470 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 471 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 472 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 473 PetscFunctionReturn(0); 474 } 475 476 /*@C 477 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 478 479 Collective 480 481 Input Parameters: 482 + comm - MPI communicator 483 . name - name of file 484 - type - type of file 485 $ FILE_MODE_WRITE - create new file for binary output 486 $ FILE_MODE_READ - open existing file for binary input 487 $ FILE_MODE_APPEND - open existing file for binary output 488 489 Output Parameter: 490 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 491 492 Options Database: 493 + -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 494 - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 495 496 Level: beginner 497 498 Note: 499 This PetscViewer should be destroyed with PetscViewerDestroy(). 500 501 502 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 503 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 504 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 505 @*/ 506 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 507 { 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 512 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 513 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 514 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 515 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 516 PetscFunctionReturn(0); 517 } 518 519 /*@C 520 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 521 522 Not collective 523 524 Input Parameter: 525 . viewer - the PetscViewer 526 527 Output Parameter: 528 . file_id - The file id 529 530 Level: intermediate 531 532 .seealso: PetscViewerHDF5Open() 533 @*/ 534 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 535 { 536 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 540 if (file_id) *file_id = hdf5->file_id; 541 PetscFunctionReturn(0); 542 } 543 544 /*@C 545 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 546 547 Not collective 548 549 Input Parameters: 550 + viewer - the PetscViewer 551 - name - The group name 552 553 Level: intermediate 554 555 Note: The group name being NULL, empty string, or a sequence of all slashes (e.g. "///") is always internally stored as NULL and interpreted as "/". 556 557 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 558 @*/ 559 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 560 { 561 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 562 PetscViewerHDF5GroupList *groupNode; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 567 if (name) PetscValidCharPointer(name,2); 568 if (name && name[0]) { 569 size_t i,len; 570 ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 571 for (i=0; i<len; i++) if (name[i] != '/') break; 572 if (i == len) name = NULL; 573 } else name = NULL; 574 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 575 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 576 groupNode->next = hdf5->groups; 577 hdf5->groups = groupNode; 578 PetscFunctionReturn(0); 579 } 580 581 /*@ 582 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 583 584 Not collective 585 586 Input Parameter: 587 . viewer - the PetscViewer 588 589 Level: intermediate 590 591 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 592 @*/ 593 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 594 { 595 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 596 PetscViewerHDF5GroupList *groupNode; 597 PetscErrorCode ierr; 598 599 PetscFunctionBegin; 600 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 601 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 602 groupNode = hdf5->groups; 603 hdf5->groups = hdf5->groups->next; 604 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 605 ierr = PetscFree(groupNode);CHKERRQ(ierr); 606 PetscFunctionReturn(0); 607 } 608 609 /*@C 610 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 611 If none has been assigned, returns NULL. 612 613 Not collective 614 615 Input Parameter: 616 . viewer - the PetscViewer 617 618 Output Parameter: 619 . name - The group name 620 621 Level: intermediate 622 623 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 624 @*/ 625 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 626 { 627 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 628 629 PetscFunctionBegin; 630 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 631 PetscValidPointer(name,2); 632 if (hdf5->groups) *name = hdf5->groups->name; 633 else *name = NULL; 634 PetscFunctionReturn(0); 635 } 636 637 /*@ 638 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 639 and return this group's ID and file ID. 640 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 641 642 Not collective 643 644 Input Parameter: 645 . viewer - the PetscViewer 646 647 Output Parameter: 648 + fileId - The HDF5 file ID 649 - groupId - The HDF5 group ID 650 651 Notes: 652 If the viewer is writable, the group is created if it doesn't exist yet. 653 654 Level: intermediate 655 656 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 657 @*/ 658 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 659 { 660 hid_t file_id; 661 H5O_type_t type; 662 const char *groupName = NULL; 663 PetscBool create; 664 PetscErrorCode ierr; 665 666 PetscFunctionBegin; 667 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 668 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 669 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 670 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 671 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 672 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 673 *fileId = file_id; 674 PetscFunctionReturn(0); 675 } 676 677 /*@ 678 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 679 680 Not collective 681 682 Input Parameter: 683 . viewer - the PetscViewer 684 685 Level: intermediate 686 687 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 688 @*/ 689 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 690 { 691 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 692 693 PetscFunctionBegin; 694 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 695 ++hdf5->timestep; 696 PetscFunctionReturn(0); 697 } 698 699 /*@ 700 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 701 of -1 disables blocking with timesteps. 702 703 Not collective 704 705 Input Parameters: 706 + viewer - the PetscViewer 707 - timestep - The timestep number 708 709 Level: intermediate 710 711 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 712 @*/ 713 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 714 { 715 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 716 717 PetscFunctionBegin; 718 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 719 hdf5->timestep = timestep; 720 PetscFunctionReturn(0); 721 } 722 723 /*@ 724 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 725 726 Not collective 727 728 Input Parameter: 729 . viewer - the PetscViewer 730 731 Output Parameter: 732 . timestep - The timestep number 733 734 Level: intermediate 735 736 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 737 @*/ 738 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 739 { 740 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 741 742 PetscFunctionBegin; 743 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 744 PetscValidPointer(timestep,2); 745 *timestep = hdf5->timestep; 746 PetscFunctionReturn(0); 747 } 748 749 /*@C 750 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 751 752 Not collective 753 754 Input Parameter: 755 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 756 757 Output Parameter: 758 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 759 760 Level: advanced 761 762 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 763 @*/ 764 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 765 { 766 PetscFunctionBegin; 767 if (ptype == PETSC_INT) 768 #if defined(PETSC_USE_64BIT_INDICES) 769 *htype = H5T_NATIVE_LLONG; 770 #else 771 *htype = H5T_NATIVE_INT; 772 #endif 773 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 774 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 775 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 776 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 777 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 778 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 779 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 780 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 781 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 782 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 783 PetscFunctionReturn(0); 784 } 785 786 /*@C 787 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 788 789 Not collective 790 791 Input Parameter: 792 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 793 794 Output Parameter: 795 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 796 797 Level: advanced 798 799 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 800 @*/ 801 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 802 { 803 PetscFunctionBegin; 804 #if defined(PETSC_USE_64BIT_INDICES) 805 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 806 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 807 #else 808 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 809 #endif 810 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 811 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 812 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 813 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 814 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 815 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 816 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 817 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 818 PetscFunctionReturn(0); 819 } 820 821 /*@C 822 PetscViewerHDF5WriteAttribute - Write an attribute 823 824 Input Parameters: 825 + viewer - The HDF5 viewer 826 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 827 . name - The attribute name 828 . datatype - The attribute type 829 - value - The attribute value 830 831 Level: advanced 832 833 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 834 @*/ 835 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 836 { 837 char *parent; 838 hid_t h5, dataspace, obj, attribute, dtype; 839 PetscBool has; 840 PetscErrorCode ierr; 841 842 PetscFunctionBegin; 843 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 844 if (dataset) PetscValidCharPointer(dataset, 2); 845 PetscValidCharPointer(name, 3); 846 PetscValidPointer(value, 5); 847 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 848 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 849 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 850 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 851 if (datatype == PETSC_STRING) { 852 size_t len; 853 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 854 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 855 } 856 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 857 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 858 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 859 if (has) { 860 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 861 } else { 862 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 863 } 864 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 865 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 866 PetscStackCallHDF5(H5Aclose,(attribute)); 867 PetscStackCallHDF5(H5Oclose,(obj)); 868 PetscStackCallHDF5(H5Sclose,(dataspace)); 869 ierr = PetscFree(parent);CHKERRQ(ierr); 870 PetscFunctionReturn(0); 871 } 872 873 /*@C 874 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 875 876 Input Parameters: 877 + viewer - The HDF5 viewer 878 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 879 . name - The attribute name 880 . datatype - The attribute type 881 - value - The attribute value 882 883 Notes: 884 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 885 You might want to check first if it does using PetscViewerHDF5HasObject(). 886 887 Level: advanced 888 889 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 890 @*/ 891 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 892 { 893 PetscErrorCode ierr; 894 895 PetscFunctionBegin; 896 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 897 PetscValidHeader(obj,2); 898 PetscValidCharPointer(name,3); 899 PetscValidPointer(value,5); 900 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 901 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 902 PetscFunctionReturn(0); 903 } 904 905 /*@C 906 PetscViewerHDF5ReadAttribute - Read an attribute 907 908 Input Parameters: 909 + viewer - The HDF5 viewer 910 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 911 . name - The attribute name 912 - datatype - The attribute type 913 914 Output Parameter: 915 . value - The attribute value 916 917 Notes: If the datatype is PETSC_STRING one must PetscFree() the obtained value when it is no longer needed. 918 919 Level: advanced 920 921 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 922 @*/ 923 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 924 { 925 char *parent; 926 hid_t h5, obj, attribute, atype, dtype; 927 PetscBool has; 928 PetscErrorCode ierr; 929 930 PetscFunctionBegin; 931 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 932 if (dataset) PetscValidCharPointer(dataset, 2); 933 PetscValidCharPointer(name, 3); 934 PetscValidPointer(value, 5); 935 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 936 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 937 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 938 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 939 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 940 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 941 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 942 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 943 if (datatype == PETSC_STRING) { 944 size_t len; 945 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 946 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 947 ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 948 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 949 PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 950 } else { 951 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 952 } 953 PetscStackCallHDF5(H5Aclose,(attribute)); 954 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 955 PetscStackCallHDF5(H5Oclose,(obj)); 956 ierr = PetscFree(parent);CHKERRQ(ierr); 957 PetscFunctionReturn(0); 958 } 959 960 /*@C 961 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 962 963 Input Parameters: 964 + viewer - The HDF5 viewer 965 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 966 . name - The attribute name 967 - datatype - The attribute type 968 969 Output Parameter: 970 . value - The attribute value 971 972 Notes: 973 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 974 You might want to check first if it does using PetscViewerHDF5HasObject(). 975 976 Level: advanced 977 978 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 979 @*/ 980 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 981 { 982 PetscErrorCode ierr; 983 984 PetscFunctionBegin; 985 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 986 PetscValidHeader(obj,2); 987 PetscValidCharPointer(name,3); 988 PetscValidPointer(value, 5); 989 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 990 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 991 PetscFunctionReturn(0); 992 } 993 994 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 995 { 996 htri_t exists; 997 hid_t group; 998 999 PetscFunctionBegin; 1000 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1001 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1002 if (!exists && createGroup) { 1003 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1004 PetscStackCallHDF5(H5Gclose,(group)); 1005 exists = PETSC_TRUE; 1006 } 1007 *exists_ = (PetscBool) exists; 1008 PetscFunctionReturn(0); 1009 } 1010 1011 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1012 { 1013 const char rootGroupName[] = "/"; 1014 hid_t h5; 1015 PetscBool exists=PETSC_FALSE; 1016 PetscInt i; 1017 int n; 1018 char **hierarchy; 1019 char buf[PETSC_MAX_PATH_LEN]=""; 1020 PetscErrorCode ierr; 1021 1022 PetscFunctionBegin; 1023 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1024 if (name) PetscValidCharPointer(name, 2); 1025 else name = rootGroupName; 1026 if (has) { 1027 PetscValidIntPointer(has, 3); 1028 *has = PETSC_FALSE; 1029 } 1030 if (otype) { 1031 PetscValidIntPointer(otype, 4); 1032 *otype = H5O_TYPE_UNKNOWN; 1033 } 1034 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1035 1036 /* 1037 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1038 Hence, each of them needs to be tested separately: 1039 1) whether it's a valid link 1040 2) whether this link resolves to an object 1041 See H5Oexists_by_name() documentation. 1042 */ 1043 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1044 if (!n) { 1045 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1046 if (has) *has = PETSC_TRUE; 1047 if (otype) *otype = H5O_TYPE_GROUP; 1048 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1049 PetscFunctionReturn(0); 1050 } 1051 for (i=0; i<n; i++) { 1052 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1053 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1054 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1055 if (!exists) break; 1056 } 1057 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1058 1059 /* If the object exists, get its type */ 1060 if (exists && otype) { 1061 H5O_info_t info; 1062 1063 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1064 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1065 *otype = info.type; 1066 } 1067 if (has) *has = exists; 1068 PetscFunctionReturn(0); 1069 } 1070 1071 /*@ 1072 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1073 1074 Input Parameters: 1075 . viewer - The HDF5 viewer 1076 1077 Output Parameter: 1078 . has - Flag for group existence 1079 1080 Notes: 1081 If the path exists but is not a group, this returns PETSC_FALSE as well. 1082 1083 Level: advanced 1084 1085 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 1086 @*/ 1087 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 1088 { 1089 H5O_type_t type; 1090 const char *name; 1091 PetscErrorCode ierr; 1092 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1095 PetscValidIntPointer(has,2); 1096 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 1097 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1098 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 1099 PetscFunctionReturn(0); 1100 } 1101 1102 /*@ 1103 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1104 1105 Input Parameters: 1106 + viewer - The HDF5 viewer 1107 - obj - The named object 1108 1109 Output Parameter: 1110 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1111 1112 Notes: 1113 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1114 1115 Level: advanced 1116 1117 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1118 @*/ 1119 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1120 { 1121 H5O_type_t type; 1122 char *path; 1123 PetscErrorCode ierr; 1124 1125 PetscFunctionBegin; 1126 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1127 PetscValidHeader(obj,2); 1128 PetscValidIntPointer(has,3); 1129 *has = PETSC_FALSE; 1130 if (!obj->name) PetscFunctionReturn(0); 1131 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1132 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1133 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1134 ierr = PetscFree(path);CHKERRQ(ierr); 1135 PetscFunctionReturn(0); 1136 } 1137 1138 /*@C 1139 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1140 1141 Input Parameters: 1142 + viewer - The HDF5 viewer 1143 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1144 - name - The attribute name 1145 1146 Output Parameter: 1147 . has - Flag for attribute existence 1148 1149 Level: advanced 1150 1151 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1152 @*/ 1153 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1154 { 1155 char *parent; 1156 PetscErrorCode ierr; 1157 1158 PetscFunctionBegin; 1159 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1160 if (dataset) PetscValidCharPointer(dataset,2); 1161 PetscValidCharPointer(name,3); 1162 PetscValidIntPointer(has,4); 1163 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1164 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1165 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1166 ierr = PetscFree(parent);CHKERRQ(ierr); 1167 PetscFunctionReturn(0); 1168 } 1169 1170 /*@C 1171 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1172 1173 Input Parameters: 1174 + viewer - The HDF5 viewer 1175 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1176 - name - The attribute name 1177 1178 Output Parameter: 1179 . has - Flag for attribute existence 1180 1181 Notes: 1182 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1183 You might want to check first if it does using PetscViewerHDF5HasObject(). 1184 1185 Level: advanced 1186 1187 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1188 @*/ 1189 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1190 { 1191 PetscErrorCode ierr; 1192 1193 PetscFunctionBegin; 1194 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1195 PetscValidHeader(obj,2); 1196 PetscValidCharPointer(name,3); 1197 PetscValidIntPointer(has,4); 1198 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1199 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1200 PetscFunctionReturn(0); 1201 } 1202 1203 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1204 { 1205 hid_t h5; 1206 htri_t hhas; 1207 PetscErrorCode ierr; 1208 1209 PetscFunctionBegin; 1210 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1211 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1212 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1213 PetscFunctionReturn(0); 1214 } 1215 1216 /* 1217 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1218 is attached to a communicator, in this case the attribute is a PetscViewer. 1219 */ 1220 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1221 1222 /*@C 1223 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1224 1225 Collective 1226 1227 Input Parameter: 1228 . comm - the MPI communicator to share the HDF5 PetscViewer 1229 1230 Level: intermediate 1231 1232 Options Database Keys: 1233 . -viewer_hdf5_filename <name> 1234 1235 Environmental variables: 1236 . PETSC_VIEWER_HDF5_FILENAME 1237 1238 Notes: 1239 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1240 an error code. The HDF5 PetscViewer is usually used in the form 1241 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1242 1243 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1244 @*/ 1245 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1246 { 1247 PetscErrorCode ierr; 1248 PetscBool flg; 1249 PetscViewer viewer; 1250 char fname[PETSC_MAX_PATH_LEN]; 1251 MPI_Comm ncomm; 1252 1253 PetscFunctionBegin; 1254 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1255 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1256 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1257 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1258 } 1259 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1260 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1261 if (!flg) { /* PetscViewer not yet created */ 1262 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1263 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1264 if (!flg) { 1265 ierr = PetscStrcpy(fname,"output.h5"); 1266 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1267 } 1268 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1269 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1270 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1271 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1272 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1273 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1274 } 1275 ierr = PetscCommDestroy(&ncomm); 1276 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1277 PetscFunctionReturn(viewer); 1278 } 1279