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