1 #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 2 3 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 4 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 5 6 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) 7 { 8 PetscBool relative = PETSC_FALSE; 9 const char *group; 10 char buf[PETSC_MAX_PATH_LEN] = ""; 11 PetscErrorCode ierr; 12 13 PetscFunctionBegin; 14 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 15 ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr); 16 if (relative) { 17 ierr = PetscStrcpy(buf, group);CHKERRQ(ierr); 18 ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 19 ierr = PetscStrcat(buf, path);CHKERRQ(ierr); 20 ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr); 21 } else { 22 ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr); 23 } 24 PetscFunctionReturn(0); 25 } 26 27 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 28 { 29 PetscBool has; 30 PetscErrorCode ierr; 31 32 PetscFunctionBegin; 33 ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 34 if (!has) { 35 const char *group; 36 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 37 SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/"); 38 } 39 PetscFunctionReturn(0); 40 } 41 42 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 43 { 44 PetscErrorCode ierr; 45 PetscBool flg = PETSC_FALSE, set; 46 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 47 48 PetscFunctionBegin; 49 ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 50 ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 51 ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 52 ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 53 if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 54 ierr = PetscOptionsTail();CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 58 static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 59 { 60 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 61 PetscBool flg; 62 PetscErrorCode ierr; 63 64 PetscFunctionBegin; 65 if (hdf5->filename) { 66 ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 67 } 68 ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 69 ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 70 ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 71 ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 72 PetscFunctionReturn(0); 73 } 74 75 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 76 { 77 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 78 PetscErrorCode ierr; 79 80 PetscFunctionBegin; 81 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 82 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 83 PetscFunctionReturn(0); 84 } 85 86 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 87 { 88 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 89 PetscErrorCode ierr; 90 91 PetscFunctionBegin; 92 PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 93 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 94 while (hdf5->groups) { 95 PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 96 97 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 98 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 99 hdf5->groups = tmp; 100 } 101 ierr = PetscFree(hdf5);CHKERRQ(ierr); 102 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 103 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 104 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 105 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 106 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 107 PetscFunctionReturn(0); 108 } 109 110 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 111 { 112 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 113 114 PetscFunctionBegin; 115 hdf5->btype = type; 116 PetscFunctionReturn(0); 117 } 118 119 static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 120 { 121 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 122 123 PetscFunctionBegin; 124 *type = hdf5->btype; 125 PetscFunctionReturn(0); 126 } 127 128 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 129 { 130 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 131 132 PetscFunctionBegin; 133 hdf5->basedimension2 = flg; 134 PetscFunctionReturn(0); 135 } 136 137 /*@ 138 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 139 dimension of 2. 140 141 Logically Collective on PetscViewer 142 143 Input Parameters: 144 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 145 - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 146 147 Options Database: 148 . -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 149 150 Notes: 151 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 152 of one when the dimension is lower. Others think the option is crazy. 153 154 Level: intermediate 155 156 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 157 158 @*/ 159 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 160 { 161 PetscErrorCode ierr; 162 163 PetscFunctionBegin; 164 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 165 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 166 PetscFunctionReturn(0); 167 } 168 169 /*@ 170 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 171 dimension of 2. 172 173 Logically Collective on PetscViewer 174 175 Input Parameter: 176 . viewer - the PetscViewer, must be of type HDF5 177 178 Output Parameter: 179 . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 180 181 Notes: 182 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 183 of one when the dimension is lower. Others think the option is crazy. 184 185 Level: intermediate 186 187 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 188 189 @*/ 190 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 191 { 192 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 193 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 196 *flg = hdf5->basedimension2; 197 PetscFunctionReturn(0); 198 } 199 200 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 201 { 202 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 203 204 PetscFunctionBegin; 205 hdf5->spoutput = flg; 206 PetscFunctionReturn(0); 207 } 208 209 /*@ 210 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 211 compiled with double precision PetscReal. 212 213 Logically Collective on PetscViewer 214 215 Input Parameters: 216 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 217 - flg - if PETSC_TRUE the data will be written to disk with single precision 218 219 Options Database: 220 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 221 222 Notes: 223 Setting this option does not make any difference if PETSc is compiled with single precision 224 in the first place. It does not affect reading datasets (HDF5 handle this internally). 225 226 Level: intermediate 227 228 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 229 PetscReal 230 231 @*/ 232 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 233 { 234 PetscErrorCode ierr; 235 236 PetscFunctionBegin; 237 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 238 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 239 PetscFunctionReturn(0); 240 } 241 242 /*@ 243 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 244 compiled with double precision PetscReal. 245 246 Logically Collective on PetscViewer 247 248 Input Parameter: 249 . viewer - the PetscViewer, must be of type HDF5 250 251 Output Parameter: 252 . flg - if PETSC_TRUE the data will be written to disk with single precision 253 254 Notes: 255 Setting this option does not make any difference if PETSc is compiled with single precision 256 in the first place. It does not affect reading datasets (HDF5 handle this internally). 257 258 Level: intermediate 259 260 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 261 PetscReal 262 263 @*/ 264 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 265 { 266 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 267 268 PetscFunctionBegin; 269 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 270 *flg = hdf5->spoutput; 271 PetscFunctionReturn(0); 272 } 273 274 static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 275 { 276 PetscFunctionBegin; 277 /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 278 - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 279 #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 280 { 281 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 282 PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 283 } 284 #else 285 if (flg) { 286 PetscErrorCode ierr; 287 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); 288 } 289 #endif 290 PetscFunctionReturn(0); 291 } 292 293 /*@ 294 PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 295 296 Logically Collective; flg must contain common value 297 298 Input Parameters: 299 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 300 - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 301 302 Options Database: 303 . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 304 305 Notes: 306 Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 307 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. 308 309 Developer notes: 310 In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 311 This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 312 See HDF5 documentation and MPI-IO documentation for details. 313 314 Level: intermediate 315 316 .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 317 318 @*/ 319 PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 320 { 321 PetscErrorCode ierr; 322 323 PetscFunctionBegin; 324 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 325 PetscValidLogicalCollectiveBool(viewer,flg,2); 326 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 331 { 332 #if defined(H5_HAVE_PARALLEL) 333 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 334 H5FD_mpio_xfer_t mode; 335 #endif 336 337 PetscFunctionBegin; 338 #if !defined(H5_HAVE_PARALLEL) 339 *flg = PETSC_FALSE; 340 #else 341 PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 342 *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 343 #endif 344 PetscFunctionReturn(0); 345 } 346 347 /*@ 348 PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 349 350 Not Collective 351 352 Input Parameters: 353 . viewer - the HDF5 PetscViewer 354 355 Output Parameters: 356 . flg - the flag 357 358 Level: intermediate 359 360 Notes: 361 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. 362 For more details, see PetscViewerHDF5SetCollective(). 363 364 .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 365 366 @*/ 367 PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 368 { 369 PetscErrorCode ierr; 370 371 PetscFunctionBegin; 372 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 373 PetscValidBoolPointer(flg,2); 374 375 ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 376 PetscFunctionReturn(0); 377 } 378 379 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 380 { 381 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 382 hid_t plist_id; 383 PetscErrorCode ierr; 384 385 PetscFunctionBegin; 386 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 387 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 388 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 389 /* Set up file access property list with parallel I/O access */ 390 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 391 #if defined(H5_HAVE_PARALLEL) 392 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 393 #endif 394 /* Create or open the file collectively */ 395 switch (hdf5->btype) { 396 case FILE_MODE_READ: 397 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 398 break; 399 case FILE_MODE_APPEND: 400 case FILE_MODE_UPDATE: 401 { 402 PetscBool flg; 403 ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 404 if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 405 else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 406 break; 407 } 408 case FILE_MODE_WRITE: 409 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 410 break; 411 case FILE_MODE_UNDEFINED: 412 SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 413 default: 414 SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 415 } 416 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 417 PetscStackCallHDF5(H5Pclose,(plist_id)); 418 PetscFunctionReturn(0); 419 } 420 421 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 422 { 423 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 424 425 PetscFunctionBegin; 426 *name = vhdf5->filename; 427 PetscFunctionReturn(0); 428 } 429 430 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 431 { 432 /* 433 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 434 PetscErrorCode ierr; 435 */ 436 437 PetscFunctionBegin; 438 PetscFunctionReturn(0); 439 } 440 441 /*MC 442 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 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 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 515 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 516 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 517 @*/ 518 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 519 { 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 524 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 525 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 526 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 527 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 528 PetscFunctionReturn(0); 529 } 530 531 /*@C 532 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 533 534 Not collective 535 536 Input Parameter: 537 . viewer - the PetscViewer 538 539 Output Parameter: 540 . file_id - The file id 541 542 Level: intermediate 543 544 .seealso: PetscViewerHDF5Open() 545 @*/ 546 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 547 { 548 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 549 550 PetscFunctionBegin; 551 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 552 if (file_id) *file_id = hdf5->file_id; 553 PetscFunctionReturn(0); 554 } 555 556 /*@C 557 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 558 559 Not collective 560 561 Input Parameters: 562 + viewer - the PetscViewer 563 - name - The group name 564 565 Level: intermediate 566 567 Notes: 568 This is designed to mnemonically resemble the Unix cd command. 569 + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group. 570 . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 571 - "." means the current group is pushed again. 572 573 Example: 574 Suppose the current group is "/a". 575 + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 576 . If name is ".", then the new group will be "/a". 577 . If name is "b", then the new group will be "/a/b". 578 - If name is "/b", then the new group will be "/b". 579 580 Developer Notes: 581 The root group "/" is internally stored as NULL. 582 583 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 584 @*/ 585 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 586 { 587 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 588 PetscViewerHDF5GroupList *groupNode; 589 size_t i,len; 590 char buf[PETSC_MAX_PATH_LEN]; 591 const char *gname; 592 PetscErrorCode ierr; 593 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 596 if (name) PetscValidCharPointer(name,2); 597 ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 598 gname = NULL; 599 if (len) { 600 if (len == 1 && name[0] == '.') { 601 /* use current name */ 602 gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 603 } else if (name[0] == '/') { 604 /* absolute */ 605 for (i=1; i<len; i++) { 606 if (name[i] != '/') { 607 gname = name; 608 break; 609 } 610 } 611 } else { 612 /* relative */ 613 const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 614 ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 615 gname = buf; 616 } 617 } 618 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 619 ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 620 groupNode->next = hdf5->groups; 621 hdf5->groups = groupNode; 622 PetscFunctionReturn(0); 623 } 624 625 /*@ 626 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 627 628 Not collective 629 630 Input Parameter: 631 . viewer - the PetscViewer 632 633 Level: intermediate 634 635 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 636 @*/ 637 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 638 { 639 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 640 PetscViewerHDF5GroupList *groupNode; 641 PetscErrorCode ierr; 642 643 PetscFunctionBegin; 644 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 645 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 646 groupNode = hdf5->groups; 647 hdf5->groups = hdf5->groups->next; 648 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 649 ierr = PetscFree(groupNode);CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 /*@C 654 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 655 If none has been assigned, returns NULL. 656 657 Not collective 658 659 Input Parameter: 660 . viewer - the PetscViewer 661 662 Output Parameter: 663 . name - The group name 664 665 Level: intermediate 666 667 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 668 @*/ 669 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 670 { 671 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 672 673 PetscFunctionBegin; 674 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 675 PetscValidPointer(name,2); 676 if (hdf5->groups) *name = hdf5->groups->name; 677 else *name = NULL; 678 PetscFunctionReturn(0); 679 } 680 681 /*@ 682 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 683 and return this group's ID and file ID. 684 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 685 686 Not collective 687 688 Input Parameter: 689 . viewer - the PetscViewer 690 691 Output Parameter: 692 + fileId - The HDF5 file ID 693 - groupId - The HDF5 group ID 694 695 Notes: 696 If the viewer is writable, the group is created if it doesn't exist yet. 697 698 Level: intermediate 699 700 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 701 @*/ 702 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 703 { 704 hid_t file_id; 705 H5O_type_t type; 706 const char *groupName = NULL; 707 PetscBool create; 708 PetscErrorCode ierr; 709 710 PetscFunctionBegin; 711 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 712 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 713 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 714 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 715 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 716 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 717 *fileId = file_id; 718 PetscFunctionReturn(0); 719 } 720 721 /*@ 722 PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 723 724 Not collective 725 726 Input Parameter: 727 . viewer - the PetscViewer 728 729 Level: intermediate 730 731 Notes: 732 On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 733 Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 734 Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()]. 735 Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 736 Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 737 738 If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 739 Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 740 741 Developer notes: 742 Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 743 744 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 745 @*/ 746 PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 747 { 748 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 749 750 PetscFunctionBegin; 751 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 752 if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 753 hdf5->timestepping = PETSC_TRUE; 754 if (hdf5->timestep < 0) hdf5->timestep = 0; 755 PetscFunctionReturn(0); 756 } 757 758 /*@ 759 PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 760 761 Not collective 762 763 Input Parameter: 764 . viewer - the PetscViewer 765 766 Level: intermediate 767 768 Notes: 769 See PetscViewerHDF5PushTimestepping() for details. 770 771 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 772 @*/ 773 PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 774 { 775 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 776 777 PetscFunctionBegin; 778 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 779 if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 780 hdf5->timestepping = PETSC_FALSE; 781 PetscFunctionReturn(0); 782 } 783 784 /*@ 785 PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 786 787 Not collective 788 789 Input Parameter: 790 . viewer - the PetscViewer 791 792 Output Parameter: 793 . flg - is timestepping active? 794 795 Level: intermediate 796 797 Notes: 798 See PetscViewerHDF5PushTimestepping() for details. 799 800 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 801 @*/ 802 PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 803 { 804 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 805 806 PetscFunctionBegin; 807 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 808 *flg = hdf5->timestepping; 809 PetscFunctionReturn(0); 810 } 811 812 /*@ 813 PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 814 815 Not collective 816 817 Input Parameter: 818 . viewer - the PetscViewer 819 820 Level: intermediate 821 822 Notes: 823 This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 824 825 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 826 @*/ 827 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 828 { 829 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 830 831 PetscFunctionBegin; 832 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 833 if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 834 ++hdf5->timestep; 835 PetscFunctionReturn(0); 836 } 837 838 /*@ 839 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 840 841 Logically collective 842 843 Input Parameters: 844 + viewer - the PetscViewer 845 - timestep - The timestep 846 847 Level: intermediate 848 849 Notes: 850 This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 851 852 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 853 @*/ 854 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 855 { 856 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 857 858 PetscFunctionBegin; 859 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 860 PetscValidLogicalCollectiveInt(viewer, timestep, 2); 861 if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep); 862 if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 863 hdf5->timestep = timestep; 864 PetscFunctionReturn(0); 865 } 866 867 /*@ 868 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 869 870 Not collective 871 872 Input Parameter: 873 . viewer - the PetscViewer 874 875 Output Parameter: 876 . timestep - The timestep 877 878 Level: intermediate 879 880 Notes: 881 This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 882 883 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 884 @*/ 885 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 886 { 887 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 888 889 PetscFunctionBegin; 890 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 891 PetscValidIntPointer(timestep,2); 892 if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 893 *timestep = hdf5->timestep; 894 PetscFunctionReturn(0); 895 } 896 897 /*@C 898 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 899 900 Not collective 901 902 Input Parameter: 903 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 904 905 Output Parameter: 906 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 907 908 Level: advanced 909 910 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 911 @*/ 912 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 913 { 914 PetscFunctionBegin; 915 if (ptype == PETSC_INT) 916 #if defined(PETSC_USE_64BIT_INDICES) 917 *htype = H5T_NATIVE_LLONG; 918 #else 919 *htype = H5T_NATIVE_INT; 920 #endif 921 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 922 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 923 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 924 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 925 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 926 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 927 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 928 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 929 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 930 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 931 PetscFunctionReturn(0); 932 } 933 934 /*@C 935 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 936 937 Not collective 938 939 Input Parameter: 940 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 941 942 Output Parameter: 943 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 944 945 Level: advanced 946 947 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 948 @*/ 949 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 950 { 951 PetscFunctionBegin; 952 #if defined(PETSC_USE_64BIT_INDICES) 953 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 954 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 955 #else 956 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 957 #endif 958 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 959 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 960 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 961 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 962 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 963 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 964 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 965 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 966 PetscFunctionReturn(0); 967 } 968 969 /*@C 970 PetscViewerHDF5WriteAttribute - Write an attribute 971 972 Collective 973 974 Input Parameters: 975 + viewer - The HDF5 viewer 976 . parent - The parent dataset/group name 977 . name - The attribute name 978 . datatype - The attribute type 979 - value - The attribute value 980 981 Level: advanced 982 983 Notes: 984 If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 985 986 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 987 @*/ 988 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 989 { 990 char *parentAbsPath; 991 hid_t h5, dataspace, obj, attribute, dtype; 992 PetscBool has; 993 PetscErrorCode ierr; 994 995 PetscFunctionBegin; 996 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 997 if (parent) PetscValidCharPointer(parent, 2); 998 PetscValidCharPointer(name, 3); 999 PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1000 PetscValidPointer(value, 5); 1001 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1002 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 1003 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 1004 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 1005 if (datatype == PETSC_STRING) { 1006 size_t len; 1007 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1008 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1009 } 1010 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1011 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 1012 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1013 if (has) { 1014 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1015 } else { 1016 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1017 } 1018 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1019 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1020 PetscStackCallHDF5(H5Aclose,(attribute)); 1021 PetscStackCallHDF5(H5Oclose,(obj)); 1022 PetscStackCallHDF5(H5Sclose,(dataspace)); 1023 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1024 PetscFunctionReturn(0); 1025 } 1026 1027 /*@C 1028 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1029 1030 Collective 1031 1032 Input Parameters: 1033 + viewer - The HDF5 viewer 1034 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1035 . name - The attribute name 1036 . datatype - The attribute type 1037 - value - The attribute value 1038 1039 Notes: 1040 This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1041 You might want to check first if it does using PetscViewerHDF5HasObject(). 1042 1043 Level: advanced 1044 1045 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1046 @*/ 1047 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1048 { 1049 PetscErrorCode ierr; 1050 1051 PetscFunctionBegin; 1052 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1053 PetscValidHeader(obj,2); 1054 PetscValidCharPointer(name,3); 1055 PetscValidPointer(value,5); 1056 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1057 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1058 PetscFunctionReturn(0); 1059 } 1060 1061 /*@C 1062 PetscViewerHDF5ReadAttribute - Read an attribute 1063 1064 Collective 1065 1066 Input Parameters: 1067 + viewer - The HDF5 viewer 1068 . parent - The parent dataset/group name 1069 . name - The attribute name 1070 . datatype - The attribute type 1071 - defaultValue - The pointer to the default value 1072 1073 Output Parameter: 1074 . value - The pointer to the read HDF5 attribute value 1075 1076 Notes: 1077 If defaultValue is NULL and the attribute is not found, an error occurs. 1078 If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1079 The pointers defaultValue and value can be the same; for instance 1080 $ flg = PETSC_FALSE; 1081 $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1082 is valid, but make sure the default value is initialized. 1083 1084 If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1085 1086 If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 1087 1088 Level: advanced 1089 1090 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1091 @*/ 1092 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1093 { 1094 char *parentAbsPath; 1095 hid_t h5, obj, attribute, dtype; 1096 PetscBool has; 1097 PetscErrorCode ierr; 1098 1099 PetscFunctionBegin; 1100 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1101 if (parent) PetscValidCharPointer(parent, 2); 1102 PetscValidCharPointer(name, 3); 1103 if (defaultValue) PetscValidPointer(defaultValue, 5); 1104 PetscValidPointer(value, 6); 1105 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 1106 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1107 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 1108 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1109 if (!has) { 1110 if (defaultValue) { 1111 if (defaultValue != value) { 1112 if (datatype == PETSC_STRING) { 1113 ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1114 } else { 1115 size_t len; 1116 PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1117 ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1118 } 1119 } 1120 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1121 PetscFunctionReturn(0); 1122 } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1123 } 1124 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1125 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1126 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1127 if (datatype == PETSC_STRING) { 1128 size_t len; 1129 hid_t atype; 1130 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1131 PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1132 ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1133 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1134 PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1135 } else { 1136 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1137 } 1138 PetscStackCallHDF5(H5Aclose,(attribute)); 1139 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1140 PetscStackCallHDF5(H5Oclose,(obj)); 1141 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1142 PetscFunctionReturn(0); 1143 } 1144 1145 /*@C 1146 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1147 1148 Collective 1149 1150 Input Parameters: 1151 + viewer - The HDF5 viewer 1152 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1153 . name - The attribute name 1154 - datatype - The attribute type 1155 1156 Output Parameter: 1157 . value - The attribute value 1158 1159 Notes: 1160 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1161 You might want to check first if it does using PetscViewerHDF5HasObject(). 1162 1163 Level: advanced 1164 1165 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1166 @*/ 1167 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1168 { 1169 PetscErrorCode ierr; 1170 1171 PetscFunctionBegin; 1172 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1173 PetscValidHeader(obj,2); 1174 PetscValidCharPointer(name,3); 1175 PetscValidPointer(value, 6); 1176 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1177 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1178 PetscFunctionReturn(0); 1179 } 1180 1181 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1182 { 1183 htri_t exists; 1184 hid_t group; 1185 1186 PetscFunctionBegin; 1187 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1188 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1189 if (!exists && createGroup) { 1190 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1191 PetscStackCallHDF5(H5Gclose,(group)); 1192 exists = PETSC_TRUE; 1193 } 1194 *exists_ = (PetscBool) exists; 1195 PetscFunctionReturn(0); 1196 } 1197 1198 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1199 { 1200 const char rootGroupName[] = "/"; 1201 hid_t h5; 1202 PetscBool exists=PETSC_FALSE; 1203 PetscInt i; 1204 int n; 1205 char **hierarchy; 1206 char buf[PETSC_MAX_PATH_LEN]=""; 1207 PetscErrorCode ierr; 1208 1209 PetscFunctionBegin; 1210 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1211 if (name) PetscValidCharPointer(name, 2); 1212 else name = rootGroupName; 1213 if (has) { 1214 PetscValidBoolPointer(has, 4); 1215 *has = PETSC_FALSE; 1216 } 1217 if (otype) { 1218 PetscValidIntPointer(otype, 5); 1219 *otype = H5O_TYPE_UNKNOWN; 1220 } 1221 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1222 1223 /* 1224 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1225 Hence, each of them needs to be tested separately: 1226 1) whether it's a valid link 1227 2) whether this link resolves to an object 1228 See H5Oexists_by_name() documentation. 1229 */ 1230 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1231 if (!n) { 1232 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1233 if (has) *has = PETSC_TRUE; 1234 if (otype) *otype = H5O_TYPE_GROUP; 1235 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1236 PetscFunctionReturn(0); 1237 } 1238 for (i=0; i<n; i++) { 1239 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1240 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1241 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1242 if (!exists) break; 1243 } 1244 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1245 1246 /* If the object exists, get its type */ 1247 if (exists && otype) { 1248 H5O_info_t info; 1249 1250 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1251 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1252 *otype = info.type; 1253 } 1254 if (has) *has = exists; 1255 PetscFunctionReturn(0); 1256 } 1257 1258 /*@C 1259 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1260 1261 Collective 1262 1263 Input Parameters: 1264 + viewer - The HDF5 viewer 1265 - path - The group path 1266 1267 Output Parameter: 1268 . has - Flag for group existence 1269 1270 Level: advanced 1271 1272 Notes: 1273 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1274 So NULL or empty path means the current pushed group. 1275 If path exists but is not a group, PETSC_FALSE is returned. 1276 1277 .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 1278 @*/ 1279 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 1280 { 1281 H5O_type_t type; 1282 char *abspath; 1283 PetscErrorCode ierr; 1284 1285 PetscFunctionBegin; 1286 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1287 if (path) PetscValidCharPointer(path,2); 1288 PetscValidBoolPointer(has,3); 1289 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 1290 ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 1291 *has = (PetscBool)(type == H5O_TYPE_GROUP); 1292 ierr = PetscFree(abspath);CHKERRQ(ierr); 1293 PetscFunctionReturn(0); 1294 } 1295 1296 /*@C 1297 PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 1298 1299 Collective 1300 1301 Input Parameters: 1302 + viewer - The HDF5 viewer 1303 - path - The dataset path 1304 1305 Output Parameter: 1306 . has - Flag whether dataset exists 1307 1308 Level: advanced 1309 1310 Notes: 1311 If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1312 If path is NULL or empty, has is set to PETSC_FALSE. 1313 If path exists but is not a dataset, has is set to PETSC_FALSE as well. 1314 1315 .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1316 @*/ 1317 PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 1318 { 1319 H5O_type_t type; 1320 char *abspath; 1321 PetscErrorCode ierr; 1322 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1325 if (path) PetscValidCharPointer(path,2); 1326 PetscValidBoolPointer(has,3); 1327 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 1328 ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 1329 *has = (PetscBool)(type == H5O_TYPE_DATASET); 1330 ierr = PetscFree(abspath);CHKERRQ(ierr); 1331 PetscFunctionReturn(0); 1332 } 1333 1334 /*@ 1335 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1336 1337 Collective 1338 1339 Input Parameters: 1340 + viewer - The HDF5 viewer 1341 - obj - The named object 1342 1343 Output Parameter: 1344 . has - Flag for dataset existence 1345 1346 Notes: 1347 If the object is unnamed, an error occurs. 1348 If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1349 1350 Level: advanced 1351 1352 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1353 @*/ 1354 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1355 { 1356 size_t len; 1357 PetscErrorCode ierr; 1358 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1361 PetscValidHeader(obj,2); 1362 PetscValidBoolPointer(has,3); 1363 ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 1364 if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 1365 ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1366 PetscFunctionReturn(0); 1367 } 1368 1369 /*@C 1370 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1371 1372 Collective 1373 1374 Input Parameters: 1375 + viewer - The HDF5 viewer 1376 . parent - The parent dataset/group name 1377 - name - The attribute name 1378 1379 Output Parameter: 1380 . has - Flag for attribute existence 1381 1382 Level: advanced 1383 1384 Notes: 1385 If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 1386 1387 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1388 @*/ 1389 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1390 { 1391 char *parentAbsPath; 1392 PetscErrorCode ierr; 1393 1394 PetscFunctionBegin; 1395 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1396 if (parent) PetscValidCharPointer(parent,2); 1397 PetscValidCharPointer(name,3); 1398 PetscValidBoolPointer(has,4); 1399 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 1400 ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1401 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 1402 ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1403 PetscFunctionReturn(0); 1404 } 1405 1406 /*@C 1407 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1408 1409 Collective 1410 1411 Input Parameters: 1412 + viewer - The HDF5 viewer 1413 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1414 - name - The attribute name 1415 1416 Output Parameter: 1417 . has - Flag for attribute existence 1418 1419 Notes: 1420 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1421 You might want to check first if it does using PetscViewerHDF5HasObject(). 1422 1423 Level: advanced 1424 1425 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1426 @*/ 1427 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1428 { 1429 PetscErrorCode ierr; 1430 1431 PetscFunctionBegin; 1432 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1433 PetscValidHeader(obj,2); 1434 PetscValidCharPointer(name,3); 1435 PetscValidBoolPointer(has,4); 1436 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1437 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1438 PetscFunctionReturn(0); 1439 } 1440 1441 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1442 { 1443 hid_t h5; 1444 htri_t hhas; 1445 PetscErrorCode ierr; 1446 1447 PetscFunctionBegin; 1448 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1449 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1450 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1451 PetscFunctionReturn(0); 1452 } 1453 1454 /* 1455 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1456 is attached to a communicator, in this case the attribute is a PetscViewer. 1457 */ 1458 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1459 1460 /*@C 1461 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1462 1463 Collective 1464 1465 Input Parameter: 1466 . comm - the MPI communicator to share the HDF5 PetscViewer 1467 1468 Level: intermediate 1469 1470 Options Database Keys: 1471 . -viewer_hdf5_filename <name> 1472 1473 Environmental variables: 1474 . PETSC_VIEWER_HDF5_FILENAME 1475 1476 Notes: 1477 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1478 an error code. The HDF5 PetscViewer is usually used in the form 1479 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1480 1481 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1482 @*/ 1483 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1484 { 1485 PetscErrorCode ierr; 1486 PetscBool flg; 1487 PetscViewer viewer; 1488 char fname[PETSC_MAX_PATH_LEN]; 1489 MPI_Comm ncomm; 1490 1491 PetscFunctionBegin; 1492 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1493 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1494 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1495 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1496 } 1497 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1498 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1499 if (!flg) { /* PetscViewer not yet created */ 1500 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1501 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1502 if (!flg) { 1503 ierr = PetscStrcpy(fname,"output.h5"); 1504 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1505 } 1506 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1507 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1508 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1509 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1510 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1511 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1512 } 1513 ierr = PetscCommDestroy(&ncomm); 1514 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1515 PetscFunctionReturn(viewer); 1516 } 1517