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