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