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