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