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