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