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