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