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