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