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 PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 16 PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 17 } PetscViewer_HDF5; 18 19 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 20 { 21 PetscErrorCode ierr; 22 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 23 24 PetscFunctionBegin; 25 ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 26 ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 27 ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 28 ierr = PetscOptionsTail();CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 33 { 34 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 35 PetscErrorCode ierr; 36 37 PetscFunctionBegin; 38 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 39 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 40 PetscFunctionReturn(0); 41 } 42 43 PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 44 { 45 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 50 while (hdf5->groups) { 51 GroupList *tmp = hdf5->groups->next; 52 53 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 54 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 55 hdf5->groups = tmp; 56 } 57 ierr = PetscFree(hdf5);CHKERRQ(ierr); 58 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 59 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 60 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 61 PetscFunctionReturn(0); 62 } 63 64 PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 65 { 66 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 67 68 PetscFunctionBegin; 69 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 70 hdf5->btype = type; 71 PetscFunctionReturn(0); 72 } 73 74 PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 75 { 76 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 77 78 PetscFunctionBegin; 79 hdf5->basedimension2 = flg; 80 PetscFunctionReturn(0); 81 } 82 83 /*@ 84 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 85 dimension of 2. 86 87 Logically Collective on PetscViewer 88 89 Input Parameters: 90 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 91 - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 92 93 Options Database: 94 . -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 95 96 97 Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 98 of one when the dimension is lower. Others think the option is crazy. 99 100 Level: intermediate 101 102 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 103 104 @*/ 105 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 106 { 107 PetscErrorCode ierr; 108 109 PetscFunctionBegin; 110 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 111 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 112 PetscFunctionReturn(0); 113 } 114 115 /*@ 116 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 117 dimension of 2. 118 119 Logically Collective on PetscViewer 120 121 Input Parameter: 122 . viewer - the PetscViewer, must be of type HDF5 123 124 Output Parameter: 125 . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 126 127 Notes: Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 128 of one when the dimension is lower. Others think the option is crazy. 129 130 Level: intermediate 131 132 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 133 134 @*/ 135 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 136 { 137 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 138 139 PetscFunctionBegin; 140 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 141 *flg = hdf5->basedimension2; 142 PetscFunctionReturn(0); 143 } 144 145 PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 146 { 147 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 148 149 PetscFunctionBegin; 150 hdf5->spoutput = flg; 151 PetscFunctionReturn(0); 152 } 153 154 /*@ 155 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 156 compiled with double precision PetscReal. 157 158 Logically Collective on PetscViewer 159 160 Input Parameters: 161 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 162 - flg - if PETSC_TRUE the data will be written to disk with single precision 163 164 Options Database: 165 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 166 167 168 Notes: Setting this option does not make any difference if PETSc is compiled with single precision 169 in the first place. It does not affect reading datasets (HDF5 handle this internally). 170 171 Level: intermediate 172 173 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 174 PetscReal 175 176 @*/ 177 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 178 { 179 PetscErrorCode ierr; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 183 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 /*@ 188 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 189 compiled with double precision PetscReal. 190 191 Logically Collective on PetscViewer 192 193 Input Parameter: 194 . viewer - the PetscViewer, must be of type HDF5 195 196 Output Parameter: 197 . flg - if PETSC_TRUE the data will be written to disk with single precision 198 199 Notes: Setting this option does not make any difference if PETSc is compiled with single precision 200 in the first place. It does not affect reading datasets (HDF5 handle this internally). 201 202 Level: intermediate 203 204 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 205 PetscReal 206 207 @*/ 208 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 209 { 210 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 211 212 PetscFunctionBegin; 213 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 214 *flg = hdf5->spoutput; 215 PetscFunctionReturn(0); 216 } 217 218 PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 219 { 220 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 221 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 222 MPI_Info info = MPI_INFO_NULL; 223 #endif 224 hid_t plist_id; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 229 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 230 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 231 /* Set up file access property list with parallel I/O access */ 232 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 233 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 234 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 235 #endif 236 /* Create or open the file collectively */ 237 switch (hdf5->btype) { 238 case FILE_MODE_READ: 239 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 240 break; 241 case FILE_MODE_APPEND: 242 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 243 break; 244 case FILE_MODE_WRITE: 245 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 246 break; 247 default: 248 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 249 } 250 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 251 PetscStackCallHDF5(H5Pclose,(plist_id)); 252 PetscFunctionReturn(0); 253 } 254 255 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 256 { 257 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 258 259 PetscFunctionBegin; 260 *name = vhdf5->filename; 261 PetscFunctionReturn(0); 262 } 263 264 /*MC 265 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 266 267 268 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 269 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 270 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 271 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 272 273 Level: beginner 274 M*/ 275 276 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 277 { 278 PetscViewer_HDF5 *hdf5; 279 PetscErrorCode ierr; 280 281 PetscFunctionBegin; 282 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 283 284 v->data = (void*) hdf5; 285 v->ops->destroy = PetscViewerDestroy_HDF5; 286 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 287 v->ops->flush = 0; 288 hdf5->btype = (PetscFileMode) -1; 289 hdf5->filename = 0; 290 hdf5->timestep = -1; 291 hdf5->groups = NULL; 292 293 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 294 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 295 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 296 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 297 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 298 PetscFunctionReturn(0); 299 } 300 301 /*@C 302 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 303 304 Collective on MPI_Comm 305 306 Input Parameters: 307 + comm - MPI communicator 308 . name - name of file 309 - type - type of file 310 $ FILE_MODE_WRITE - create new file for binary output 311 $ FILE_MODE_READ - open existing file for binary input 312 $ FILE_MODE_APPEND - open existing file for binary output 313 314 Output Parameter: 315 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 316 317 Options Database: 318 . -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 319 . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 320 321 Level: beginner 322 323 Note: 324 This PetscViewer should be destroyed with PetscViewerDestroy(). 325 326 Concepts: HDF5 files 327 Concepts: PetscViewerHDF5^creating 328 329 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 330 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 331 MatLoad(), PetscFileMode, PetscViewer 332 @*/ 333 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 334 { 335 PetscErrorCode ierr; 336 337 PetscFunctionBegin; 338 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 339 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 340 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 341 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 342 PetscFunctionReturn(0); 343 } 344 345 /*@C 346 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 347 348 Not collective 349 350 Input Parameter: 351 . viewer - the PetscViewer 352 353 Output Parameter: 354 . file_id - The file id 355 356 Level: intermediate 357 358 .seealso: PetscViewerHDF5Open() 359 @*/ 360 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 361 { 362 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 363 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 366 if (file_id) *file_id = hdf5->file_id; 367 PetscFunctionReturn(0); 368 } 369 370 /*@C 371 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 372 373 Not collective 374 375 Input Parameters: 376 + viewer - the PetscViewer 377 - name - The group name 378 379 Level: intermediate 380 381 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 382 @*/ 383 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 384 { 385 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 386 GroupList *groupNode; 387 PetscErrorCode ierr; 388 389 PetscFunctionBegin; 390 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 391 PetscValidCharPointer(name,2); 392 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 393 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 394 395 groupNode->next = hdf5->groups; 396 hdf5->groups = groupNode; 397 PetscFunctionReturn(0); 398 } 399 400 /*@ 401 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 402 403 Not collective 404 405 Input Parameter: 406 . viewer - the PetscViewer 407 408 Level: intermediate 409 410 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup() 411 @*/ 412 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 413 { 414 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 415 GroupList *groupNode; 416 PetscErrorCode ierr; 417 418 PetscFunctionBegin; 419 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 420 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 421 groupNode = hdf5->groups; 422 hdf5->groups = hdf5->groups->next; 423 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 424 ierr = PetscFree(groupNode);CHKERRQ(ierr); 425 PetscFunctionReturn(0); 426 } 427 428 /*@C 429 PetscViewerHDF5GetGroup - Get the current HDF5 group for output. If none has been assigned, returns NULL. 430 431 Not collective 432 433 Input Parameter: 434 . viewer - the PetscViewer 435 436 Output Parameter: 437 . name - The group name 438 439 Level: intermediate 440 441 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup() 442 @*/ 443 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 444 { 445 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 446 447 PetscFunctionBegin; 448 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 449 PetscValidPointer(name,2); 450 if (hdf5->groups) *name = hdf5->groups->name; 451 else *name = NULL; 452 PetscFunctionReturn(0); 453 } 454 455 /*@ 456 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 457 458 Not collective 459 460 Input Parameter: 461 . viewer - the PetscViewer 462 463 Level: intermediate 464 465 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 466 @*/ 467 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 468 { 469 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 470 471 PetscFunctionBegin; 472 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 473 ++hdf5->timestep; 474 PetscFunctionReturn(0); 475 } 476 477 /*@ 478 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 479 of -1 disables blocking with timesteps. 480 481 Not collective 482 483 Input Parameters: 484 + viewer - the PetscViewer 485 - timestep - The timestep number 486 487 Level: intermediate 488 489 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 490 @*/ 491 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 492 { 493 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 494 495 PetscFunctionBegin; 496 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 497 hdf5->timestep = timestep; 498 PetscFunctionReturn(0); 499 } 500 501 /*@ 502 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 503 504 Not collective 505 506 Input Parameter: 507 . viewer - the PetscViewer 508 509 Output Parameter: 510 . timestep - The timestep number 511 512 Level: intermediate 513 514 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 515 @*/ 516 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 517 { 518 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 519 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 522 PetscValidPointer(timestep,2); 523 *timestep = hdf5->timestep; 524 PetscFunctionReturn(0); 525 } 526 527 /*@C 528 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 529 530 Not collective 531 532 Input Parameter: 533 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 534 535 Output Parameter: 536 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 537 538 Level: advanced 539 540 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 541 @*/ 542 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 543 { 544 PetscFunctionBegin; 545 if (ptype == PETSC_INT) 546 #if defined(PETSC_USE_64BIT_INDICES) 547 *htype = H5T_NATIVE_LLONG; 548 #else 549 *htype = H5T_NATIVE_INT; 550 #endif 551 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 552 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 553 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 554 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 555 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 556 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 557 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 558 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 559 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 560 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 561 PetscFunctionReturn(0); 562 } 563 564 /*@C 565 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 566 567 Not collective 568 569 Input Parameter: 570 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 571 572 Output Parameter: 573 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 574 575 Level: advanced 576 577 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 578 @*/ 579 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 580 { 581 PetscFunctionBegin; 582 #if defined(PETSC_USE_64BIT_INDICES) 583 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 584 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 585 #else 586 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 587 #endif 588 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 589 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 590 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 591 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 592 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 593 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 594 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 595 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 596 PetscFunctionReturn(0); 597 } 598 599 /*@C 600 PetscViewerHDF5WriteAttribute - Write a scalar attribute 601 602 Input Parameters: 603 + viewer - The HDF5 viewer 604 . parent - The parent name 605 . name - The attribute name 606 . datatype - The attribute type 607 - value - The attribute value 608 609 Level: advanced 610 611 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 612 @*/ 613 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 614 { 615 hid_t h5, dataspace, dataset, attribute, dtype; 616 PetscErrorCode ierr; 617 618 PetscFunctionBegin; 619 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 620 PetscValidPointer(parent, 2); 621 PetscValidPointer(name, 3); 622 PetscValidPointer(value, 4); 623 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 624 if (datatype == PETSC_STRING) { 625 size_t len; 626 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 627 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 628 } 629 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 630 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 631 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 632 PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 633 PetscStackCallHDF5Return(attribute,H5Acreate2,(dataset, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 634 #else 635 PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 636 PetscStackCallHDF5Return(attribute,H5Acreate,(dataset, name, dtype, dataspace, H5P_DEFAULT)); 637 #endif 638 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 639 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 640 PetscStackCallHDF5(H5Aclose,(attribute)); 641 PetscStackCallHDF5(H5Dclose,(dataset)); 642 PetscStackCallHDF5(H5Sclose,(dataspace)); 643 PetscFunctionReturn(0); 644 } 645 646 /*@C 647 PetscViewerHDF5ReadAttribute - Read a scalar attribute 648 649 Input Parameters: 650 + viewer - The HDF5 viewer 651 . parent - The parent name 652 . name - The attribute name 653 - datatype - The attribute type 654 655 Output Parameter: 656 . value - The attribute value 657 658 Level: advanced 659 660 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 661 @*/ 662 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 663 { 664 hid_t h5, dataspace, dataset, attribute, dtype; 665 PetscErrorCode ierr; 666 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 669 PetscValidPointer(parent, 2); 670 PetscValidPointer(name, 3); 671 PetscValidPointer(value, 4); 672 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 673 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 674 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 675 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 676 PetscStackCallHDF5Return(dataset,H5Dopen2,(h5, parent, H5P_DEFAULT)); 677 #else 678 PetscStackCallHDF5Return(dataset,H5Dopen,(h5, parent)); 679 #endif 680 PetscStackCallHDF5Return(attribute,H5Aopen_name,(dataset, name)); 681 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 682 PetscStackCallHDF5(H5Aclose,(attribute)); 683 PetscStackCallHDF5(H5Dclose,(dataset)); 684 PetscStackCallHDF5(H5Sclose,(dataspace)); 685 PetscFunctionReturn(0); 686 } 687 688 static PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has) 689 { 690 hid_t h5; 691 PetscErrorCode ierr; 692 693 PetscFunctionBegin; 694 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 695 PetscValidPointer(name, 2); 696 PetscValidPointer(has, 3); 697 *has = PETSC_FALSE; 698 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 699 if (H5Lexists(h5, name, H5P_DEFAULT)) { 700 H5O_info_t info; 701 hid_t obj; 702 703 PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT)); 704 PetscStackCallHDF5(H5Oget_info,(obj, &info)); 705 if (otype == info.type) *has = PETSC_TRUE; 706 PetscStackCallHDF5(H5Oclose,(obj)); 707 } 708 PetscFunctionReturn(0); 709 } 710 711 /*@C 712 PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists 713 714 Input Parameters: 715 + viewer - The HDF5 viewer 716 . parent - The parent name 717 - name - The attribute name 718 719 Output Parameter: 720 . has - Flag for attribute existence 721 722 Level: advanced 723 724 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute() 725 @*/ 726 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 727 { 728 hid_t h5, dataset; 729 htri_t hhas; 730 PetscBool exists; 731 PetscErrorCode ierr; 732 733 PetscFunctionBegin; 734 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 735 PetscValidPointer(parent, 2); 736 PetscValidPointer(name, 3); 737 PetscValidPointer(has, 4); 738 *has = PETSC_FALSE; 739 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 740 ierr = PetscViewerHDF5HasObject(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr); 741 if (exists) { 742 #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 743 PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT)); 744 #else 745 PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent)); 746 #endif 747 if (dataset < 0) PetscFunctionReturn(0); 748 PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name)); 749 if (hhas < 0) { 750 PetscStackCallHDF5(H5Dclose,(dataset)); 751 PetscFunctionReturn(0); 752 } 753 PetscStackCallHDF5(H5Dclose,(dataset)); 754 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 755 } 756 PetscFunctionReturn(0); 757 } 758 759 /* 760 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 761 is attached to a communicator, in this case the attribute is a PetscViewer. 762 */ 763 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 764 765 /*@C 766 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 767 768 Collective on MPI_Comm 769 770 Input Parameter: 771 . comm - the MPI communicator to share the HDF5 PetscViewer 772 773 Level: intermediate 774 775 Options Database Keys: 776 . -viewer_hdf5_filename <name> 777 778 Environmental variables: 779 . PETSC_VIEWER_HDF5_FILENAME 780 781 Notes: 782 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 783 an error code. The HDF5 PetscViewer is usually used in the form 784 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 785 786 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 787 @*/ 788 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 789 { 790 PetscErrorCode ierr; 791 PetscBool flg; 792 PetscViewer viewer; 793 char fname[PETSC_MAX_PATH_LEN]; 794 MPI_Comm ncomm; 795 796 PetscFunctionBegin; 797 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 798 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 799 ierr = MPI_Keyval_create(MPI_NULL_COPY_FN,MPI_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 800 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 801 } 802 ierr = MPI_Attr_get(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 803 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 804 if (!flg) { /* PetscViewer not yet created */ 805 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 806 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 807 if (!flg) { 808 ierr = PetscStrcpy(fname,"output.h5"); 809 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 810 } 811 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 812 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 813 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 814 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 815 ierr = MPI_Attr_put(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 816 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 817 } 818 ierr = PetscCommDestroy(&ncomm); 819 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 820 PetscFunctionReturn(viewer); 821 } 822