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