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 on MPI_Comm 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 Concepts: HDF5 files 406 Concepts: PetscViewerHDF5^creating 407 408 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 409 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 410 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 411 @*/ 412 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 413 { 414 PetscErrorCode ierr; 415 416 PetscFunctionBegin; 417 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 418 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 419 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 420 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 421 ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 422 PetscFunctionReturn(0); 423 } 424 425 /*@C 426 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 427 428 Not collective 429 430 Input Parameter: 431 . viewer - the PetscViewer 432 433 Output Parameter: 434 . file_id - The file id 435 436 Level: intermediate 437 438 .seealso: PetscViewerHDF5Open() 439 @*/ 440 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 441 { 442 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 443 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 446 if (file_id) *file_id = hdf5->file_id; 447 PetscFunctionReturn(0); 448 } 449 450 /*@C 451 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 452 453 Not collective 454 455 Input Parameters: 456 + viewer - the PetscViewer 457 - name - The group name 458 459 Level: intermediate 460 461 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 462 @*/ 463 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 464 { 465 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 466 GroupList *groupNode; 467 PetscErrorCode ierr; 468 469 PetscFunctionBegin; 470 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 471 PetscValidCharPointer(name,2); 472 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 473 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 474 475 groupNode->next = hdf5->groups; 476 hdf5->groups = groupNode; 477 PetscFunctionReturn(0); 478 } 479 480 /*@ 481 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 482 483 Not collective 484 485 Input Parameter: 486 . viewer - the PetscViewer 487 488 Level: intermediate 489 490 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 491 @*/ 492 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 493 { 494 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 495 GroupList *groupNode; 496 PetscErrorCode ierr; 497 498 PetscFunctionBegin; 499 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 500 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 501 groupNode = hdf5->groups; 502 hdf5->groups = hdf5->groups->next; 503 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 504 ierr = PetscFree(groupNode);CHKERRQ(ierr); 505 PetscFunctionReturn(0); 506 } 507 508 /*@C 509 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 510 If none has been assigned, returns NULL. 511 512 Not collective 513 514 Input Parameter: 515 . viewer - the PetscViewer 516 517 Output Parameter: 518 . name - The group name 519 520 Level: intermediate 521 522 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 523 @*/ 524 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 525 { 526 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 527 528 PetscFunctionBegin; 529 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 530 PetscValidPointer(name,2); 531 if (hdf5->groups) *name = hdf5->groups->name; 532 else *name = NULL; 533 PetscFunctionReturn(0); 534 } 535 536 /*@ 537 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 538 and return this group's ID and file ID. 539 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 540 541 Not collective 542 543 Input Parameter: 544 . viewer - the PetscViewer 545 546 Output Parameter: 547 + fileId - The HDF5 file ID 548 - groupId - The HDF5 group ID 549 550 Notes: 551 If the viewer is writable, the group is created if it doesn't exist yet. 552 553 Level: intermediate 554 555 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 556 @*/ 557 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 558 { 559 hid_t file_id; 560 H5O_type_t type; 561 const char *groupName = NULL; 562 PetscBool create; 563 PetscErrorCode ierr; 564 565 PetscFunctionBegin; 566 ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 567 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 568 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 569 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 570 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 571 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 572 *fileId = file_id; 573 PetscFunctionReturn(0); 574 } 575 576 /*@ 577 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 578 579 Not collective 580 581 Input Parameter: 582 . viewer - the PetscViewer 583 584 Level: intermediate 585 586 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 587 @*/ 588 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 589 { 590 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 591 592 PetscFunctionBegin; 593 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 594 ++hdf5->timestep; 595 PetscFunctionReturn(0); 596 } 597 598 /*@ 599 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 600 of -1 disables blocking with timesteps. 601 602 Not collective 603 604 Input Parameters: 605 + viewer - the PetscViewer 606 - timestep - The timestep number 607 608 Level: intermediate 609 610 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 611 @*/ 612 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 613 { 614 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 615 616 PetscFunctionBegin; 617 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 618 hdf5->timestep = timestep; 619 PetscFunctionReturn(0); 620 } 621 622 /*@ 623 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 624 625 Not collective 626 627 Input Parameter: 628 . viewer - the PetscViewer 629 630 Output Parameter: 631 . timestep - The timestep number 632 633 Level: intermediate 634 635 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 636 @*/ 637 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 638 { 639 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 640 641 PetscFunctionBegin; 642 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 643 PetscValidPointer(timestep,2); 644 *timestep = hdf5->timestep; 645 PetscFunctionReturn(0); 646 } 647 648 /*@C 649 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 650 651 Not collective 652 653 Input Parameter: 654 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 655 656 Output Parameter: 657 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 658 659 Level: advanced 660 661 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 662 @*/ 663 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 664 { 665 PetscFunctionBegin; 666 if (ptype == PETSC_INT) 667 #if defined(PETSC_USE_64BIT_INDICES) 668 *htype = H5T_NATIVE_LLONG; 669 #else 670 *htype = H5T_NATIVE_INT; 671 #endif 672 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 673 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 674 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 675 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 676 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 677 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 678 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 679 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 680 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 681 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 682 PetscFunctionReturn(0); 683 } 684 685 /*@C 686 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 687 688 Not collective 689 690 Input Parameter: 691 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 692 693 Output Parameter: 694 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 695 696 Level: advanced 697 698 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 699 @*/ 700 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 701 { 702 PetscFunctionBegin; 703 #if defined(PETSC_USE_64BIT_INDICES) 704 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 705 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 706 #else 707 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 708 #endif 709 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 710 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 711 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 712 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 713 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 714 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 715 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 716 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 717 PetscFunctionReturn(0); 718 } 719 720 /*@C 721 PetscViewerHDF5WriteAttribute - Write an attribute 722 723 Input Parameters: 724 + viewer - The HDF5 viewer 725 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 726 . name - The attribute name 727 . datatype - The attribute type 728 - value - The attribute value 729 730 Level: advanced 731 732 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 733 @*/ 734 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 735 { 736 char *parent; 737 hid_t h5, dataspace, obj, attribute, dtype; 738 PetscBool has; 739 PetscErrorCode ierr; 740 741 PetscFunctionBegin; 742 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 743 if (dataset) PetscValidCharPointer(dataset, 2); 744 PetscValidCharPointer(name, 3); 745 PetscValidPointer(value, 5); 746 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 747 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 748 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 749 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 750 if (datatype == PETSC_STRING) { 751 size_t len; 752 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 753 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 754 } 755 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 756 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 757 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 758 if (has) { 759 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 760 } else { 761 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 762 } 763 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 764 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 765 PetscStackCallHDF5(H5Aclose,(attribute)); 766 PetscStackCallHDF5(H5Oclose,(obj)); 767 PetscStackCallHDF5(H5Sclose,(dataspace)); 768 ierr = PetscFree(parent);CHKERRQ(ierr); 769 PetscFunctionReturn(0); 770 } 771 772 /*@C 773 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 774 775 Input Parameters: 776 + viewer - The HDF5 viewer 777 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 778 . name - The attribute name 779 . datatype - The attribute type 780 - value - The attribute value 781 782 Notes: 783 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 784 You might want to check first if it does using PetscViewerHDF5HasObject(). 785 786 Level: advanced 787 788 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 789 @*/ 790 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 791 { 792 PetscErrorCode ierr; 793 794 PetscFunctionBegin; 795 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 796 PetscValidHeader(obj,2); 797 PetscValidCharPointer(name,3); 798 PetscValidPointer(value,5); 799 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 800 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 801 PetscFunctionReturn(0); 802 } 803 804 /*@C 805 PetscViewerHDF5ReadAttribute - Read an attribute 806 807 Input Parameters: 808 + viewer - The HDF5 viewer 809 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 810 . name - The attribute name 811 - datatype - The attribute type 812 813 Output Parameter: 814 . value - The attribute value 815 816 Level: advanced 817 818 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 819 @*/ 820 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 821 { 822 char *parent; 823 hid_t h5, obj, attribute, atype, dtype; 824 PetscBool has; 825 PetscErrorCode ierr; 826 827 PetscFunctionBegin; 828 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 829 if (dataset) PetscValidCharPointer(dataset, 2); 830 PetscValidCharPointer(name, 3); 831 PetscValidPointer(value, 5); 832 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 833 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 834 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 835 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 836 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 837 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 838 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 839 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 840 if (datatype == PETSC_STRING) { 841 size_t len; 842 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 843 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 844 PetscStackCallHDF5(H5Tclose,(atype)); 845 ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 846 } 847 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 848 PetscStackCallHDF5(H5Aclose,(attribute)); 849 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 850 PetscStackCallHDF5(H5Oclose,(obj)); 851 ierr = PetscFree(parent);CHKERRQ(ierr); 852 PetscFunctionReturn(0); 853 } 854 855 /*@C 856 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 857 858 Input Parameters: 859 + viewer - The HDF5 viewer 860 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 861 . name - The attribute name 862 - datatype - The attribute type 863 864 Output Parameter: 865 . value - The attribute value 866 867 Notes: 868 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 869 You might want to check first if it does using PetscViewerHDF5HasObject(). 870 871 Level: advanced 872 873 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 874 @*/ 875 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 876 { 877 PetscErrorCode ierr; 878 879 PetscFunctionBegin; 880 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 881 PetscValidHeader(obj,2); 882 PetscValidCharPointer(name,3); 883 PetscValidPointer(value, 5); 884 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 885 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 886 PetscFunctionReturn(0); 887 } 888 889 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 890 { 891 htri_t exists; 892 hid_t group; 893 894 PetscFunctionBegin; 895 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 896 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 897 if (!exists && createGroup) { 898 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 899 PetscStackCallHDF5(H5Gclose,(group)); 900 exists = PETSC_TRUE; 901 } 902 *exists_ = (PetscBool) exists; 903 PetscFunctionReturn(0); 904 } 905 906 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 907 { 908 const char rootGroupName[] = "/"; 909 hid_t h5; 910 PetscBool exists=PETSC_FALSE; 911 PetscInt i; 912 int n; 913 char **hierarchy; 914 char buf[PETSC_MAX_PATH_LEN]=""; 915 PetscErrorCode ierr; 916 917 PetscFunctionBegin; 918 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 919 if (name) PetscValidCharPointer(name, 2); 920 else name = rootGroupName; 921 if (has) { 922 PetscValidIntPointer(has, 3); 923 *has = PETSC_FALSE; 924 } 925 if (otype) { 926 PetscValidIntPointer(otype, 4); 927 *otype = H5O_TYPE_UNKNOWN; 928 } 929 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 930 931 /* 932 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 933 Hence, each of them needs to be tested separately: 934 1) whether it's a valid link 935 2) whether this link resolves to an object 936 See H5Oexists_by_name() documentation. 937 */ 938 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 939 if (!n) { 940 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 941 if (has) *has = PETSC_TRUE; 942 if (otype) *otype = H5O_TYPE_GROUP; 943 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 944 PetscFunctionReturn(0); 945 } 946 for (i=0; i<n; i++) { 947 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 948 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 949 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 950 if (!exists) break; 951 } 952 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 953 954 /* If the object exists, get its type */ 955 if (exists && otype) { 956 H5O_info_t info; 957 958 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 959 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 960 *otype = info.type; 961 } 962 if (has) *has = exists; 963 PetscFunctionReturn(0); 964 } 965 966 /*@ 967 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 968 969 Input Parameters: 970 . viewer - The HDF5 viewer 971 972 Output Parameter: 973 . has - Flag for group existence 974 975 Notes: 976 If the path exists but is not a group, this returns PETSC_FALSE as well. 977 978 Level: advanced 979 980 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 981 @*/ 982 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 983 { 984 H5O_type_t type; 985 const char *name; 986 PetscErrorCode ierr; 987 988 PetscFunctionBegin; 989 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 990 PetscValidIntPointer(has,2); 991 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 992 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 993 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 994 PetscFunctionReturn(0); 995 } 996 997 /*@ 998 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 999 1000 Input Parameters: 1001 + viewer - The HDF5 viewer 1002 - obj - The named object 1003 1004 Output Parameter: 1005 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1006 1007 Notes: 1008 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1009 1010 Level: advanced 1011 1012 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1013 @*/ 1014 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1015 { 1016 H5O_type_t type; 1017 char *path; 1018 PetscErrorCode ierr; 1019 1020 PetscFunctionBegin; 1021 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1022 PetscValidHeader(obj,2); 1023 PetscValidIntPointer(has,3); 1024 *has = PETSC_FALSE; 1025 if (!obj->name) PetscFunctionReturn(0); 1026 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1027 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1028 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1029 ierr = PetscFree(path);CHKERRQ(ierr); 1030 PetscFunctionReturn(0); 1031 } 1032 1033 /*@C 1034 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1035 1036 Input Parameters: 1037 + viewer - The HDF5 viewer 1038 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1039 - name - The attribute name 1040 1041 Output Parameter: 1042 . has - Flag for attribute existence 1043 1044 Level: advanced 1045 1046 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1047 @*/ 1048 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1049 { 1050 char *parent; 1051 PetscErrorCode ierr; 1052 1053 PetscFunctionBegin; 1054 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1055 if (dataset) PetscValidCharPointer(dataset,2); 1056 PetscValidCharPointer(name,3); 1057 PetscValidIntPointer(has,4); 1058 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1059 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1060 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1061 ierr = PetscFree(parent);CHKERRQ(ierr); 1062 PetscFunctionReturn(0); 1063 } 1064 1065 /*@C 1066 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1067 1068 Input Parameters: 1069 + viewer - The HDF5 viewer 1070 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1071 - name - The attribute name 1072 1073 Output Parameter: 1074 . has - Flag for attribute existence 1075 1076 Notes: 1077 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1078 You might want to check first if it does using PetscViewerHDF5HasObject(). 1079 1080 Level: advanced 1081 1082 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1083 @*/ 1084 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1085 { 1086 PetscErrorCode ierr; 1087 1088 PetscFunctionBegin; 1089 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1090 PetscValidHeader(obj,2); 1091 PetscValidCharPointer(name,3); 1092 PetscValidIntPointer(has,4); 1093 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1094 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1095 PetscFunctionReturn(0); 1096 } 1097 1098 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1099 { 1100 hid_t h5; 1101 htri_t hhas; 1102 PetscErrorCode ierr; 1103 1104 PetscFunctionBegin; 1105 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1106 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1107 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1108 PetscFunctionReturn(0); 1109 } 1110 1111 /* 1112 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1113 is attached to a communicator, in this case the attribute is a PetscViewer. 1114 */ 1115 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1116 1117 /*@C 1118 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1119 1120 Collective on MPI_Comm 1121 1122 Input Parameter: 1123 . comm - the MPI communicator to share the HDF5 PetscViewer 1124 1125 Level: intermediate 1126 1127 Options Database Keys: 1128 . -viewer_hdf5_filename <name> 1129 1130 Environmental variables: 1131 . PETSC_VIEWER_HDF5_FILENAME 1132 1133 Notes: 1134 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1135 an error code. The HDF5 PetscViewer is usually used in the form 1136 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1137 1138 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1139 @*/ 1140 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1141 { 1142 PetscErrorCode ierr; 1143 PetscBool flg; 1144 PetscViewer viewer; 1145 char fname[PETSC_MAX_PATH_LEN]; 1146 MPI_Comm ncomm; 1147 1148 PetscFunctionBegin; 1149 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1150 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1151 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1152 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1153 } 1154 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1155 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1156 if (!flg) { /* PetscViewer not yet created */ 1157 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1158 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1159 if (!flg) { 1160 ierr = PetscStrcpy(fname,"output.h5"); 1161 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1162 } 1163 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1164 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1165 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1166 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1167 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1168 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1169 } 1170 ierr = PetscCommDestroy(&ncomm); 1171 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1172 PetscFunctionReturn(viewer); 1173 } 1174