1 #include <petsc-private/viewerimpl.h> /*I "petscsys.h" I*/ 2 #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 3 4 typedef struct GroupList { 5 const char *name; 6 struct GroupList *next; 7 } GroupList; 8 9 typedef struct { 10 char *filename; 11 PetscFileMode btype; 12 hid_t file_id; 13 PetscInt timestep; 14 GroupList *groups; 15 } PetscViewer_HDF5; 16 17 #undef __FUNCT__ 18 #define __FUNCT__ "PetscViewerFileClose_HDF5" 19 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 20 { 21 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 22 PetscErrorCode ierr; 23 24 PetscFunctionBegin; 25 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 26 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 27 PetscFunctionReturn(0); 28 } 29 30 #undef __FUNCT__ 31 #define __FUNCT__ "PetscViewerDestroy_HDF5" 32 PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 33 { 34 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 39 while (hdf5->groups) { 40 GroupList *tmp = hdf5->groups->next; 41 42 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 43 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 44 hdf5->groups = tmp; 45 } 46 ierr = PetscFree(hdf5);CHKERRQ(ierr); 47 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 48 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "PetscViewerFileSetMode_HDF5" 54 PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 55 { 56 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 57 58 PetscFunctionBegin; 59 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 60 hdf5->btype = type; 61 PetscFunctionReturn(0); 62 } 63 64 #undef __FUNCT__ 65 #define __FUNCT__ "PetscViewerFileSetName_HDF5" 66 PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 67 { 68 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 69 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 70 MPI_Info info = MPI_INFO_NULL; 71 #endif 72 hid_t plist_id; 73 PetscErrorCode ierr; 74 75 PetscFunctionBegin; 76 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 77 /* Set up file access property list with parallel I/O access */ 78 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 79 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 80 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 81 #endif 82 /* Create or open the file collectively */ 83 switch (hdf5->btype) { 84 case FILE_MODE_READ: 85 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 86 break; 87 case FILE_MODE_APPEND: 88 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 89 break; 90 case FILE_MODE_WRITE: 91 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 92 break; 93 default: 94 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 95 } 96 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 97 PetscStackCallHDF5(H5Pclose,(plist_id)); 98 PetscFunctionReturn(0); 99 } 100 101 #undef __FUNCT__ 102 #define __FUNCT__ "PetscViewerCreate_HDF5" 103 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 104 { 105 PetscViewer_HDF5 *hdf5; 106 PetscErrorCode ierr; 107 108 PetscFunctionBegin; 109 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 110 111 v->data = (void*) hdf5; 112 v->ops->destroy = PetscViewerDestroy_HDF5; 113 v->ops->flush = 0; 114 hdf5->btype = (PetscFileMode) -1; 115 hdf5->filename = 0; 116 hdf5->timestep = -1; 117 hdf5->groups = NULL; 118 119 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 120 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 124 #undef __FUNCT__ 125 #define __FUNCT__ "PetscViewerHDF5Open" 126 /*@C 127 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 128 129 Collective on MPI_Comm 130 131 Input Parameters: 132 + comm - MPI communicator 133 . name - name of file 134 - type - type of file 135 $ FILE_MODE_WRITE - create new file for binary output 136 $ FILE_MODE_READ - open existing file for binary input 137 $ FILE_MODE_APPEND - open existing file for binary output 138 139 Output Parameter: 140 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 141 142 Level: beginner 143 144 Note: 145 This PetscViewer should be destroyed with PetscViewerDestroy(). 146 147 Concepts: HDF5 files 148 Concepts: PetscViewerHDF5^creating 149 150 .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(), 151 VecView(), MatView(), VecLoad(), MatLoad(), 152 PetscFileMode, PetscViewer 153 @*/ 154 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 155 { 156 PetscErrorCode ierr; 157 158 PetscFunctionBegin; 159 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 160 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 161 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 162 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 163 PetscFunctionReturn(0); 164 } 165 166 #undef __FUNCT__ 167 #define __FUNCT__ "PetscViewerHDF5GetFileId" 168 /*@C 169 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 170 171 Not collective 172 173 Input Parameter: 174 . viewer - the PetscViewer 175 176 Output Parameter: 177 . file_id - The file id 178 179 Level: intermediate 180 181 .seealso: PetscViewerHDF5Open() 182 @*/ 183 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 184 { 185 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 186 187 PetscFunctionBegin; 188 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 189 if (file_id) *file_id = hdf5->file_id; 190 PetscFunctionReturn(0); 191 } 192 193 #undef __FUNCT__ 194 #define __FUNCT__ "PetscViewerHDF5PushGroup" 195 /*@C 196 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 197 198 Not collective 199 200 Input Parameters: 201 + viewer - the PetscViewer 202 - name - The group name 203 204 Level: intermediate 205 206 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 207 @*/ 208 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 209 { 210 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 211 GroupList *groupNode; 212 PetscErrorCode ierr; 213 214 PetscFunctionBegin; 215 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 216 PetscValidCharPointer(name,2); 217 ierr = PetscMalloc(sizeof(GroupList), &groupNode);CHKERRQ(ierr); 218 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 219 220 groupNode->next = hdf5->groups; 221 hdf5->groups = groupNode; 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "PetscViewerHDF5PopGroup" 227 /*@ 228 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 229 230 Not collective 231 232 Input Parameter: 233 . viewer - the PetscViewer 234 235 Level: intermediate 236 237 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup() 238 @*/ 239 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 240 { 241 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 242 GroupList *groupNode; 243 PetscErrorCode ierr; 244 245 PetscFunctionBegin; 246 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 247 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 248 groupNode = hdf5->groups; 249 hdf5->groups = hdf5->groups->next; 250 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 251 ierr = PetscFree(groupNode);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "PetscViewerHDF5GetGroup" 257 /*@C 258 PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL. 259 260 Not collective 261 262 Input Parameter: 263 . viewer - the PetscViewer 264 265 Output Parameter: 266 . name - The group name 267 268 Level: intermediate 269 270 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup() 271 @*/ 272 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 273 { 274 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 275 276 PetscFunctionBegin; 277 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 278 PetscValidPointer(name,2); 279 if (hdf5->groups) *name = hdf5->groups->name; 280 else *name = NULL; 281 PetscFunctionReturn(0); 282 } 283 284 #undef __FUNCT__ 285 #define __FUNCT__ "PetscViewerHDF5IncrementTimestep" 286 /*@ 287 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 288 289 Not collective 290 291 Input Parameter: 292 . viewer - the PetscViewer 293 294 Level: intermediate 295 296 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 297 @*/ 298 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 299 { 300 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 301 302 PetscFunctionBegin; 303 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 304 ++hdf5->timestep; 305 PetscFunctionReturn(0); 306 } 307 308 #undef __FUNCT__ 309 #define __FUNCT__ "PetscViewerHDF5SetTimestep" 310 /*@ 311 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 312 of -1 disables blocking with timesteps. 313 314 Not collective 315 316 Input Parameters: 317 + viewer - the PetscViewer 318 - timestep - The timestep number 319 320 Level: intermediate 321 322 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 323 @*/ 324 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 325 { 326 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 327 328 PetscFunctionBegin; 329 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 330 hdf5->timestep = timestep; 331 PetscFunctionReturn(0); 332 } 333 334 #undef __FUNCT__ 335 #define __FUNCT__ "PetscViewerHDF5GetTimestep" 336 /*@ 337 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 338 339 Not collective 340 341 Input Parameter: 342 . viewer - the PetscViewer 343 344 Output Parameter: 345 . timestep - The timestep number 346 347 Level: intermediate 348 349 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 350 @*/ 351 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 352 { 353 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 354 355 PetscFunctionBegin; 356 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 357 PetscValidPointer(timestep,2); 358 *timestep = hdf5->timestep; 359 PetscFunctionReturn(0); 360 } 361 362 #undef __FUNCT__ 363 #define __FUNCT__ "PetscDataTypeToHDF5DataType" 364 /*@C 365 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 366 367 Not collective 368 369 Input Parameter: 370 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 371 372 Output Parameter: 373 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 374 375 Level: advanced 376 377 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 378 @*/ 379 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 380 { 381 PetscFunctionBegin; 382 if (ptype == PETSC_INT) 383 #if defined(PETSC_USE_64BIT_INDICES) 384 *htype = H5T_NATIVE_LLONG; 385 #else 386 *htype = H5T_NATIVE_INT; 387 #endif 388 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 389 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 390 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 391 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 392 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 393 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 394 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 395 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 396 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 397 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 398 PetscFunctionReturn(0); 399 } 400 401 #undef __FUNCT__ 402 #define __FUNCT__ "PetscHDF5DataTypeToPetscDataType" 403 /*@C 404 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 405 406 Not collective 407 408 Input Parameter: 409 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 410 411 Output Parameter: 412 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 413 414 Level: advanced 415 416 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 417 @*/ 418 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 419 { 420 PetscFunctionBegin; 421 #if defined(PETSC_USE_64BIT_INDICES) 422 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 423 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 424 #else 425 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 426 #endif 427 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 428 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 429 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 430 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 431 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 432 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 433 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 434 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 435 PetscFunctionReturn(0); 436 } 437 438 #undef __FUNCT__ 439 #define __FUNCT__ "PetscViewerHDF5WriteAttribute" 440 /*@ 441 PetscViewerHDF5WriteAttribute - Write a scalar attribute 442 443 Input Parameters: 444 + viewer - The HDF5 viewer 445 . parent - The parent name 446 . name - The attribute name 447 . datatype - The attribute type 448 - value - The attribute value 449 450 Level: advanced 451 452 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 453 @*/ 454 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 455 { 456 hid_t h5, dataspace, dataset, attribute, dtype; 457 PetscErrorCode ierr; 458 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 461 PetscValidPointer(parent, 2); 462 PetscValidPointer(name, 3); 463 PetscValidPointer(value, 4); 464 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 465 if (datatype == PETSC_STRING) { 466 size_t len; 467 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 468 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 469 } 470 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 471 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 472 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 473 PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 474 PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 475 #else 476 PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 477 PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT)); 478 #endif 479 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 480 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 481 PetscStackCallHDF5(H5Aclose,(attribute)); 482 PetscStackCallHDF5(H5Dclose,(dataset)); 483 PetscStackCallHDF5(H5Sclose,(dataspace)); 484 PetscFunctionReturn(0); 485 } 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "PetscViewerHDF5ReadAttribute" 489 /*@ 490 PetscViewerHDF5ReadAttribute - Read a scalar attribute 491 492 Input Parameters: 493 + viewer - The HDF5 viewer 494 . parent - The parent name 495 . name - The attribute name 496 - datatype - The attribute type 497 498 Output Parameter: 499 . value - The attribute value 500 501 Level: advanced 502 503 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 504 @*/ 505 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 506 { 507 hid_t h5, dataspace, dataset, attribute, dtype; 508 PetscErrorCode ierr; 509 510 PetscFunctionBegin; 511 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 512 PetscValidPointer(parent, 2); 513 PetscValidPointer(name, 3); 514 PetscValidPointer(value, 4); 515 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 516 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 517 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 518 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 519 PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 520 #else 521 PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 522 #endif 523 PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name)); 524 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 525 PetscStackCallHDF5(H5Aclose,(attribute)); 526 PetscStackCallHDF5(H5Dclose,(dataset)); 527 PetscStackCallHDF5(H5Sclose,(dataspace)); 528 PetscFunctionReturn(0); 529 } 530 531 #undef __FUNCT__ 532 #define __FUNCT__ "PetscViewerHDF5HasObject" 533 static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has) 534 { 535 hid_t h5; 536 PetscErrorCode ierr; 537 538 PetscFunctionBegin; 539 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 540 PetscValidPointer(name, 2); 541 PetscValidPointer(has, 3); 542 *has = PETSC_FALSE; 543 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 544 if (H5Lexists(h5, name, H5P_DEFAULT)) { 545 H5O_info_t info; 546 hid_t obj; 547 548 PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT)); 549 PetscStackCallHDF5(H5Oget_info,(obj, &info)); 550 if (otype == info.type) *has = PETSC_TRUE; 551 PetscStackCallHDF5(H5Oclose,(obj)); 552 } 553 PetscFunctionReturn(0); 554 } 555 556 #undef __FUNCT__ 557 #define __FUNCT__ "PetscViewerHDF5HasAttribute" 558 /*@ 559 PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists 560 561 Input Parameters: 562 + viewer - The HDF5 viewer 563 . parent - The parent name 564 - name - The attribute name 565 566 Output Parameter: 567 . has - Flag for attribute existence 568 569 Level: advanced 570 571 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute() 572 @*/ 573 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 574 { 575 hid_t h5, dataset; 576 htri_t hhas; 577 PetscBool exists; 578 PetscErrorCode ierr; 579 580 PetscFunctionBegin; 581 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 582 PetscValidPointer(parent, 2); 583 PetscValidPointer(name, 3); 584 PetscValidPointer(has, 4); 585 *has = PETSC_FALSE; 586 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 587 ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr); 588 if (exists) { 589 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 590 PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT)); 591 #else 592 PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent)); 593 #endif 594 if (dataset < 0) PetscFunctionReturn(0); 595 PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name)); 596 if (hhas < 0) { 597 PetscStackCallHDF5(H5Dclose,(dataset)); 598 PetscFunctionReturn(0); 599 } 600 PetscStackCallHDF5(H5Dclose,(dataset)); 601 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 602 } 603 PetscFunctionReturn(0); 604 } 605 606 /* 607 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 608 is attached to a communicator, in this case the attribute is a PetscViewer. 609 */ 610 static int Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 611 612 #undef __FUNCT__ 613 #define __FUNCT__ "PETSC_VIEWER_HDF5_" 614 /*@C 615 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 616 617 Collective on MPI_Comm 618 619 Input Parameter: 620 . comm - the MPI communicator to share the HDF5 PetscViewer 621 622 Level: intermediate 623 624 Options Database Keys: 625 . -viewer_hdf5_filename <name> 626 627 Environmental variables: 628 . PETSC_VIEWER_HDF5_FILENAME 629 630 Notes: 631 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 632 an error code. The HDF5 PetscViewer is usually used in the form 633 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 634 635 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 636 @*/ 637 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 638 { 639 PetscErrorCode ierr; 640 PetscBool flg; 641 PetscViewer viewer; 642 char fname[PETSC_MAX_PATH_LEN]; 643 MPI_Comm ncomm; 644 645 PetscFunctionBegin; 646 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 647 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 648 ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 649 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 650 } 651 ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 652 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 653 if (!flg) { /* PetscViewer not yet created */ 654 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 655 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 656 if (!flg) { 657 ierr = PetscStrcpy(fname,"output.h5"); 658 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 659 } 660 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 661 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 662 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 663 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 664 ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 665 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 666 } 667 ierr = PetscCommDestroy(&ncomm); 668 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 669 PetscFunctionReturn(viewer); 670 } 671