1 #include <petsc/private/viewerimpl.h> 2 #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 3 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800) 4 #error "PETSc needs HDF5 version >= 1.8.0" 5 #endif 6 7 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 8 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 9 10 typedef struct GroupList { 11 const char *name; 12 struct GroupList *next; 13 } GroupList; 14 15 typedef struct { 16 char *filename; 17 PetscFileMode btype; 18 hid_t file_id; 19 PetscInt timestep; 20 GroupList *groups; 21 PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 22 PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 23 } PetscViewer_HDF5; 24 25 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 26 { 27 const char *group; 28 char buf[PETSC_MAX_PATH_LEN]=""; 29 PetscErrorCode ierr; 30 31 PetscFunctionBegin; 32 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 33 ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 34 ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 35 ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 36 ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 37 PetscFunctionReturn(0); 38 } 39 40 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 41 { 42 PetscBool has; 43 const char *group; 44 PetscErrorCode ierr; 45 46 PetscFunctionBegin; 47 if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 48 ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 49 if (!has) { 50 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 51 SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 52 } 53 PetscFunctionReturn(0); 54 } 55 56 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 57 { 58 PetscErrorCode ierr; 59 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 60 61 PetscFunctionBegin; 62 ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 63 ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 64 ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 65 ierr = PetscOptionsTail();CHKERRQ(ierr); 66 PetscFunctionReturn(0); 67 } 68 69 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 70 { 71 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 72 PetscErrorCode ierr; 73 74 PetscFunctionBegin; 75 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 76 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 77 PetscFunctionReturn(0); 78 } 79 80 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 81 { 82 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 83 PetscErrorCode ierr; 84 85 PetscFunctionBegin; 86 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 87 while (hdf5->groups) { 88 GroupList *tmp = hdf5->groups->next; 89 90 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 91 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 92 hdf5->groups = tmp; 93 } 94 ierr = PetscFree(hdf5);CHKERRQ(ierr); 95 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 96 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 97 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 98 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 99 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 100 PetscFunctionReturn(0); 101 } 102 103 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 104 { 105 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 106 107 PetscFunctionBegin; 108 hdf5->btype = type; 109 PetscFunctionReturn(0); 110 } 111 112 static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 113 { 114 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 115 116 PetscFunctionBegin; 117 *type = hdf5->btype; 118 PetscFunctionReturn(0); 119 } 120 121 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 122 { 123 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 124 125 PetscFunctionBegin; 126 hdf5->basedimension2 = flg; 127 PetscFunctionReturn(0); 128 } 129 130 /*@ 131 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 132 dimension of 2. 133 134 Logically Collective on PetscViewer 135 136 Input Parameters: 137 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 138 - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 139 140 Options Database: 141 . -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 142 143 144 Notes: 145 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 146 of one when the dimension is lower. Others think the option is crazy. 147 148 Level: intermediate 149 150 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 151 152 @*/ 153 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 154 { 155 PetscErrorCode ierr; 156 157 PetscFunctionBegin; 158 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 159 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 160 PetscFunctionReturn(0); 161 } 162 163 /*@ 164 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 165 dimension of 2. 166 167 Logically Collective on PetscViewer 168 169 Input Parameter: 170 . viewer - the PetscViewer, must be of type HDF5 171 172 Output Parameter: 173 . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 174 175 Notes: 176 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 177 of one when the dimension is lower. Others think the option is crazy. 178 179 Level: intermediate 180 181 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 182 183 @*/ 184 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 185 { 186 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 187 188 PetscFunctionBegin; 189 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 190 *flg = hdf5->basedimension2; 191 PetscFunctionReturn(0); 192 } 193 194 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 195 { 196 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 197 198 PetscFunctionBegin; 199 hdf5->spoutput = flg; 200 PetscFunctionReturn(0); 201 } 202 203 /*@ 204 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 205 compiled with double precision PetscReal. 206 207 Logically Collective on PetscViewer 208 209 Input Parameters: 210 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 211 - flg - if PETSC_TRUE the data will be written to disk with single precision 212 213 Options Database: 214 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 215 216 217 Notes: 218 Setting this option does not make any difference if PETSc is compiled with single precision 219 in the first place. It does not affect reading datasets (HDF5 handle this internally). 220 221 Level: intermediate 222 223 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 224 PetscReal 225 226 @*/ 227 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 228 { 229 PetscErrorCode ierr; 230 231 PetscFunctionBegin; 232 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 233 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 234 PetscFunctionReturn(0); 235 } 236 237 /*@ 238 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 239 compiled with double precision PetscReal. 240 241 Logically Collective on PetscViewer 242 243 Input Parameter: 244 . viewer - the PetscViewer, must be of type HDF5 245 246 Output Parameter: 247 . flg - if PETSC_TRUE the data will be written to disk with single precision 248 249 Notes: 250 Setting this option does not make any difference if PETSc is compiled with single precision 251 in the first place. It does not affect reading datasets (HDF5 handle this internally). 252 253 Level: intermediate 254 255 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 256 PetscReal 257 258 @*/ 259 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 260 { 261 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 262 263 PetscFunctionBegin; 264 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 265 *flg = hdf5->spoutput; 266 PetscFunctionReturn(0); 267 } 268 269 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 270 { 271 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 272 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 273 MPI_Info info = MPI_INFO_NULL; 274 #endif 275 hid_t plist_id; 276 PetscErrorCode ierr; 277 278 PetscFunctionBegin; 279 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 280 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 281 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 282 /* Set up file access property list with parallel I/O access */ 283 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 284 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 285 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 286 #endif 287 /* Create or open the file collectively */ 288 switch (hdf5->btype) { 289 case FILE_MODE_READ: 290 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 291 break; 292 case FILE_MODE_APPEND: 293 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 294 break; 295 case FILE_MODE_WRITE: 296 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 297 break; 298 default: 299 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 300 } 301 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 302 PetscStackCallHDF5(H5Pclose,(plist_id)); 303 PetscFunctionReturn(0); 304 } 305 306 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 307 { 308 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 309 310 PetscFunctionBegin; 311 *name = vhdf5->filename; 312 PetscFunctionReturn(0); 313 } 314 315 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 316 { 317 /* 318 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 319 PetscErrorCode ierr; 320 */ 321 322 PetscFunctionBegin; 323 PetscFunctionReturn(0); 324 } 325 326 /*MC 327 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 328 329 330 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 331 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 332 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 333 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 334 335 Level: beginner 336 M*/ 337 338 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 339 { 340 PetscViewer_HDF5 *hdf5; 341 PetscErrorCode ierr; 342 343 PetscFunctionBegin; 344 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 345 346 v->data = (void*) hdf5; 347 v->ops->destroy = PetscViewerDestroy_HDF5; 348 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 349 v->ops->setup = PetscViewerSetUp_HDF5; 350 v->ops->flush = 0; 351 hdf5->btype = (PetscFileMode) -1; 352 hdf5->filename = 0; 353 hdf5->timestep = -1; 354 hdf5->groups = NULL; 355 356 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 357 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 358 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 359 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 360 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 361 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 /*@C 366 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 367 368 Collective on MPI_Comm 369 370 Input Parameters: 371 + comm - MPI communicator 372 . name - name of file 373 - type - type of file 374 $ FILE_MODE_WRITE - create new file for binary output 375 $ FILE_MODE_READ - open existing file for binary input 376 $ FILE_MODE_APPEND - open existing file for binary output 377 378 Output Parameter: 379 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 380 381 Options Database: 382 . -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 383 . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 384 385 Level: beginner 386 387 Note: 388 This PetscViewer should be destroyed with PetscViewerDestroy(). 389 390 Concepts: HDF5 files 391 Concepts: PetscViewerHDF5^creating 392 393 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 394 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 395 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 396 @*/ 397 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 398 { 399 PetscErrorCode ierr; 400 401 PetscFunctionBegin; 402 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 403 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 404 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 405 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 406 PetscFunctionReturn(0); 407 } 408 409 /*@C 410 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 411 412 Not collective 413 414 Input Parameter: 415 . viewer - the PetscViewer 416 417 Output Parameter: 418 . file_id - The file id 419 420 Level: intermediate 421 422 .seealso: PetscViewerHDF5Open() 423 @*/ 424 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 425 { 426 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 427 428 PetscFunctionBegin; 429 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 430 if (file_id) *file_id = hdf5->file_id; 431 PetscFunctionReturn(0); 432 } 433 434 /*@C 435 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 436 437 Not collective 438 439 Input Parameters: 440 + viewer - the PetscViewer 441 - name - The group name 442 443 Level: intermediate 444 445 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 446 @*/ 447 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 448 { 449 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 450 GroupList *groupNode; 451 PetscErrorCode ierr; 452 453 PetscFunctionBegin; 454 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 455 PetscValidCharPointer(name,2); 456 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 457 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 458 459 groupNode->next = hdf5->groups; 460 hdf5->groups = groupNode; 461 PetscFunctionReturn(0); 462 } 463 464 /*@ 465 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 466 467 Not collective 468 469 Input Parameter: 470 . viewer - the PetscViewer 471 472 Level: intermediate 473 474 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 475 @*/ 476 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 477 { 478 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 479 GroupList *groupNode; 480 PetscErrorCode ierr; 481 482 PetscFunctionBegin; 483 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 484 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 485 groupNode = hdf5->groups; 486 hdf5->groups = hdf5->groups->next; 487 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 488 ierr = PetscFree(groupNode);CHKERRQ(ierr); 489 PetscFunctionReturn(0); 490 } 491 492 /*@C 493 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 494 If none has been assigned, returns NULL. 495 496 Not collective 497 498 Input Parameter: 499 . viewer - the PetscViewer 500 501 Output Parameter: 502 . name - The group name 503 504 Level: intermediate 505 506 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 507 @*/ 508 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 509 { 510 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 511 512 PetscFunctionBegin; 513 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 514 PetscValidPointer(name,2); 515 if (hdf5->groups) *name = hdf5->groups->name; 516 else *name = NULL; 517 PetscFunctionReturn(0); 518 } 519 520 /*@ 521 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 522 and return this group's ID and file ID. 523 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 524 525 Not collective 526 527 Input Parameter: 528 . viewer - the PetscViewer 529 530 Output Parameter: 531 + fileId - The HDF5 file ID 532 - groupId - The HDF5 group ID 533 534 Notes: 535 If the viewer is writable, the group is created if it doesn't exist yet. 536 537 Level: intermediate 538 539 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 540 @*/ 541 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 542 { 543 hid_t file_id; 544 H5O_type_t type; 545 const char *groupName = NULL; 546 PetscBool create; 547 PetscErrorCode ierr; 548 549 PetscFunctionBegin; 550 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 551 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 552 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 553 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 554 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 555 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 556 *fileId = file_id; 557 PetscFunctionReturn(0); 558 } 559 560 /*@ 561 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 562 563 Not collective 564 565 Input Parameter: 566 . viewer - the PetscViewer 567 568 Level: intermediate 569 570 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 571 @*/ 572 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 573 { 574 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 575 576 PetscFunctionBegin; 577 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 578 ++hdf5->timestep; 579 PetscFunctionReturn(0); 580 } 581 582 /*@ 583 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 584 of -1 disables blocking with timesteps. 585 586 Not collective 587 588 Input Parameters: 589 + viewer - the PetscViewer 590 - timestep - The timestep number 591 592 Level: intermediate 593 594 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 595 @*/ 596 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 597 { 598 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 599 600 PetscFunctionBegin; 601 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 602 hdf5->timestep = timestep; 603 PetscFunctionReturn(0); 604 } 605 606 /*@ 607 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 608 609 Not collective 610 611 Input Parameter: 612 . viewer - the PetscViewer 613 614 Output Parameter: 615 . timestep - The timestep number 616 617 Level: intermediate 618 619 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 620 @*/ 621 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 622 { 623 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 624 625 PetscFunctionBegin; 626 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 627 PetscValidPointer(timestep,2); 628 *timestep = hdf5->timestep; 629 PetscFunctionReturn(0); 630 } 631 632 /*@C 633 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 634 635 Not collective 636 637 Input Parameter: 638 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 639 640 Output Parameter: 641 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 642 643 Level: advanced 644 645 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 646 @*/ 647 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 648 { 649 PetscFunctionBegin; 650 if (ptype == PETSC_INT) 651 #if defined(PETSC_USE_64BIT_INDICES) 652 *htype = H5T_NATIVE_LLONG; 653 #else 654 *htype = H5T_NATIVE_INT; 655 #endif 656 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 657 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 658 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 659 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 660 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 661 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 662 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 663 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 664 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 665 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 666 PetscFunctionReturn(0); 667 } 668 669 /*@C 670 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 671 672 Not collective 673 674 Input Parameter: 675 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 676 677 Output Parameter: 678 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 679 680 Level: advanced 681 682 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 683 @*/ 684 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 685 { 686 PetscFunctionBegin; 687 #if defined(PETSC_USE_64BIT_INDICES) 688 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 689 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 690 #else 691 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 692 #endif 693 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 694 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 695 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 696 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 697 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 698 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 699 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 700 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 701 PetscFunctionReturn(0); 702 } 703 704 /*@C 705 PetscViewerHDF5WriteAttribute - Write an attribute 706 707 Input Parameters: 708 + viewer - The HDF5 viewer 709 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 710 . name - The attribute name 711 . datatype - The attribute type 712 - value - The attribute value 713 714 Level: advanced 715 716 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 717 @*/ 718 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 719 { 720 char *parent; 721 hid_t h5, dataspace, obj, attribute, dtype; 722 PetscBool has; 723 PetscErrorCode ierr; 724 725 PetscFunctionBegin; 726 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 727 if (dataset) PetscValidCharPointer(dataset, 2); 728 PetscValidCharPointer(name, 3); 729 PetscValidPointer(value, 5); 730 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 731 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 732 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 733 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 734 if (datatype == PETSC_STRING) { 735 size_t len; 736 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 737 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 738 } 739 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 740 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 741 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 742 if (has) { 743 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 744 } else { 745 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 746 } 747 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 748 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 749 PetscStackCallHDF5(H5Aclose,(attribute)); 750 PetscStackCallHDF5(H5Oclose,(obj)); 751 PetscStackCallHDF5(H5Sclose,(dataspace)); 752 ierr = PetscFree(parent);CHKERRQ(ierr); 753 PetscFunctionReturn(0); 754 } 755 756 /*@C 757 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 758 759 Input Parameters: 760 + viewer - The HDF5 viewer 761 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 762 . name - The attribute name 763 . datatype - The attribute type 764 - value - The attribute value 765 766 Notes: 767 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 768 You might want to check first if it does using PetscViewerHDF5HasObject(). 769 770 Level: advanced 771 772 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 773 @*/ 774 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 775 { 776 PetscErrorCode ierr; 777 778 PetscFunctionBegin; 779 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 780 PetscValidHeader(obj,2); 781 PetscValidCharPointer(name,3); 782 PetscValidPointer(value,5); 783 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 784 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 785 PetscFunctionReturn(0); 786 } 787 788 /*@C 789 PetscViewerHDF5ReadAttribute - Read an attribute 790 791 Input Parameters: 792 + viewer - The HDF5 viewer 793 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 794 . name - The attribute name 795 - datatype - The attribute type 796 797 Output Parameter: 798 . value - The attribute value 799 800 Level: advanced 801 802 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 803 @*/ 804 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 805 { 806 char *parent; 807 hid_t h5, obj, attribute, atype, dtype; 808 PetscBool has; 809 PetscErrorCode ierr; 810 811 PetscFunctionBegin; 812 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 813 if (dataset) PetscValidCharPointer(dataset, 2); 814 PetscValidCharPointer(name, 3); 815 PetscValidPointer(value, 5); 816 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 817 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 818 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 819 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 820 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 821 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 822 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 823 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 824 if (datatype == PETSC_STRING) { 825 size_t len; 826 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 827 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 828 PetscStackCallHDF5(H5Tclose,(atype)); 829 ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 830 } 831 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 832 PetscStackCallHDF5(H5Aclose,(attribute)); 833 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 834 PetscStackCallHDF5(H5Oclose,(obj)); 835 ierr = PetscFree(parent);CHKERRQ(ierr); 836 PetscFunctionReturn(0); 837 } 838 839 /*@C 840 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 841 842 Input Parameters: 843 + viewer - The HDF5 viewer 844 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 845 . name - The attribute name 846 - datatype - The attribute type 847 848 Output Parameter: 849 . value - The attribute value 850 851 Notes: 852 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 853 You might want to check first if it does using PetscViewerHDF5HasObject(). 854 855 Level: advanced 856 857 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 858 @*/ 859 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 860 { 861 PetscErrorCode ierr; 862 863 PetscFunctionBegin; 864 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 865 PetscValidHeader(obj,2); 866 PetscValidCharPointer(name,3); 867 PetscValidPointer(value, 5); 868 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 869 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 870 PetscFunctionReturn(0); 871 } 872 873 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 874 { 875 htri_t exists; 876 hid_t group; 877 878 PetscFunctionBegin; 879 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 880 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 881 if (!exists && createGroup) { 882 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 883 PetscStackCallHDF5(H5Gclose,(group)); 884 exists = PETSC_TRUE; 885 } 886 *exists_ = (PetscBool) exists; 887 PetscFunctionReturn(0); 888 } 889 890 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 891 { 892 const char rootGroupName[] = "/"; 893 hid_t h5; 894 PetscBool exists=PETSC_FALSE; 895 PetscInt i; 896 int n; 897 char **hierarchy; 898 char buf[PETSC_MAX_PATH_LEN]=""; 899 PetscErrorCode ierr; 900 901 PetscFunctionBegin; 902 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 903 if (name) PetscValidCharPointer(name, 2); 904 else name = rootGroupName; 905 if (has) { 906 PetscValidIntPointer(has, 3); 907 *has = PETSC_FALSE; 908 } 909 if (otype) { 910 PetscValidIntPointer(otype, 4); 911 *otype = H5O_TYPE_UNKNOWN; 912 } 913 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 914 915 /* 916 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 917 Hence, each of them needs to be tested separately: 918 1) whether it's a valid link 919 2) whether this link resolves to an object 920 See H5Oexists_by_name() documentation. 921 */ 922 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 923 if (!n) { 924 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 925 if (has) *has = PETSC_TRUE; 926 if (otype) *otype = H5O_TYPE_GROUP; 927 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 928 PetscFunctionReturn(0); 929 } 930 for (i=0; i<n; i++) { 931 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 932 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 933 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 934 if (!exists) break; 935 } 936 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 937 938 /* If the object exists, get its type */ 939 if (exists && otype) { 940 H5O_info_t info; 941 942 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 943 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 944 *otype = info.type; 945 } 946 if (has) *has = exists; 947 PetscFunctionReturn(0); 948 } 949 950 /*@ 951 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 952 953 Input Parameters: 954 . viewer - The HDF5 viewer 955 956 Output Parameter: 957 . has - Flag for group existence 958 959 Notes: 960 If the path exists but is not a group, this returns PETSC_FALSE as well. 961 962 Level: advanced 963 964 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 965 @*/ 966 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 967 { 968 H5O_type_t type; 969 const char *name; 970 PetscErrorCode ierr; 971 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 974 PetscValidIntPointer(has,2); 975 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 976 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 977 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 978 PetscFunctionReturn(0); 979 } 980 981 /*@ 982 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 983 984 Input Parameters: 985 + viewer - The HDF5 viewer 986 - obj - The named object 987 988 Output Parameter: 989 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 990 991 Notes: 992 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 993 994 Level: advanced 995 996 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 997 @*/ 998 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 999 { 1000 H5O_type_t type; 1001 char *path; 1002 PetscErrorCode ierr; 1003 1004 PetscFunctionBegin; 1005 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1006 PetscValidHeader(obj,2); 1007 PetscValidIntPointer(has,3); 1008 *has = PETSC_FALSE; 1009 if (!obj->name) PetscFunctionReturn(0); 1010 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1011 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1012 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1013 ierr = PetscFree(path);CHKERRQ(ierr); 1014 PetscFunctionReturn(0); 1015 } 1016 1017 /*@C 1018 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1019 1020 Input Parameters: 1021 + viewer - The HDF5 viewer 1022 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1023 - name - The attribute name 1024 1025 Output Parameter: 1026 . has - Flag for attribute existence 1027 1028 Level: advanced 1029 1030 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1031 @*/ 1032 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1033 { 1034 char *parent; 1035 PetscErrorCode ierr; 1036 1037 PetscFunctionBegin; 1038 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1039 if (dataset) PetscValidCharPointer(dataset,2); 1040 PetscValidCharPointer(name,3); 1041 PetscValidIntPointer(has,4); 1042 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1043 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1044 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1045 ierr = PetscFree(parent);CHKERRQ(ierr); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 /*@C 1050 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1051 1052 Input Parameters: 1053 + viewer - The HDF5 viewer 1054 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1055 - name - The attribute name 1056 1057 Output Parameter: 1058 . has - Flag for attribute existence 1059 1060 Notes: 1061 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1062 You might want to check first if it does using PetscViewerHDF5HasObject(). 1063 1064 Level: advanced 1065 1066 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1067 @*/ 1068 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1069 { 1070 PetscErrorCode ierr; 1071 1072 PetscFunctionBegin; 1073 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1074 PetscValidHeader(obj,2); 1075 PetscValidCharPointer(name,3); 1076 PetscValidIntPointer(has,4); 1077 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1078 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1079 PetscFunctionReturn(0); 1080 } 1081 1082 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1083 { 1084 hid_t h5; 1085 htri_t hhas; 1086 PetscErrorCode ierr; 1087 1088 PetscFunctionBegin; 1089 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1090 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1091 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1092 PetscFunctionReturn(0); 1093 } 1094 1095 /* 1096 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1097 is attached to a communicator, in this case the attribute is a PetscViewer. 1098 */ 1099 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1100 1101 /*@C 1102 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1103 1104 Collective on MPI_Comm 1105 1106 Input Parameter: 1107 . comm - the MPI communicator to share the HDF5 PetscViewer 1108 1109 Level: intermediate 1110 1111 Options Database Keys: 1112 . -viewer_hdf5_filename <name> 1113 1114 Environmental variables: 1115 . PETSC_VIEWER_HDF5_FILENAME 1116 1117 Notes: 1118 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1119 an error code. The HDF5 PetscViewer is usually used in the form 1120 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1121 1122 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1123 @*/ 1124 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1125 { 1126 PetscErrorCode ierr; 1127 PetscBool flg; 1128 PetscViewer viewer; 1129 char fname[PETSC_MAX_PATH_LEN]; 1130 MPI_Comm ncomm; 1131 1132 PetscFunctionBegin; 1133 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1134 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1135 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1136 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1137 } 1138 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1139 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1140 if (!flg) { /* PetscViewer not yet created */ 1141 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1142 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1143 if (!flg) { 1144 ierr = PetscStrcpy(fname,"output.h5"); 1145 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1146 } 1147 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1148 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1149 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1150 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1151 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1152 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1153 } 1154 ierr = PetscCommDestroy(&ncomm); 1155 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1156 PetscFunctionReturn(viewer); 1157 } 1158