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