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