1 #include <petsc/private/viewerimpl.h> /*I "petscsys.h" I*/ 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 char *mataij_iname; 24 char *mataij_jname; 25 char *mataij_aname; 26 char *mataij_cname; 27 PetscBool mataij_names_set; 28 } PetscViewer_HDF5; 29 30 struct _n_HDF5ReadCtx { 31 hid_t file, group, dataset, dataspace, plist; 32 PetscInt timestep; 33 PetscBool complexVal, dim2, horizontal; 34 }; 35 typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 36 37 static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 38 { 39 const char *group; 40 char buf[PETSC_MAX_PATH_LEN]=""; 41 PetscErrorCode ierr; 42 43 PetscFunctionBegin; 44 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 45 ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 46 ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 47 ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 48 ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 49 PetscFunctionReturn(0); 50 } 51 52 static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 53 { 54 PetscBool has; 55 const char *group; 56 PetscErrorCode ierr; 57 58 PetscFunctionBegin; 59 if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 60 ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 61 if (!has) { 62 ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 63 SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 64 } 65 PetscFunctionReturn(0); 66 } 67 68 static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 69 { 70 PetscErrorCode ierr; 71 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 72 73 PetscFunctionBegin; 74 ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 75 ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 76 ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 77 ierr = PetscOptionsTail();CHKERRQ(ierr); 78 PetscFunctionReturn(0); 79 } 80 81 static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 82 { 83 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 84 PetscErrorCode ierr; 85 86 PetscFunctionBegin; 87 ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 88 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 89 PetscFunctionReturn(0); 90 } 91 92 static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 93 { 94 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 95 PetscErrorCode ierr; 96 97 PetscFunctionBegin; 98 ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 99 while (hdf5->groups) { 100 GroupList *tmp = hdf5->groups->next; 101 102 ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 103 ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 104 hdf5->groups = tmp; 105 } 106 ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 107 ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 108 ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 109 ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 110 ierr = PetscFree(hdf5);CHKERRQ(ierr); 111 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 112 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 113 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 114 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 115 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 116 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 117 ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 118 PetscFunctionReturn(0); 119 } 120 121 static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 122 { 123 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 124 125 PetscFunctionBegin; 126 hdf5->btype = type; 127 PetscFunctionReturn(0); 128 } 129 130 static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 131 { 132 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 133 134 PetscFunctionBegin; 135 hdf5->basedimension2 = flg; 136 PetscFunctionReturn(0); 137 } 138 139 /*@ 140 PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 141 dimension of 2. 142 143 Logically Collective on PetscViewer 144 145 Input Parameters: 146 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 147 - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 148 149 Options Database: 150 . -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 151 152 153 Notes: 154 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 155 of one when the dimension is lower. Others think the option is crazy. 156 157 Level: intermediate 158 159 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 160 161 @*/ 162 PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 163 { 164 PetscErrorCode ierr; 165 166 PetscFunctionBegin; 167 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 168 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 169 PetscFunctionReturn(0); 170 } 171 172 /*@ 173 PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 174 dimension of 2. 175 176 Logically Collective on PetscViewer 177 178 Input Parameter: 179 . viewer - the PetscViewer, must be of type HDF5 180 181 Output Parameter: 182 . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 183 184 Notes: 185 Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 186 of one when the dimension is lower. Others think the option is crazy. 187 188 Level: intermediate 189 190 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 191 192 @*/ 193 PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 194 { 195 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 196 197 PetscFunctionBegin; 198 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 199 *flg = hdf5->basedimension2; 200 PetscFunctionReturn(0); 201 } 202 203 static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 204 { 205 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 206 207 PetscFunctionBegin; 208 hdf5->spoutput = flg; 209 PetscFunctionReturn(0); 210 } 211 212 /*@ 213 PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 214 compiled with double precision PetscReal. 215 216 Logically Collective on PetscViewer 217 218 Input Parameters: 219 + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 220 - flg - if PETSC_TRUE the data will be written to disk with single precision 221 222 Options Database: 223 . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 224 225 226 Notes: 227 Setting this option does not make any difference if PETSc is compiled with single precision 228 in the first place. It does not affect reading datasets (HDF5 handle this internally). 229 230 Level: intermediate 231 232 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 233 PetscReal 234 235 @*/ 236 PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 237 { 238 PetscErrorCode ierr; 239 240 PetscFunctionBegin; 241 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 242 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 243 PetscFunctionReturn(0); 244 } 245 246 /*@ 247 PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 248 compiled with double precision PetscReal. 249 250 Logically Collective on PetscViewer 251 252 Input Parameter: 253 . viewer - the PetscViewer, must be of type HDF5 254 255 Output Parameter: 256 . flg - if PETSC_TRUE the data will be written to disk with single precision 257 258 Notes: 259 Setting this option does not make any difference if PETSc is compiled with single precision 260 in the first place. It does not affect reading datasets (HDF5 handle this internally). 261 262 Level: intermediate 263 264 .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 265 PetscReal 266 267 @*/ 268 PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 269 { 270 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 271 272 PetscFunctionBegin; 273 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 274 *flg = hdf5->spoutput; 275 PetscFunctionReturn(0); 276 } 277 278 static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 279 { 280 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 281 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 282 MPI_Info info = MPI_INFO_NULL; 283 #endif 284 hid_t plist_id; 285 PetscErrorCode ierr; 286 287 PetscFunctionBegin; 288 if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 289 if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 290 ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 291 /* Set up file access property list with parallel I/O access */ 292 PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 293 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 294 PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 295 #endif 296 /* Create or open the file collectively */ 297 switch (hdf5->btype) { 298 case FILE_MODE_READ: 299 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 300 break; 301 case FILE_MODE_APPEND: 302 PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 303 break; 304 case FILE_MODE_WRITE: 305 PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 306 break; 307 default: 308 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 309 } 310 if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 311 PetscStackCallHDF5(H5Pclose,(plist_id)); 312 PetscFunctionReturn(0); 313 } 314 315 static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 316 { 317 PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 318 319 PetscFunctionBegin; 320 *name = vhdf5->filename; 321 PetscFunctionReturn(0); 322 } 323 324 static PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 325 { 326 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 327 PetscErrorCode ierr; 328 329 PetscFunctionBegin; 330 ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 331 ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 332 ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 333 ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 334 ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 335 ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 336 ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 337 ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 338 hdf5->mataij_names_set = PETSC_TRUE; 339 PetscFunctionReturn(0); 340 } 341 342 /*@C 343 PetscViewerHDF5SetAIJNames - Set the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file. 344 345 Collective on PetscViewer 346 347 Input Parameters: 348 + viewer - the PetscViewer; either ASCII or binary 349 . iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix 350 . jname - name of dataset j representing column indices 351 . aname - name of dataset a representing matrix values 352 - cname - name of attribute stoting column count 353 354 Level: advanced 355 356 Notes: 357 Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols"). 358 For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded. 359 360 .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 361 @*/ 362 PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 363 { 364 PetscErrorCode ierr; 365 366 PetscFunctionBegin; 367 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 368 PetscValidCharPointer(iname,2); 369 PetscValidCharPointer(jname,3); 370 PetscValidCharPointer(aname,4); 371 PetscValidCharPointer(cname,5); 372 ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 373 PetscFunctionReturn(0); 374 } 375 376 static PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 377 { 378 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 379 380 PetscFunctionBegin; 381 *iname = hdf5->mataij_iname; 382 *jname = hdf5->mataij_jname; 383 *aname = hdf5->mataij_aname; 384 *cname = hdf5->mataij_cname; 385 PetscFunctionReturn(0); 386 } 387 388 /*@C 389 PetscViewerHDF5GetAIJNames - Get the names of the datasets representing the three AIJ (CRS) arrays and the name of the attribute storing the number of columns within the HDF5 file. 390 391 Collective on PetscViewer 392 393 Input Parameters: 394 . viewer - the PetscViewer; either ASCII or binary 395 396 Output Parameters: 397 + iname - name of dataset i representing row pointers; that is i[0] = 0, i[row] = i[row-1] + number of elements in that row of the matrix 398 . jname - name of dataset j representing column indices 399 . aname - name of dataset a representing matrix values 400 - cname - name of attribute stoting column count 401 402 Level: advanced 403 404 Notes: 405 Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols"). 406 For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded. 407 408 .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 409 @*/ 410 PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 411 { 412 PetscErrorCode ierr; 413 414 PetscFunctionBegin; 415 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 416 PetscValidPointer(iname,2); 417 PetscValidPointer(jname,3); 418 PetscValidPointer(aname,4); 419 PetscValidPointer(cname,5); 420 ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 421 PetscFunctionReturn(0); 422 } 423 424 static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 425 { 426 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 427 PetscErrorCode ierr; 428 429 PetscFunctionBegin; 430 if (!hdf5->mataij_names_set) { 431 if (viewer->format == PETSC_VIEWER_HDF5_MAT) { 432 ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"jc","ir","data","MATLAB_sparse");CHKERRQ(ierr); 433 } else { 434 ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"i","j","a","ncols");CHKERRQ(ierr); 435 } 436 } 437 PetscFunctionReturn(0); 438 } 439 440 /*MC 441 PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 442 443 444 .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 445 PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 446 PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 447 PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 448 449 Level: beginner 450 M*/ 451 452 PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 453 { 454 PetscViewer_HDF5 *hdf5; 455 PetscErrorCode ierr; 456 457 PetscFunctionBegin; 458 ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 459 460 v->data = (void*) hdf5; 461 v->ops->destroy = PetscViewerDestroy_HDF5; 462 v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 463 v->ops->setup = PetscViewerSetUp_HDF5; 464 v->ops->flush = 0; 465 hdf5->btype = (PetscFileMode) -1; 466 hdf5->filename = 0; 467 hdf5->timestep = -1; 468 hdf5->groups = NULL; 469 470 hdf5->mataij_iname = NULL; 471 hdf5->mataij_jname = NULL; 472 hdf5->mataij_aname = NULL; 473 hdf5->mataij_cname = NULL; 474 hdf5->mataij_names_set = PETSC_FALSE; 475 476 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 477 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 478 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 479 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 480 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 481 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 482 ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 483 PetscFunctionReturn(0); 484 } 485 486 /*@C 487 PetscViewerHDF5Open - Opens a file for HDF5 input/output. 488 489 Collective on MPI_Comm 490 491 Input Parameters: 492 + comm - MPI communicator 493 . name - name of file 494 - type - type of file 495 $ FILE_MODE_WRITE - create new file for binary output 496 $ FILE_MODE_READ - open existing file for binary input 497 $ FILE_MODE_APPEND - open existing file for binary output 498 499 Output Parameter: 500 . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 501 502 Options Database: 503 . -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 504 . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 505 506 Level: beginner 507 508 Note: 509 This PetscViewer should be destroyed with PetscViewerDestroy(). 510 511 Concepts: HDF5 files 512 Concepts: PetscViewerHDF5^creating 513 514 .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 515 PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 516 MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 517 @*/ 518 PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 519 { 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 524 ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 525 ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 526 ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 527 PetscFunctionReturn(0); 528 } 529 530 /*@C 531 PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 532 533 Not collective 534 535 Input Parameter: 536 . viewer - the PetscViewer 537 538 Output Parameter: 539 . file_id - The file id 540 541 Level: intermediate 542 543 .seealso: PetscViewerHDF5Open() 544 @*/ 545 PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 546 { 547 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 548 549 PetscFunctionBegin; 550 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 551 if (file_id) *file_id = hdf5->file_id; 552 PetscFunctionReturn(0); 553 } 554 555 /*@C 556 PetscViewerHDF5PushGroup - Set the current HDF5 group for output 557 558 Not collective 559 560 Input Parameters: 561 + viewer - the PetscViewer 562 - name - The group name 563 564 Level: intermediate 565 566 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 567 @*/ 568 PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 569 { 570 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 571 GroupList *groupNode; 572 PetscErrorCode ierr; 573 574 PetscFunctionBegin; 575 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 576 PetscValidCharPointer(name,2); 577 ierr = PetscNew(&groupNode);CHKERRQ(ierr); 578 ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 579 580 groupNode->next = hdf5->groups; 581 hdf5->groups = groupNode; 582 PetscFunctionReturn(0); 583 } 584 585 /*@ 586 PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 587 588 Not collective 589 590 Input Parameter: 591 . viewer - the PetscViewer 592 593 Level: intermediate 594 595 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 596 @*/ 597 PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 598 { 599 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 600 GroupList *groupNode; 601 PetscErrorCode ierr; 602 603 PetscFunctionBegin; 604 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 605 if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 606 groupNode = hdf5->groups; 607 hdf5->groups = hdf5->groups->next; 608 ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 609 ierr = PetscFree(groupNode);CHKERRQ(ierr); 610 PetscFunctionReturn(0); 611 } 612 613 /*@C 614 PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 615 If none has been assigned, returns NULL. 616 617 Not collective 618 619 Input Parameter: 620 . viewer - the PetscViewer 621 622 Output Parameter: 623 . name - The group name 624 625 Level: intermediate 626 627 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 628 @*/ 629 PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 630 { 631 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 632 633 PetscFunctionBegin; 634 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 635 PetscValidPointer(name,2); 636 if (hdf5->groups) *name = hdf5->groups->name; 637 else *name = NULL; 638 PetscFunctionReturn(0); 639 } 640 641 /*@ 642 PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 643 and return this group's ID and file ID. 644 If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 645 646 Not collective 647 648 Input Parameter: 649 . viewer - the PetscViewer 650 651 Output Parameter: 652 + fileId - The HDF5 file ID 653 - groupId - The HDF5 group ID 654 655 Level: intermediate 656 657 .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 658 @*/ 659 PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 660 { 661 hid_t file_id; 662 H5O_type_t type; 663 const char *groupName = NULL; 664 PetscErrorCode ierr; 665 666 PetscFunctionBegin; 667 ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 668 ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 669 ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, PETSC_TRUE, NULL, &type);CHKERRQ(ierr); 670 if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 671 PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 672 *fileId = file_id; 673 PetscFunctionReturn(0); 674 } 675 676 /*@ 677 PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 678 679 Not collective 680 681 Input Parameter: 682 . viewer - the PetscViewer 683 684 Level: intermediate 685 686 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 687 @*/ 688 PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 689 { 690 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 691 692 PetscFunctionBegin; 693 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 694 ++hdf5->timestep; 695 PetscFunctionReturn(0); 696 } 697 698 /*@ 699 PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 700 of -1 disables blocking with timesteps. 701 702 Not collective 703 704 Input Parameters: 705 + viewer - the PetscViewer 706 - timestep - The timestep number 707 708 Level: intermediate 709 710 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 711 @*/ 712 PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 713 { 714 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 715 716 PetscFunctionBegin; 717 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 718 hdf5->timestep = timestep; 719 PetscFunctionReturn(0); 720 } 721 722 /*@ 723 PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 724 725 Not collective 726 727 Input Parameter: 728 . viewer - the PetscViewer 729 730 Output Parameter: 731 . timestep - The timestep number 732 733 Level: intermediate 734 735 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 736 @*/ 737 PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 738 { 739 PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 740 741 PetscFunctionBegin; 742 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 743 PetscValidPointer(timestep,2); 744 *timestep = hdf5->timestep; 745 PetscFunctionReturn(0); 746 } 747 748 /*@C 749 PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 750 751 Not collective 752 753 Input Parameter: 754 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 755 756 Output Parameter: 757 . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 758 759 Level: advanced 760 761 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 762 @*/ 763 PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 764 { 765 PetscFunctionBegin; 766 if (ptype == PETSC_INT) 767 #if defined(PETSC_USE_64BIT_INDICES) 768 *htype = H5T_NATIVE_LLONG; 769 #else 770 *htype = H5T_NATIVE_INT; 771 #endif 772 else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 773 else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 774 else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 775 else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 776 else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 777 else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 778 else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 779 else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 780 else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 781 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 782 PetscFunctionReturn(0); 783 } 784 785 /*@C 786 PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 787 788 Not collective 789 790 Input Parameter: 791 . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 792 793 Output Parameter: 794 . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 795 796 Level: advanced 797 798 .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 799 @*/ 800 PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 801 { 802 PetscFunctionBegin; 803 #if defined(PETSC_USE_64BIT_INDICES) 804 if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 805 else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 806 #else 807 if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 808 #endif 809 else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 810 else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 811 else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 812 else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 813 else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 814 else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 815 else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 816 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 817 PetscFunctionReturn(0); 818 } 819 820 /*@C 821 PetscViewerHDF5WriteAttribute - Write an attribute 822 823 Input Parameters: 824 + viewer - The HDF5 viewer 825 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 826 . name - The attribute name 827 . datatype - The attribute type 828 - value - The attribute value 829 830 Level: advanced 831 832 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 833 @*/ 834 PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 835 { 836 char *parent; 837 hid_t h5, dataspace, obj, attribute, dtype; 838 PetscBool has; 839 PetscErrorCode ierr; 840 841 PetscFunctionBegin; 842 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 843 if (dataset) PetscValidCharPointer(dataset, 2); 844 PetscValidCharPointer(name, 3); 845 PetscValidPointer(value, 5); 846 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 847 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 848 ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 849 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 850 if (datatype == PETSC_STRING) { 851 size_t len; 852 ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 853 PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 854 } 855 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 856 PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 857 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 858 if (has) { 859 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 860 } else { 861 PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 862 } 863 PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 864 if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 865 PetscStackCallHDF5(H5Aclose,(attribute)); 866 PetscStackCallHDF5(H5Oclose,(obj)); 867 PetscStackCallHDF5(H5Sclose,(dataspace)); 868 ierr = PetscFree(parent);CHKERRQ(ierr); 869 PetscFunctionReturn(0); 870 } 871 872 /*@C 873 PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 874 875 Input Parameters: 876 + viewer - The HDF5 viewer 877 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 878 . name - The attribute name 879 . datatype - The attribute type 880 - value - The attribute value 881 882 Notes: 883 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 884 You might want to check first if it does using PetscViewerHDF5HasObject(). 885 886 Level: advanced 887 888 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 889 @*/ 890 PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 891 { 892 PetscErrorCode ierr; 893 894 PetscFunctionBegin; 895 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 896 PetscValidHeader(obj,2); 897 PetscValidCharPointer(name,3); 898 PetscValidPointer(value,5); 899 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 900 ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 901 PetscFunctionReturn(0); 902 } 903 904 /*@C 905 PetscViewerHDF5ReadAttribute - Read an attribute 906 907 Input Parameters: 908 + viewer - The HDF5 viewer 909 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 910 . name - The attribute name 911 - datatype - The attribute type 912 913 Output Parameter: 914 . value - The attribute value 915 916 Level: advanced 917 918 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 919 @*/ 920 PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 921 { 922 char *parent; 923 hid_t h5, obj, attribute, atype, dtype; 924 PetscBool has; 925 PetscErrorCode ierr; 926 927 PetscFunctionBegin; 928 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 929 if (dataset) PetscValidCharPointer(dataset, 2); 930 PetscValidCharPointer(name, 3); 931 PetscValidPointer(value, 5); 932 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 933 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 934 if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 935 if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 936 ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 937 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 938 PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 939 PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 940 if (datatype == PETSC_STRING) { 941 size_t len; 942 PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 943 PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 944 PetscStackCallHDF5(H5Tclose,(atype)); 945 ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 946 } 947 PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 948 PetscStackCallHDF5(H5Aclose,(attribute)); 949 /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 950 PetscStackCallHDF5(H5Oclose,(obj)); 951 ierr = PetscFree(parent);CHKERRQ(ierr); 952 PetscFunctionReturn(0); 953 } 954 955 /*@C 956 PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 957 958 Input Parameters: 959 + viewer - The HDF5 viewer 960 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 961 . name - The attribute name 962 - datatype - The attribute type 963 964 Output Parameter: 965 . value - The attribute value 966 967 Notes: 968 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 969 You might want to check first if it does using PetscViewerHDF5HasObject(). 970 971 Level: advanced 972 973 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 974 @*/ 975 PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 976 { 977 PetscErrorCode ierr; 978 979 PetscFunctionBegin; 980 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 981 PetscValidHeader(obj,2); 982 PetscValidCharPointer(name,3); 983 PetscValidPointer(value, 5); 984 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 985 ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 986 PetscFunctionReturn(0); 987 } 988 989 PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 990 { 991 htri_t exists; 992 hid_t group; 993 994 PetscFunctionBegin; 995 PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 996 if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 997 if (!exists && createGroup) { 998 PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 999 PetscStackCallHDF5(H5Gclose,(group)); 1000 exists = PETSC_TRUE; 1001 } 1002 *exists_ = (PetscBool) exists; 1003 PetscFunctionReturn(0); 1004 } 1005 1006 static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 1007 { 1008 const char rootGroupName[] = "/"; 1009 hid_t h5; 1010 PetscBool exists=PETSC_FALSE; 1011 PetscInt i; 1012 int n; 1013 char **hierarchy; 1014 char buf[PETSC_MAX_PATH_LEN]=""; 1015 PetscErrorCode ierr; 1016 1017 PetscFunctionBegin; 1018 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1019 if (name) PetscValidCharPointer(name, 2); 1020 else name = rootGroupName; 1021 if (has) { 1022 PetscValidIntPointer(has, 3); 1023 *has = PETSC_FALSE; 1024 } 1025 if (otype) { 1026 PetscValidIntPointer(otype, 4); 1027 *otype = H5O_TYPE_UNKNOWN; 1028 } 1029 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1030 1031 /* 1032 Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 1033 Hence, each of them needs to be tested separately: 1034 1) whether it's a valid link 1035 2) whether this link resolves to an object 1036 See H5Oexists_by_name() documentation. 1037 */ 1038 ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 1039 if (!n) { 1040 /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1041 if (has) *has = PETSC_TRUE; 1042 if (otype) *otype = H5O_TYPE_GROUP; 1043 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1044 PetscFunctionReturn(0); 1045 } 1046 for (i=0; i<n; i++) { 1047 ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 1048 ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1049 ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1050 if (!exists) break; 1051 } 1052 ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 1053 1054 /* If the object exists, get its type */ 1055 if (exists && otype) { 1056 H5O_info_t info; 1057 1058 /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 1059 PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 1060 *otype = info.type; 1061 } 1062 if (has) *has = exists; 1063 PetscFunctionReturn(0); 1064 } 1065 1066 /*@ 1067 PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 1068 1069 Input Parameters: 1070 . viewer - The HDF5 viewer 1071 1072 Output Parameter: 1073 . has - Flag for group existence 1074 1075 Notes: 1076 If the path exists but is not a group, this returns PETSC_FALSE as well. 1077 1078 Level: advanced 1079 1080 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 1081 @*/ 1082 PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 1083 { 1084 H5O_type_t type; 1085 const char *name; 1086 PetscErrorCode ierr; 1087 1088 PetscFunctionBegin; 1089 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1090 PetscValidIntPointer(has,2); 1091 ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 1092 ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1093 *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 1094 PetscFunctionReturn(0); 1095 } 1096 1097 /*@ 1098 PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1099 1100 Input Parameters: 1101 + viewer - The HDF5 viewer 1102 - obj - The named object 1103 1104 Output Parameter: 1105 . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1106 1107 Notes: 1108 If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1109 1110 Level: advanced 1111 1112 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1113 @*/ 1114 PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1115 { 1116 H5O_type_t type; 1117 char *path; 1118 PetscErrorCode ierr; 1119 1120 PetscFunctionBegin; 1121 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1122 PetscValidHeader(obj,2); 1123 PetscValidIntPointer(has,3); 1124 *has = PETSC_FALSE; 1125 if (!obj->name) PetscFunctionReturn(0); 1126 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 1127 ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 1128 *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 1129 ierr = PetscFree(path);CHKERRQ(ierr); 1130 PetscFunctionReturn(0); 1131 } 1132 1133 /*@C 1134 PetscViewerHDF5HasAttribute - Check whether an attribute exists 1135 1136 Input Parameters: 1137 + viewer - The HDF5 viewer 1138 . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1139 - name - The attribute name 1140 1141 Output Parameter: 1142 . has - Flag for attribute existence 1143 1144 Level: advanced 1145 1146 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1147 @*/ 1148 PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1149 { 1150 char *parent; 1151 PetscErrorCode ierr; 1152 1153 PetscFunctionBegin; 1154 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1155 if (dataset) PetscValidCharPointer(dataset,2); 1156 PetscValidCharPointer(name,3); 1157 PetscValidIntPointer(has,4); 1158 ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1159 ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 1160 if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 1161 ierr = PetscFree(parent);CHKERRQ(ierr); 1162 PetscFunctionReturn(0); 1163 } 1164 1165 /*@C 1166 PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1167 1168 Input Parameters: 1169 + viewer - The HDF5 viewer 1170 . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1171 - name - The attribute name 1172 1173 Output Parameter: 1174 . has - Flag for attribute existence 1175 1176 Notes: 1177 This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1178 You might want to check first if it does using PetscViewerHDF5HasObject(). 1179 1180 Level: advanced 1181 1182 .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1183 @*/ 1184 PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1185 { 1186 PetscErrorCode ierr; 1187 1188 PetscFunctionBegin; 1189 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1190 PetscValidHeader(obj,2); 1191 PetscValidCharPointer(name,3); 1192 PetscValidIntPointer(has,4); 1193 ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1194 ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1195 PetscFunctionReturn(0); 1196 } 1197 1198 static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1199 { 1200 hid_t h5; 1201 htri_t hhas; 1202 PetscErrorCode ierr; 1203 1204 PetscFunctionBegin; 1205 ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1206 PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 1207 *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1208 PetscFunctionReturn(0); 1209 } 1210 1211 static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1212 { 1213 HDF5ReadCtx h=NULL; 1214 PetscErrorCode ierr; 1215 1216 PetscFunctionBegin; 1217 ierr = PetscNew(&h);CHKERRQ(ierr); 1218 ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 1219 PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 1220 PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 1221 ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 1222 ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr); 1223 if (h->complexVal) {ierr = PetscViewerHDF5ReadAttribute(viewer,name,"complex",PETSC_BOOL,&h->complexVal);CHKERRQ(ierr);} 1224 /* MATLAB stores column vectors horizontally */ 1225 ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1226 /* Create property list for collective dataset read */ 1227 PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1228 #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1229 PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1230 #endif 1231 /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 1232 *ctx = h; 1233 PetscFunctionReturn(0); 1234 } 1235 1236 static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 1237 { 1238 HDF5ReadCtx h; 1239 PetscErrorCode ierr; 1240 1241 PetscFunctionBegin; 1242 h = *ctx; 1243 PetscStackCallHDF5(H5Pclose,(h->plist)); 1244 PetscStackCallHDF5(H5Gclose,(h->group)); 1245 PetscStackCallHDF5(H5Sclose,(h->dataspace)); 1246 PetscStackCallHDF5(H5Dclose,(h->dataset)); 1247 ierr = PetscFree(*ctx);CHKERRQ(ierr); 1248 PetscFunctionReturn(0); 1249 } 1250 1251 static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 1252 { 1253 int rdim, dim; 1254 hsize_t dims[4]; 1255 PetscInt bsInd, lenInd, bs, len, N; 1256 PetscLayout map; 1257 PetscErrorCode ierr; 1258 1259 PetscFunctionBegin; 1260 if (!(*map_)) { 1261 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 1262 } 1263 map = *map_; 1264 /* calculate expected number of dimensions */ 1265 dim = 0; 1266 if (ctx->timestep >= 0) ++dim; 1267 ++dim; /* length in blocks */ 1268 if (ctx->complexVal) ++dim; 1269 /* get actual number of dimensions in dataset */ 1270 PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 1271 /* calculate expected dimension indices */ 1272 lenInd = 0; 1273 if (ctx->timestep >= 0) ++lenInd; 1274 bsInd = lenInd + 1; 1275 ctx->dim2 = PETSC_FALSE; 1276 if (rdim == dim) { 1277 bs = 1; /* support vectors stored as 1D array */ 1278 } else if (rdim == dim+1) { 1279 bs = (PetscInt) dims[bsInd]; 1280 if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1281 } else { 1282 SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim); 1283 } 1284 len = dims[lenInd]; 1285 if (ctx->horizontal) { 1286 if (len != 1) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Cannot have horizontal array with number of rows > 1. In case of MATLAB MAT-file, vectors must be saved as column vectors."); 1287 len = bs; 1288 bs = 1; 1289 } 1290 N = (PetscInt) len*bs; 1291 1292 /* Set Vec sizes,blocksize,and type if not already set */ 1293 if (map->bs < 0) { 1294 ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 1295 } else if (map->bs != bs) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_FILE_UNEXPECTED, "Block size of array in file is %D, not %D as expected",bs,map->bs); 1296 if (map->N < 0) { 1297 ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 1298 } else if (map->N != N) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %D, not %D as expected",N,map->N); 1299 PetscFunctionReturn(0); 1300 } 1301 1302 static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 1303 { 1304 hsize_t count[4], offset[4]; 1305 int dim; 1306 PetscInt bs, n, low; 1307 PetscErrorCode ierr; 1308 1309 PetscFunctionBegin; 1310 /* Compute local size and ownership range */ 1311 ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1312 ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 1313 ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 1314 ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 1315 1316 /* Each process defines a dataset and reads it from the hyperslab in the file */ 1317 dim = 0; 1318 if (ctx->timestep >= 0) { 1319 count[dim] = 1; 1320 offset[dim] = ctx->timestep; 1321 ++dim; 1322 } 1323 if (ctx->horizontal) { 1324 count[dim] = 1; 1325 offset[dim] = 0; 1326 ++dim; 1327 } 1328 { 1329 ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 1330 ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 1331 ++dim; 1332 } 1333 if (bs > 1 || ctx->dim2) { 1334 if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 1335 count[dim] = bs; 1336 offset[dim] = 0; 1337 ++dim; 1338 } 1339 if (ctx->complexVal) { 1340 count[dim] = 2; 1341 offset[dim] = 0; 1342 ++dim; 1343 } 1344 PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 1345 PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 1346 PetscFunctionReturn(0); 1347 } 1348 1349 static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1350 { 1351 PetscFunctionBegin; 1352 PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1353 PetscFunctionReturn(0); 1354 } 1355 1356 PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1357 { 1358 HDF5ReadCtx h=NULL; 1359 hid_t memspace=0; 1360 size_t unitsize; 1361 void *arr; 1362 PetscErrorCode ierr; 1363 1364 PetscFunctionBegin; 1365 ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1366 #if defined(PETSC_USE_COMPLEX) 1367 if (!h->complexVal) { 1368 H5T_class_t clazz = H5Tget_class(datatype); 1369 if (clazz == H5T_FLOAT) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains real numbers but PETSc is configured for complex. The conversion is not yet implemented. Configure with --with-scalar-type=real."); 1370 } 1371 #else 1372 if (h->complexVal) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"File contains complex numbers but PETSc not configured for them. Configure with --with-scalar-type=complex."); 1373 #endif 1374 1375 ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1376 ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1377 ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 1378 1379 unitsize = H5Tget_size(datatype); 1380 if (h->complexVal) unitsize *= 2; 1381 if (unitsize <= 0 || unitsize > PetscMax(sizeof(PetscInt),sizeof(PetscScalar))) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"Sanity check failed: HDF5 function H5Tget_size(datatype) returned suspicious value %D",unitsize); 1382 ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 1383 1384 ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1385 PetscStackCallHDF5(H5Sclose,(memspace)); 1386 ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1387 *newarr = arr; 1388 PetscFunctionReturn(0); 1389 } 1390 1391 /*@C 1392 PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1393 1394 Input Parameters: 1395 + viewer - The HDF5 viewer 1396 - name - The vector name 1397 1398 Output Parameter: 1399 + bs - block size 1400 - N - global size 1401 1402 Note: 1403 A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1404 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1405 1406 A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1407 1408 Level: advanced 1409 1410 .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1411 @*/ 1412 PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 1413 { 1414 HDF5ReadCtx h=NULL; 1415 PetscLayout map=NULL; 1416 PetscErrorCode ierr; 1417 1418 PetscFunctionBegin; 1419 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1420 ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1421 ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1422 ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1423 if (bs) *bs = map->bs; 1424 if (N) *N = map->N; 1425 ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1426 PetscFunctionReturn(0); 1427 } 1428 1429 /* 1430 The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1431 is attached to a communicator, in this case the attribute is a PetscViewer. 1432 */ 1433 PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1434 1435 /*@C 1436 PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1437 1438 Collective on MPI_Comm 1439 1440 Input Parameter: 1441 . comm - the MPI communicator to share the HDF5 PetscViewer 1442 1443 Level: intermediate 1444 1445 Options Database Keys: 1446 . -viewer_hdf5_filename <name> 1447 1448 Environmental variables: 1449 . PETSC_VIEWER_HDF5_FILENAME 1450 1451 Notes: 1452 Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1453 an error code. The HDF5 PetscViewer is usually used in the form 1454 $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1455 1456 .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1457 @*/ 1458 PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1459 { 1460 PetscErrorCode ierr; 1461 PetscBool flg; 1462 PetscViewer viewer; 1463 char fname[PETSC_MAX_PATH_LEN]; 1464 MPI_Comm ncomm; 1465 1466 PetscFunctionBegin; 1467 ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1468 if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1469 ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1470 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1471 } 1472 ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1473 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1474 if (!flg) { /* PetscViewer not yet created */ 1475 ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1476 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1477 if (!flg) { 1478 ierr = PetscStrcpy(fname,"output.h5"); 1479 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1480 } 1481 ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1482 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1483 ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1484 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1485 ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1486 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1487 } 1488 ierr = PetscCommDestroy(&ncomm); 1489 if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1490 PetscFunctionReturn(viewer); 1491 } 1492