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