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