1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I "petscsys.h" I*/ 2d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 3ded06c3fSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800) 4ded06c3fSVaclav Hapla #error "PETSc needs HDF5 version >= 1.8.0" 5ded06c3fSVaclav Hapla #endif 65c6c1daeSBarry Smith 7bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 906db490cSVaclav Hapla 105c6c1daeSBarry Smith typedef struct GroupList { 115c6c1daeSBarry Smith const char *name; 125c6c1daeSBarry Smith struct GroupList *next; 135c6c1daeSBarry Smith } GroupList; 145c6c1daeSBarry Smith 155c6c1daeSBarry Smith typedef struct { 165c6c1daeSBarry Smith char *filename; 175c6c1daeSBarry Smith PetscFileMode btype; 185c6c1daeSBarry Smith hid_t file_id; 195c6c1daeSBarry Smith PetscInt timestep; 205c6c1daeSBarry Smith GroupList *groups; 2182ea9b62SBarry Smith PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 229a0502c6SHåkon Strandenes PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 23058bd781SVaclav Hapla char *mataij_iname; 24058bd781SVaclav Hapla char *mataij_jname; 25058bd781SVaclav Hapla char *mataij_aname; 26058bd781SVaclav Hapla char *mataij_cname; 275c6c1daeSBarry Smith } PetscViewer_HDF5; 285c6c1daeSBarry Smith 29eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx { 30eb5a92b4SVaclav Hapla hid_t file, group, dataset, dataspace, plist; 31eb5a92b4SVaclav Hapla PetscInt timestep; 3209dabeb0SVaclav Hapla PetscBool complexVal, dim2, horizontal; 33eb5a92b4SVaclav Hapla }; 34eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 35eb5a92b4SVaclav Hapla 364416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 3782ea9b62SBarry Smith { 3882ea9b62SBarry Smith PetscErrorCode ierr; 3982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4082ea9b62SBarry Smith 4182ea9b62SBarry Smith PetscFunctionBegin; 4282ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 4382ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 449a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 4582ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4682ea9b62SBarry Smith PetscFunctionReturn(0); 4782ea9b62SBarry Smith } 4882ea9b62SBarry Smith 495c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 505c6c1daeSBarry Smith { 515c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 525c6c1daeSBarry Smith PetscErrorCode ierr; 535c6c1daeSBarry Smith 545c6c1daeSBarry Smith PetscFunctionBegin; 555c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 56729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 575c6c1daeSBarry Smith PetscFunctionReturn(0); 585c6c1daeSBarry Smith } 595c6c1daeSBarry Smith 605c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 615c6c1daeSBarry Smith { 625c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 635c6c1daeSBarry Smith PetscErrorCode ierr; 645c6c1daeSBarry Smith 655c6c1daeSBarry Smith PetscFunctionBegin; 665c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 675c6c1daeSBarry Smith while (hdf5->groups) { 685c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 695c6c1daeSBarry Smith 705c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 715c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 725c6c1daeSBarry Smith hdf5->groups = tmp; 735c6c1daeSBarry Smith } 7461912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 7561912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 7661912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 7761912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 785c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 790b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 80d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 810b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 82058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 83058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 84058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 85058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 865c6c1daeSBarry Smith PetscFunctionReturn(0); 875c6c1daeSBarry Smith } 885c6c1daeSBarry Smith 895c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 905c6c1daeSBarry Smith { 915c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 925c6c1daeSBarry Smith 935c6c1daeSBarry Smith PetscFunctionBegin; 945c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 955c6c1daeSBarry Smith hdf5->btype = type; 965c6c1daeSBarry Smith PetscFunctionReturn(0); 975c6c1daeSBarry Smith } 985c6c1daeSBarry Smith 9982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 10082ea9b62SBarry Smith { 10182ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 10282ea9b62SBarry Smith 10382ea9b62SBarry Smith PetscFunctionBegin; 10482ea9b62SBarry Smith hdf5->basedimension2 = flg; 10582ea9b62SBarry Smith PetscFunctionReturn(0); 10682ea9b62SBarry Smith } 10782ea9b62SBarry Smith 108df863907SAlex Fikl /*@ 10982ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 11082ea9b62SBarry Smith dimension of 2. 11182ea9b62SBarry Smith 11282ea9b62SBarry Smith Logically Collective on PetscViewer 11382ea9b62SBarry Smith 11482ea9b62SBarry Smith Input Parameters: 11582ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 11682ea9b62SBarry Smith - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 11782ea9b62SBarry Smith 11882ea9b62SBarry Smith Options Database: 11982ea9b62SBarry Smith . -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 12082ea9b62SBarry Smith 12182ea9b62SBarry Smith 12295452b02SPatrick Sanan Notes: 12395452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 12482ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 12582ea9b62SBarry Smith 12682ea9b62SBarry Smith Level: intermediate 12782ea9b62SBarry Smith 12882ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 12982ea9b62SBarry Smith 13082ea9b62SBarry Smith @*/ 13182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 13282ea9b62SBarry Smith { 13382ea9b62SBarry Smith PetscErrorCode ierr; 13482ea9b62SBarry Smith 13582ea9b62SBarry Smith PetscFunctionBegin; 13682ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 13782ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 13882ea9b62SBarry Smith PetscFunctionReturn(0); 13982ea9b62SBarry Smith } 14082ea9b62SBarry Smith 141df863907SAlex Fikl /*@ 14282ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14382ea9b62SBarry Smith dimension of 2. 14482ea9b62SBarry Smith 14582ea9b62SBarry Smith Logically Collective on PetscViewer 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith Input Parameter: 14882ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith Output Parameter: 15182ea9b62SBarry Smith . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 15282ea9b62SBarry Smith 15395452b02SPatrick Sanan Notes: 15495452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 15582ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15682ea9b62SBarry Smith 15782ea9b62SBarry Smith Level: intermediate 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16082ea9b62SBarry Smith 16182ea9b62SBarry Smith @*/ 16282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 16382ea9b62SBarry Smith { 16482ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 16582ea9b62SBarry Smith 16682ea9b62SBarry Smith PetscFunctionBegin; 16782ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16882ea9b62SBarry Smith *flg = hdf5->basedimension2; 16982ea9b62SBarry Smith PetscFunctionReturn(0); 17082ea9b62SBarry Smith } 17182ea9b62SBarry Smith 1729a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1739a0502c6SHåkon Strandenes { 1749a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1759a0502c6SHåkon Strandenes 1769a0502c6SHåkon Strandenes PetscFunctionBegin; 1779a0502c6SHåkon Strandenes hdf5->spoutput = flg; 1789a0502c6SHåkon Strandenes PetscFunctionReturn(0); 1799a0502c6SHåkon Strandenes } 1809a0502c6SHåkon Strandenes 181df863907SAlex Fikl /*@ 1829a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 1839a0502c6SHåkon Strandenes compiled with double precision PetscReal. 1849a0502c6SHåkon Strandenes 1859a0502c6SHåkon Strandenes Logically Collective on PetscViewer 1869a0502c6SHåkon Strandenes 1879a0502c6SHåkon Strandenes Input Parameters: 1889a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 1899a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 1909a0502c6SHåkon Strandenes 1919a0502c6SHåkon Strandenes Options Database: 1929a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 1939a0502c6SHåkon Strandenes 1949a0502c6SHåkon Strandenes 19595452b02SPatrick Sanan Notes: 19695452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 1979a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 1989a0502c6SHåkon Strandenes 1999a0502c6SHåkon Strandenes Level: intermediate 2009a0502c6SHåkon Strandenes 2019a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2029a0502c6SHåkon Strandenes PetscReal 2039a0502c6SHåkon Strandenes 2049a0502c6SHåkon Strandenes @*/ 2059a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2069a0502c6SHåkon Strandenes { 2079a0502c6SHåkon Strandenes PetscErrorCode ierr; 2089a0502c6SHåkon Strandenes 2099a0502c6SHåkon Strandenes PetscFunctionBegin; 2109a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2119a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2129a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2139a0502c6SHåkon Strandenes } 2149a0502c6SHåkon Strandenes 215df863907SAlex Fikl /*@ 2169a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2179a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2189a0502c6SHåkon Strandenes 2199a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2209a0502c6SHåkon Strandenes 2219a0502c6SHåkon Strandenes Input Parameter: 2229a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes Output Parameter: 2259a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2269a0502c6SHåkon Strandenes 22795452b02SPatrick Sanan Notes: 22895452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2299a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2309a0502c6SHåkon Strandenes 2319a0502c6SHåkon Strandenes Level: intermediate 2329a0502c6SHåkon Strandenes 2339a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2349a0502c6SHåkon Strandenes PetscReal 2359a0502c6SHåkon Strandenes 2369a0502c6SHåkon Strandenes @*/ 2379a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2389a0502c6SHåkon Strandenes { 2399a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2409a0502c6SHåkon Strandenes 2419a0502c6SHåkon Strandenes PetscFunctionBegin; 2429a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2439a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2449a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2459a0502c6SHåkon Strandenes } 2469a0502c6SHåkon Strandenes 2475c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2485c6c1daeSBarry Smith { 2495c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2505c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2515c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2525c6c1daeSBarry Smith #endif 2535c6c1daeSBarry Smith hid_t plist_id; 2545c6c1daeSBarry Smith PetscErrorCode ierr; 2555c6c1daeSBarry Smith 2565c6c1daeSBarry Smith PetscFunctionBegin; 257f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 258f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2595c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2605c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 261729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2625c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 263729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2645c6c1daeSBarry Smith #endif 2655c6c1daeSBarry Smith /* Create or open the file collectively */ 2665c6c1daeSBarry Smith switch (hdf5->btype) { 2675c6c1daeSBarry Smith case FILE_MODE_READ: 268729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 2695c6c1daeSBarry Smith break; 2705c6c1daeSBarry Smith case FILE_MODE_APPEND: 271729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 2725c6c1daeSBarry Smith break; 2735c6c1daeSBarry Smith case FILE_MODE_WRITE: 274729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 2755c6c1daeSBarry Smith break; 2765c6c1daeSBarry Smith default: 2775c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 2785c6c1daeSBarry Smith } 2795c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 280729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 2815c6c1daeSBarry Smith PetscFunctionReturn(0); 2825c6c1daeSBarry Smith } 2835c6c1daeSBarry Smith 284d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 285d1232d7fSToby Isaac { 286d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 287d1232d7fSToby Isaac 288d1232d7fSToby Isaac PetscFunctionBegin; 289d1232d7fSToby Isaac *name = vhdf5->filename; 290d1232d7fSToby Isaac PetscFunctionReturn(0); 291d1232d7fSToby Isaac } 292d1232d7fSToby Isaac 293058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 294058bd781SVaclav Hapla { 295058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 296058bd781SVaclav Hapla PetscErrorCode ierr; 297058bd781SVaclav Hapla 298058bd781SVaclav Hapla PetscFunctionBegin; 299058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 300058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 301058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 302058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 303058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 304058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 305058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 306058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 307058bd781SVaclav Hapla PetscFunctionReturn(0); 308058bd781SVaclav Hapla } 309058bd781SVaclav Hapla 310058bd781SVaclav Hapla /*@C 311058bd781SVaclav Hapla 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. 312058bd781SVaclav Hapla 313058bd781SVaclav Hapla Collective on PetscViewer 314058bd781SVaclav Hapla 315058bd781SVaclav Hapla Input Parameters: 316058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 317058bd781SVaclav Hapla . 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 318058bd781SVaclav Hapla . jname - name of dataset j representing column indices 319058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 320058bd781SVaclav Hapla - cname - name of attribute stoting column count 321058bd781SVaclav Hapla 322058bd781SVaclav Hapla Level: advanced 323058bd781SVaclav Hapla 324058bd781SVaclav Hapla Notes: 325058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 326058bd781SVaclav Hapla 327058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 328058bd781SVaclav Hapla @*/ 329058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 330058bd781SVaclav Hapla { 331058bd781SVaclav Hapla PetscErrorCode ierr; 332058bd781SVaclav Hapla 333058bd781SVaclav Hapla PetscFunctionBegin; 334058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 335058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 336058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 337058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 338058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 339058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 340058bd781SVaclav Hapla PetscFunctionReturn(0); 341058bd781SVaclav Hapla } 342058bd781SVaclav Hapla 343058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 344058bd781SVaclav Hapla { 345058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 346058bd781SVaclav Hapla 347058bd781SVaclav Hapla PetscFunctionBegin; 348058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 349058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 350058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 351058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 352058bd781SVaclav Hapla PetscFunctionReturn(0); 353058bd781SVaclav Hapla } 354058bd781SVaclav Hapla 355058bd781SVaclav Hapla /*@C 356058bd781SVaclav Hapla 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. 357058bd781SVaclav Hapla 358058bd781SVaclav Hapla Collective on PetscViewer 359058bd781SVaclav Hapla 360058bd781SVaclav Hapla Input Parameters: 361058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 362058bd781SVaclav Hapla 363058bd781SVaclav Hapla Output Parameters: 364058bd781SVaclav Hapla + 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 365058bd781SVaclav Hapla . jname - name of dataset j representing column indices 366058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 367058bd781SVaclav Hapla - cname - name of attribute stoting column count 368058bd781SVaclav Hapla 369058bd781SVaclav Hapla Level: advanced 370058bd781SVaclav Hapla 371058bd781SVaclav Hapla Notes: 372058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 373058bd781SVaclav Hapla 374058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 375058bd781SVaclav Hapla @*/ 376058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 377058bd781SVaclav Hapla { 378058bd781SVaclav Hapla PetscErrorCode ierr; 379058bd781SVaclav Hapla 380058bd781SVaclav Hapla PetscFunctionBegin; 381058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 382058bd781SVaclav Hapla PetscValidPointer(iname,2); 383058bd781SVaclav Hapla PetscValidPointer(jname,3); 384058bd781SVaclav Hapla PetscValidPointer(aname,4); 385058bd781SVaclav Hapla PetscValidPointer(cname,5); 386058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 387058bd781SVaclav Hapla PetscFunctionReturn(0); 388058bd781SVaclav Hapla } 389058bd781SVaclav Hapla 3908556b5ebSBarry Smith /*MC 3918556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 3928556b5ebSBarry Smith 3938556b5ebSBarry Smith 3948556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 3958556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 3968556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 3978556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 3988556b5ebSBarry Smith 3991b266c99SBarry Smith Level: beginner 4008556b5ebSBarry Smith M*/ 401d1232d7fSToby Isaac 4028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4035c6c1daeSBarry Smith { 4045c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4055c6c1daeSBarry Smith PetscErrorCode ierr; 4065c6c1daeSBarry Smith 4075c6c1daeSBarry Smith PetscFunctionBegin; 408b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4095c6c1daeSBarry Smith 4105c6c1daeSBarry Smith v->data = (void*) hdf5; 4115c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 41282ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 4135c6c1daeSBarry Smith v->ops->flush = 0; 4145c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4155c6c1daeSBarry Smith hdf5->filename = 0; 4165c6c1daeSBarry Smith hdf5->timestep = -1; 4170298fd71SBarry Smith hdf5->groups = NULL; 4185c6c1daeSBarry Smith 419058bd781SVaclav Hapla /* ir and jc are deliberately swapped as MATLAB uses column-major format */ 420058bd781SVaclav Hapla ierr = PetscStrallocpy("jc", &hdf5->mataij_iname);CHKERRQ(ierr); 421058bd781SVaclav Hapla ierr = PetscStrallocpy("ir", &hdf5->mataij_jname);CHKERRQ(ierr); 422058bd781SVaclav Hapla ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr); 423058bd781SVaclav Hapla ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr); 424058bd781SVaclav Hapla 4250b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 426d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4270b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 42882ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4299a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 430058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 431058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4325c6c1daeSBarry Smith PetscFunctionReturn(0); 4335c6c1daeSBarry Smith } 4345c6c1daeSBarry Smith 4355c6c1daeSBarry Smith /*@C 4365c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4375c6c1daeSBarry Smith 4385c6c1daeSBarry Smith Collective on MPI_Comm 4395c6c1daeSBarry Smith 4405c6c1daeSBarry Smith Input Parameters: 4415c6c1daeSBarry Smith + comm - MPI communicator 4425c6c1daeSBarry Smith . name - name of file 4435c6c1daeSBarry Smith - type - type of file 4445c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4455c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4465c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4475c6c1daeSBarry Smith 4485c6c1daeSBarry Smith Output Parameter: 4495c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4505c6c1daeSBarry Smith 45182ea9b62SBarry Smith Options Database: 45282ea9b62SBarry Smith . -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 4539a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 45482ea9b62SBarry Smith 4555c6c1daeSBarry Smith Level: beginner 4565c6c1daeSBarry Smith 4575c6c1daeSBarry Smith Note: 4585c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4595c6c1daeSBarry Smith 4605c6c1daeSBarry Smith Concepts: HDF5 files 4615c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4625c6c1daeSBarry Smith 4636a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4649a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 465a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4665c6c1daeSBarry Smith @*/ 4675c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4685c6c1daeSBarry Smith { 4695c6c1daeSBarry Smith PetscErrorCode ierr; 4705c6c1daeSBarry Smith 4715c6c1daeSBarry Smith PetscFunctionBegin; 4725c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 4735c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 4745c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 4755c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 4765c6c1daeSBarry Smith PetscFunctionReturn(0); 4775c6c1daeSBarry Smith } 4785c6c1daeSBarry Smith 4795c6c1daeSBarry Smith /*@C 4805c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 4815c6c1daeSBarry Smith 4825c6c1daeSBarry Smith Not collective 4835c6c1daeSBarry Smith 4845c6c1daeSBarry Smith Input Parameter: 4855c6c1daeSBarry Smith . viewer - the PetscViewer 4865c6c1daeSBarry Smith 4875c6c1daeSBarry Smith Output Parameter: 4885c6c1daeSBarry Smith . file_id - The file id 4895c6c1daeSBarry Smith 4905c6c1daeSBarry Smith Level: intermediate 4915c6c1daeSBarry Smith 4925c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 4935c6c1daeSBarry Smith @*/ 4945c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 4955c6c1daeSBarry Smith { 4965c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4975c6c1daeSBarry Smith 4985c6c1daeSBarry Smith PetscFunctionBegin; 4995c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5005c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5015c6c1daeSBarry Smith PetscFunctionReturn(0); 5025c6c1daeSBarry Smith } 5035c6c1daeSBarry Smith 5045c6c1daeSBarry Smith /*@C 5055c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5065c6c1daeSBarry Smith 5075c6c1daeSBarry Smith Not collective 5085c6c1daeSBarry Smith 5095c6c1daeSBarry Smith Input Parameters: 5105c6c1daeSBarry Smith + viewer - the PetscViewer 5115c6c1daeSBarry Smith - name - The group name 5125c6c1daeSBarry Smith 5135c6c1daeSBarry Smith Level: intermediate 5145c6c1daeSBarry Smith 515874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5165c6c1daeSBarry Smith @*/ 5175c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5185c6c1daeSBarry Smith { 5195c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5205c6c1daeSBarry Smith GroupList *groupNode; 5215c6c1daeSBarry Smith PetscErrorCode ierr; 5225c6c1daeSBarry Smith 5235c6c1daeSBarry Smith PetscFunctionBegin; 5245c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5255c6c1daeSBarry Smith PetscValidCharPointer(name,2); 52695dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5275c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 528a297a907SKarl Rupp 5295c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5305c6c1daeSBarry Smith hdf5->groups = groupNode; 5315c6c1daeSBarry Smith PetscFunctionReturn(0); 5325c6c1daeSBarry Smith } 5335c6c1daeSBarry Smith 5343ef9c667SSatish Balay /*@ 5355c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5365c6c1daeSBarry Smith 5375c6c1daeSBarry Smith Not collective 5385c6c1daeSBarry Smith 5395c6c1daeSBarry Smith Input Parameter: 5405c6c1daeSBarry Smith . viewer - the PetscViewer 5415c6c1daeSBarry Smith 5425c6c1daeSBarry Smith Level: intermediate 5435c6c1daeSBarry Smith 544874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5455c6c1daeSBarry Smith @*/ 5465c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5475c6c1daeSBarry Smith { 5485c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5495c6c1daeSBarry Smith GroupList *groupNode; 5505c6c1daeSBarry Smith PetscErrorCode ierr; 5515c6c1daeSBarry Smith 5525c6c1daeSBarry Smith PetscFunctionBegin; 5535c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 55482f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5555c6c1daeSBarry Smith groupNode = hdf5->groups; 5565c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5575c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5585c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5595c6c1daeSBarry Smith PetscFunctionReturn(0); 5605c6c1daeSBarry Smith } 5615c6c1daeSBarry Smith 5625c6c1daeSBarry Smith /*@C 563874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 564874270d9SVaclav Hapla If none has been assigned, returns NULL. 5655c6c1daeSBarry Smith 5665c6c1daeSBarry Smith Not collective 5675c6c1daeSBarry Smith 5685c6c1daeSBarry Smith Input Parameter: 5695c6c1daeSBarry Smith . viewer - the PetscViewer 5705c6c1daeSBarry Smith 5715c6c1daeSBarry Smith Output Parameter: 5725c6c1daeSBarry Smith . name - The group name 5735c6c1daeSBarry Smith 5745c6c1daeSBarry Smith Level: intermediate 5755c6c1daeSBarry Smith 576874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 5775c6c1daeSBarry Smith @*/ 5785c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 5795c6c1daeSBarry Smith { 5805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 5815c6c1daeSBarry Smith 5825c6c1daeSBarry Smith PetscFunctionBegin; 5835c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 584c959eef4SJed Brown PetscValidPointer(name,2); 585a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 5860298fd71SBarry Smith else *name = NULL; 5875c6c1daeSBarry Smith PetscFunctionReturn(0); 5885c6c1daeSBarry Smith } 5895c6c1daeSBarry Smith 5908c557b5aSVaclav Hapla /*@ 591874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 592874270d9SVaclav Hapla and return this group's ID and file ID. 593874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 594874270d9SVaclav Hapla 595874270d9SVaclav Hapla Not collective 596874270d9SVaclav Hapla 597874270d9SVaclav Hapla Input Parameter: 598874270d9SVaclav Hapla . viewer - the PetscViewer 599874270d9SVaclav Hapla 600874270d9SVaclav Hapla Output Parameter: 601874270d9SVaclav Hapla + fileId - The HDF5 file ID 602874270d9SVaclav Hapla - groupId - The HDF5 group ID 603874270d9SVaclav Hapla 604874270d9SVaclav Hapla Level: intermediate 605874270d9SVaclav Hapla 606874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 607874270d9SVaclav Hapla @*/ 60854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 60954dbf706SVaclav Hapla { 61054dbf706SVaclav Hapla hid_t file_id, group; 61154dbf706SVaclav Hapla htri_t found; 61254dbf706SVaclav Hapla const char *groupName = NULL; 61354dbf706SVaclav Hapla PetscErrorCode ierr; 61454dbf706SVaclav Hapla 61554dbf706SVaclav Hapla PetscFunctionBegin; 61654dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 61754dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 61854dbf706SVaclav Hapla /* Open group */ 61954dbf706SVaclav Hapla if (groupName) { 62054dbf706SVaclav Hapla PetscBool root; 62154dbf706SVaclav Hapla 62254dbf706SVaclav Hapla ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr); 62354dbf706SVaclav Hapla PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT)); 62454dbf706SVaclav Hapla if (!root && (found <= 0)) { 62554dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT)); 62654dbf706SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 62754dbf706SVaclav Hapla } 62854dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT)); 62954dbf706SVaclav Hapla } else group = file_id; 63054dbf706SVaclav Hapla 63154dbf706SVaclav Hapla *fileId = file_id; 63254dbf706SVaclav Hapla *groupId = group; 63354dbf706SVaclav Hapla PetscFunctionReturn(0); 63454dbf706SVaclav Hapla } 63554dbf706SVaclav Hapla 6363ef9c667SSatish Balay /*@ 6375c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6385c6c1daeSBarry Smith 6395c6c1daeSBarry Smith Not collective 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith Input Parameter: 6425c6c1daeSBarry Smith . viewer - the PetscViewer 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith Level: intermediate 6455c6c1daeSBarry Smith 6465c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6475c6c1daeSBarry Smith @*/ 6485c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6495c6c1daeSBarry Smith { 6505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6515c6c1daeSBarry Smith 6525c6c1daeSBarry Smith PetscFunctionBegin; 6535c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6545c6c1daeSBarry Smith ++hdf5->timestep; 6555c6c1daeSBarry Smith PetscFunctionReturn(0); 6565c6c1daeSBarry Smith } 6575c6c1daeSBarry Smith 6583ef9c667SSatish Balay /*@ 6595c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6605c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6615c6c1daeSBarry Smith 6625c6c1daeSBarry Smith Not collective 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith Input Parameters: 6655c6c1daeSBarry Smith + viewer - the PetscViewer 6665c6c1daeSBarry Smith - timestep - The timestep number 6675c6c1daeSBarry Smith 6685c6c1daeSBarry Smith Level: intermediate 6695c6c1daeSBarry Smith 6705c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6715c6c1daeSBarry Smith @*/ 6725c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6735c6c1daeSBarry Smith { 6745c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6755c6c1daeSBarry Smith 6765c6c1daeSBarry Smith PetscFunctionBegin; 6775c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6785c6c1daeSBarry Smith hdf5->timestep = timestep; 6795c6c1daeSBarry Smith PetscFunctionReturn(0); 6805c6c1daeSBarry Smith } 6815c6c1daeSBarry Smith 6823ef9c667SSatish Balay /*@ 6835c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 6845c6c1daeSBarry Smith 6855c6c1daeSBarry Smith Not collective 6865c6c1daeSBarry Smith 6875c6c1daeSBarry Smith Input Parameter: 6885c6c1daeSBarry Smith . viewer - the PetscViewer 6895c6c1daeSBarry Smith 6905c6c1daeSBarry Smith Output Parameter: 6915c6c1daeSBarry Smith . timestep - The timestep number 6925c6c1daeSBarry Smith 6935c6c1daeSBarry Smith Level: intermediate 6945c6c1daeSBarry Smith 6955c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 6965c6c1daeSBarry Smith @*/ 6975c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 6985c6c1daeSBarry Smith { 6995c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7005c6c1daeSBarry Smith 7015c6c1daeSBarry Smith PetscFunctionBegin; 7025c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7035c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7045c6c1daeSBarry Smith *timestep = hdf5->timestep; 7055c6c1daeSBarry Smith PetscFunctionReturn(0); 7065c6c1daeSBarry Smith } 7075c6c1daeSBarry Smith 70836bce6e8SMatthew G. Knepley /*@C 70936bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 71036bce6e8SMatthew G. Knepley 71136bce6e8SMatthew G. Knepley Not collective 71236bce6e8SMatthew G. Knepley 71336bce6e8SMatthew G. Knepley Input Parameter: 71436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 71536bce6e8SMatthew G. Knepley 71636bce6e8SMatthew G. Knepley Output Parameter: 71736bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 71836bce6e8SMatthew G. Knepley 71936bce6e8SMatthew G. Knepley Level: advanced 72036bce6e8SMatthew G. Knepley 72136bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 72236bce6e8SMatthew G. Knepley @*/ 72336bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 72436bce6e8SMatthew G. Knepley { 72536bce6e8SMatthew G. Knepley PetscFunctionBegin; 72636bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 72736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 72836bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 72936bce6e8SMatthew G. Knepley #else 73036bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 73136bce6e8SMatthew G. Knepley #endif 73236bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 73336bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 73436bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 73536bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 73636bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 73736bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 73836bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 73936bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7407e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 74136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 74236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 74336bce6e8SMatthew G. Knepley } 74436bce6e8SMatthew G. Knepley 74536bce6e8SMatthew G. Knepley /*@C 74636bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 74736bce6e8SMatthew G. Knepley 74836bce6e8SMatthew G. Knepley Not collective 74936bce6e8SMatthew G. Knepley 75036bce6e8SMatthew G. Knepley Input Parameter: 75136bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 75236bce6e8SMatthew G. Knepley 75336bce6e8SMatthew G. Knepley Output Parameter: 75436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 75536bce6e8SMatthew G. Knepley 75636bce6e8SMatthew G. Knepley Level: advanced 75736bce6e8SMatthew G. Knepley 75836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 75936bce6e8SMatthew G. Knepley @*/ 76036bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 76136bce6e8SMatthew G. Knepley { 76236bce6e8SMatthew G. Knepley PetscFunctionBegin; 76336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 76436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 76536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 76636bce6e8SMatthew G. Knepley #else 76736bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 76836bce6e8SMatthew G. Knepley #endif 76936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 77036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 77136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 77236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 77336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 77436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7757e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 77636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 77736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 77836bce6e8SMatthew G. Knepley } 77936bce6e8SMatthew G. Knepley 780df863907SAlex Fikl /*@C 781b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 78236bce6e8SMatthew G. Knepley 78336bce6e8SMatthew G. Knepley Input Parameters: 78436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 78536bce6e8SMatthew G. Knepley . parent - The parent name 78636bce6e8SMatthew G. Knepley . name - The attribute name 78736bce6e8SMatthew G. Knepley . datatype - The attribute type 78836bce6e8SMatthew G. Knepley - value - The attribute value 78936bce6e8SMatthew G. Knepley 79036bce6e8SMatthew G. Knepley Level: advanced 79136bce6e8SMatthew G. Knepley 792e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 79336bce6e8SMatthew G. Knepley @*/ 79436bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 79536bce6e8SMatthew G. Knepley { 79660568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 797080f144cSVaclav Hapla PetscBool has; 79836bce6e8SMatthew G. Knepley PetscErrorCode ierr; 79936bce6e8SMatthew G. Knepley 80036bce6e8SMatthew G. Knepley PetscFunctionBegin; 8015cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 802c57a1a9aSVaclav Hapla PetscValidCharPointer(parent, 2); 803c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 80436bce6e8SMatthew G. Knepley PetscValidPointer(value, 4); 805080f144cSVaclav Hapla 806bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 807b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 80836bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8097e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8107e97c476SMatthew G. Knepley size_t len; 8117e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 812729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8137e97c476SMatthew G. Knepley } 81436bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 815729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 81660568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 817080f144cSVaclav Hapla if (has) { 818080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 819080f144cSVaclav Hapla } else { 82060568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 821080f144cSVaclav Hapla } 822729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 823729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 824729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 82560568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 826729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 82736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 82836bce6e8SMatthew G. Knepley } 82936bce6e8SMatthew G. Knepley 830df863907SAlex Fikl /*@C 831b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 832ad0c61feSMatthew G. Knepley 833ad0c61feSMatthew G. Knepley Input Parameters: 834ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 835ad0c61feSMatthew G. Knepley . parent - The parent name 836ad0c61feSMatthew G. Knepley . name - The attribute name 837ad0c61feSMatthew G. Knepley - datatype - The attribute type 838ad0c61feSMatthew G. Knepley 839ad0c61feSMatthew G. Knepley Output Parameter: 840ad0c61feSMatthew G. Knepley . value - The attribute value 841ad0c61feSMatthew G. Knepley 842ad0c61feSMatthew G. Knepley Level: advanced 843ad0c61feSMatthew G. Knepley 844e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 845ad0c61feSMatthew G. Knepley @*/ 846ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 847ad0c61feSMatthew G. Knepley { 848f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 849ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 850ad0c61feSMatthew G. Knepley 851ad0c61feSMatthew G. Knepley PetscFunctionBegin; 8525cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 853c57a1a9aSVaclav Hapla PetscValidCharPointer(parent, 2); 854c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 855ad0c61feSMatthew G. Knepley PetscValidPointer(value, 4); 856ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 857ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 85860568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 85960568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 860f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 861f0b58479SMatthew G. Knepley size_t len; 862e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 863f0b58479SMatthew G. Knepley PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 864f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 865f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 866f0b58479SMatthew G. Knepley } 86770efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 868729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 869e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 870e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 871ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 872ad0c61feSMatthew G. Knepley } 873ad0c61feSMatthew G. Knepley 874*a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 875*a75c3fd4SVaclav Hapla { 876*a75c3fd4SVaclav Hapla htri_t exists; 877*a75c3fd4SVaclav Hapla hid_t group; 878*a75c3fd4SVaclav Hapla 879*a75c3fd4SVaclav Hapla PetscFunctionBegin; 880*a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 881*a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 882*a75c3fd4SVaclav Hapla if (!exists && createGroup) { 883*a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 884*a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 885*a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 886*a75c3fd4SVaclav Hapla } 887*a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 888*a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 889*a75c3fd4SVaclav Hapla } 890*a75c3fd4SVaclav Hapla 891bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 8925cdeb108SMatthew G. Knepley { 8935cdeb108SMatthew G. Knepley hid_t h5; 894*a75c3fd4SVaclav Hapla PetscBool exists; 89585688ae2SVaclav Hapla PetscInt i,n; 89685688ae2SVaclav Hapla char **hierarchy; 89785688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 8985cdeb108SMatthew G. Knepley PetscErrorCode ierr; 8995cdeb108SMatthew G. Knepley 9005cdeb108SMatthew G. Knepley PetscFunctionBegin; 9015cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 90200385fc5SVaclav Hapla PetscValidCharPointer(name, 2); 903ccd66a83SVaclav Hapla if (has) { 90456cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 9055cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 906ccd66a83SVaclav Hapla } 907ccd66a83SVaclav Hapla if (otype) { 908ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 90956cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 910ccd66a83SVaclav Hapla } 911c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 91285688ae2SVaclav Hapla 91385688ae2SVaclav Hapla /* 91485688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 91585688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 91685688ae2SVaclav Hapla 1) whether it's a valid link 91785688ae2SVaclav Hapla 2) whether this link resolves to an object 91885688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 91985688ae2SVaclav Hapla */ 92085688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 92185688ae2SVaclav Hapla if (!n) { 92285688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 923ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 924ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 92585688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 92685688ae2SVaclav Hapla PetscFunctionReturn(0); 92785688ae2SVaclav Hapla } 92885688ae2SVaclav Hapla for (i=0; i<n; i++) { 92985688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 93085688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 931*a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 932*a75c3fd4SVaclav Hapla if (!exists) break; 93385688ae2SVaclav Hapla } 93485688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 93585688ae2SVaclav Hapla 93685688ae2SVaclav Hapla /* If the object exists, get its type */ 937ccd66a83SVaclav Hapla if (exists && otype) { 9385cdeb108SMatthew G. Knepley H5O_info_t info; 9395cdeb108SMatthew G. Knepley 94004633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 94156cc0592SVaclav Hapla *otype = info.type; 9425cdeb108SMatthew G. Knepley } 943ccd66a83SVaclav Hapla if (has) *has = exists; 9445cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 9455cdeb108SMatthew G. Knepley } 9465cdeb108SMatthew G. Knepley 947ecce1506SVaclav Hapla /*@ 9480a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 9490a9f38efSVaclav Hapla 9500a9f38efSVaclav Hapla Input Parameters: 9510a9f38efSVaclav Hapla . viewer - The HDF5 viewer 9520a9f38efSVaclav Hapla 9530a9f38efSVaclav Hapla Output Parameter: 9540a9f38efSVaclav Hapla . has - Flag for group existence 9550a9f38efSVaclav Hapla 9560a9f38efSVaclav Hapla Notes: 9570a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 9580a9f38efSVaclav Hapla 9590a9f38efSVaclav Hapla Level: advanced 9600a9f38efSVaclav Hapla 9610a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 9620a9f38efSVaclav Hapla @*/ 9630a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 9640a9f38efSVaclav Hapla { 9650a9f38efSVaclav Hapla H5O_type_t type; 9660a9f38efSVaclav Hapla const char *name; 9670a9f38efSVaclav Hapla PetscErrorCode ierr; 9680a9f38efSVaclav Hapla 9690a9f38efSVaclav Hapla PetscFunctionBegin; 9700a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9710a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 9720a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 9730a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 9740a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 9750a9f38efSVaclav Hapla PetscFunctionReturn(0); 9760a9f38efSVaclav Hapla } 9770a9f38efSVaclav Hapla 9780a9f38efSVaclav Hapla /*@ 979ecce1506SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file 980ecce1506SVaclav Hapla 981ecce1506SVaclav Hapla Input Parameters: 982ecce1506SVaclav Hapla + viewer - The HDF5 viewer 983ecce1506SVaclav Hapla - obj - The named object 984ecce1506SVaclav Hapla 985ecce1506SVaclav Hapla Output Parameter: 986ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 987ecce1506SVaclav Hapla 988ecce1506SVaclav Hapla Level: advanced 989ecce1506SVaclav Hapla 990ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute() 991ecce1506SVaclav Hapla @*/ 992ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 993ecce1506SVaclav Hapla { 99456cc0592SVaclav Hapla H5O_type_t type; 995ecce1506SVaclav Hapla PetscErrorCode ierr; 996ecce1506SVaclav Hapla 997ecce1506SVaclav Hapla PetscFunctionBegin; 998c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 999c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1000c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1001ecce1506SVaclav Hapla *has = PETSC_FALSE; 100256cc0592SVaclav Hapla if (obj->name) { 1003bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, obj->name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 100456cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 100556cc0592SVaclav Hapla } 1006ecce1506SVaclav Hapla PetscFunctionReturn(0); 1007ecce1506SVaclav Hapla } 1008ecce1506SVaclav Hapla 1009df863907SAlex Fikl /*@C 1010ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1011e7bdbf8eSMatthew G. Knepley 1012e7bdbf8eSMatthew G. Knepley Input Parameters: 1013e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 1014e7bdbf8eSMatthew G. Knepley . parent - The parent name 1015e7bdbf8eSMatthew G. Knepley - name - The attribute name 1016e7bdbf8eSMatthew G. Knepley 1017e7bdbf8eSMatthew G. Knepley Output Parameter: 1018e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1019e7bdbf8eSMatthew G. Knepley 1020e7bdbf8eSMatthew G. Knepley Level: advanced 1021e7bdbf8eSMatthew G. Knepley 1022ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasObject() 1023e7bdbf8eSMatthew G. Knepley @*/ 1024e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1025e7bdbf8eSMatthew G. Knepley { 1026e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1027e7bdbf8eSMatthew G. Knepley 1028e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 10295cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1030c57a1a9aSVaclav Hapla PetscValidCharPointer(parent,2); 1031c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1032c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 1033e7bdbf8eSMatthew G. Knepley *has = PETSC_FALSE; 1034bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 10352f936e54SVaclav Hapla if (!*has) PetscFunctionReturn(0); 103606db490cSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr); 103706db490cSVaclav Hapla PetscFunctionReturn(0); 103806db490cSVaclav Hapla } 103906db490cSVaclav Hapla 104006db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 104106db490cSVaclav Hapla { 104206db490cSVaclav Hapla hid_t h5; 104306db490cSVaclav Hapla htri_t hhas; 104406db490cSVaclav Hapla PetscErrorCode ierr; 104506db490cSVaclav Hapla 104606db490cSVaclav Hapla PetscFunctionBegin; 104706db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 10482f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 10495cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1050e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1051e7bdbf8eSMatthew G. Knepley } 1052e7bdbf8eSMatthew G. Knepley 1053eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1054a5e1feadSVaclav Hapla { 1055fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1056fbd37863SVaclav Hapla const char *groupname=NULL; 105769a06e7bSVaclav Hapla char vecgroup[PETSC_MAX_PATH_LEN]; 1058a5e1feadSVaclav Hapla PetscErrorCode ierr; 1059a5e1feadSVaclav Hapla 1060a5e1feadSVaclav Hapla PetscFunctionBegin; 106169a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 106269a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 106369a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 106469a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 10659568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 106669a06e7bSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr); 1067b44ade6dSVaclav Hapla ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr); 10689568d583SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr); 106909dabeb0SVaclav Hapla /* MATLAB stores column vectors horizontally */ 107009dabeb0SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1071b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 1072b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1073b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1074b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1075b8ef406cSVaclav Hapla #endif 1076b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 107769a06e7bSVaclav Hapla *ctx = h; 107869a06e7bSVaclav Hapla PetscFunctionReturn(0); 107969a06e7bSVaclav Hapla } 108069a06e7bSVaclav Hapla 1081eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 108269a06e7bSVaclav Hapla { 108369a06e7bSVaclav Hapla HDF5ReadCtx h; 108469a06e7bSVaclav Hapla PetscErrorCode ierr; 108569a06e7bSVaclav Hapla 108669a06e7bSVaclav Hapla PetscFunctionBegin; 108769a06e7bSVaclav Hapla h = *ctx; 1088b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 108969a06e7bSVaclav Hapla if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group)); 109069a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 109169a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 109269a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 109369a06e7bSVaclav Hapla PetscFunctionReturn(0); 109469a06e7bSVaclav Hapla } 109569a06e7bSVaclav Hapla 1096eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 109769a06e7bSVaclav Hapla { 109869a06e7bSVaclav Hapla int rdim, dim; 109969a06e7bSVaclav Hapla hsize_t dims[4]; 110009dabeb0SVaclav Hapla PetscInt bsInd, lenInd, bs, len, N; 11018374c777SVaclav Hapla PetscLayout map; 11028374c777SVaclav Hapla PetscErrorCode ierr; 110369a06e7bSVaclav Hapla 110469a06e7bSVaclav Hapla PetscFunctionBegin; 11058374c777SVaclav Hapla if (!(*map_)) { 11068374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 11078374c777SVaclav Hapla } 11088374c777SVaclav Hapla map = *map_; 11099787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1110a5e1feadSVaclav Hapla dim = 0; 11119568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1112a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 11139568d583SVaclav Hapla if (ctx->complexVal) ++dim; 11149787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 111569a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 11169787f08bSVaclav Hapla /* calculate expected dimension indices */ 11179787f08bSVaclav Hapla lenInd = 0; 11189568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 11199787f08bSVaclav Hapla bsInd = lenInd + 1; 11209568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 11219787f08bSVaclav Hapla if (rdim == dim) { 112245e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 11239787f08bSVaclav Hapla } else if (rdim == dim+1) { 112445e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 11259568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1126a5e1feadSVaclav Hapla } else { 11279787f08bSVaclav Hapla SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim); 1128a5e1feadSVaclav Hapla } 112909dabeb0SVaclav Hapla len = dims[lenInd]; 113009dabeb0SVaclav Hapla if (ctx->horizontal) { 113109dabeb0SVaclav Hapla 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."); 113209dabeb0SVaclav Hapla len = bs; 113309dabeb0SVaclav Hapla bs = 1; 113409dabeb0SVaclav Hapla } 113509dabeb0SVaclav Hapla N = (PetscInt) len*bs; 11368374c777SVaclav Hapla 11378374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 11388374c777SVaclav Hapla if (map->bs < 0) { 113945e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 114045e21b7fSVaclav Hapla } 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); 11418374c777SVaclav Hapla if (map->N < 0) { 114245e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 114345e21b7fSVaclav Hapla } 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); 114469a06e7bSVaclav Hapla PetscFunctionReturn(0); 114569a06e7bSVaclav Hapla } 114669a06e7bSVaclav Hapla 1147eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 114816f6402dSVaclav Hapla { 114916f6402dSVaclav Hapla hsize_t count[4], offset[4]; 115016f6402dSVaclav Hapla int dim; 115116f6402dSVaclav Hapla PetscInt bs, n, low; 115216f6402dSVaclav Hapla PetscErrorCode ierr; 115316f6402dSVaclav Hapla 115416f6402dSVaclav Hapla PetscFunctionBegin; 115516f6402dSVaclav Hapla /* Compute local size and ownership range */ 115616f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 115716f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 115816f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 115916f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 116016f6402dSVaclav Hapla 116116f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 116216f6402dSVaclav Hapla dim = 0; 116316f6402dSVaclav Hapla if (ctx->timestep >= 0) { 116416f6402dSVaclav Hapla count[dim] = 1; 116516f6402dSVaclav Hapla offset[dim] = ctx->timestep; 116616f6402dSVaclav Hapla ++dim; 116716f6402dSVaclav Hapla } 116809dabeb0SVaclav Hapla if (ctx->horizontal) { 116909dabeb0SVaclav Hapla count[dim] = 1; 117009dabeb0SVaclav Hapla offset[dim] = 0; 117109dabeb0SVaclav Hapla ++dim; 117209dabeb0SVaclav Hapla } 117316f6402dSVaclav Hapla { 117416f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 117516f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 117616f6402dSVaclav Hapla ++dim; 117716f6402dSVaclav Hapla } 117816f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 117909dabeb0SVaclav Hapla if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 118016f6402dSVaclav Hapla count[dim] = bs; 118116f6402dSVaclav Hapla offset[dim] = 0; 118216f6402dSVaclav Hapla ++dim; 118316f6402dSVaclav Hapla } 118416f6402dSVaclav Hapla if (ctx->complexVal) { 118516f6402dSVaclav Hapla count[dim] = 2; 118616f6402dSVaclav Hapla offset[dim] = 0; 118716f6402dSVaclav Hapla ++dim; 118816f6402dSVaclav Hapla } 118916f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 119016f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 119116f6402dSVaclav Hapla PetscFunctionReturn(0); 119216f6402dSVaclav Hapla } 119316f6402dSVaclav Hapla 1194eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1195ef2d82ceSVaclav Hapla { 1196ef2d82ceSVaclav Hapla PetscFunctionBegin; 1197ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1198ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1199ef2d82ceSVaclav Hapla } 1200ef2d82ceSVaclav Hapla 1201dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1202a25c73c6SVaclav Hapla { 1203fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1204fbd37863SVaclav Hapla hid_t memspace=0; 1205a25c73c6SVaclav Hapla size_t unitsize; 1206a25c73c6SVaclav Hapla void *arr; 1207a25c73c6SVaclav Hapla PetscErrorCode ierr; 1208a25c73c6SVaclav Hapla 1209a25c73c6SVaclav Hapla PetscFunctionBegin; 1210eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1211eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1212a25c73c6SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1213eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 12144fc17bcdSVaclav Hapla 12155a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 12165a89fdf4SVaclav Hapla if (!h->complexVal) { 1217c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1218c76551afSVaclav Hapla 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."); 12195a89fdf4SVaclav Hapla } 12205a89fdf4SVaclav Hapla #else 12215a89fdf4SVaclav Hapla 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."); 12225a89fdf4SVaclav Hapla #endif 12234fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 12244fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 1225dff35581SVaclav Hapla 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); 12264fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 12274fc17bcdSVaclav Hapla 1228eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1229a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1230eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1231a25c73c6SVaclav Hapla *newarr = arr; 1232a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1233a25c73c6SVaclav Hapla } 1234a25c73c6SVaclav Hapla 1235c1aaad9cSVaclav Hapla /*@C 1236c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1237c1aaad9cSVaclav Hapla 1238c1aaad9cSVaclav Hapla Input Parameters: 1239c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1240c1aaad9cSVaclav Hapla - name - The vector name 1241c1aaad9cSVaclav Hapla 1242c1aaad9cSVaclav Hapla Output Parameter: 1243c1aaad9cSVaclav Hapla + bs - block size 1244c1aaad9cSVaclav Hapla - N - global size 1245c1aaad9cSVaclav Hapla 1246c1aaad9cSVaclav Hapla Note: 1247c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1248c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1249c1aaad9cSVaclav Hapla 1250c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1251c1aaad9cSVaclav Hapla 1252c1aaad9cSVaclav Hapla Level: advanced 1253c1aaad9cSVaclav Hapla 1254c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1255c1aaad9cSVaclav Hapla @*/ 125669a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 125769a06e7bSVaclav Hapla { 1258fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 12598374c777SVaclav Hapla PetscLayout map=NULL; 126069a06e7bSVaclav Hapla PetscErrorCode ierr; 126169a06e7bSVaclav Hapla 126269a06e7bSVaclav Hapla PetscFunctionBegin; 1263c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1264eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1265eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1266eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 12678374c777SVaclav Hapla if (bs) *bs = map->bs; 12688374c777SVaclav Hapla if (N) *N = map->N; 12698374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1270a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1271a5e1feadSVaclav Hapla } 1272a5e1feadSVaclav Hapla 1273a75e6a4aSMatthew G. Knepley /* 1274a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1275a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1276a75e6a4aSMatthew G. Knepley */ 1277d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1278a75e6a4aSMatthew G. Knepley 1279a75e6a4aSMatthew G. Knepley /*@C 1280a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1281a75e6a4aSMatthew G. Knepley 1282a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1283a75e6a4aSMatthew G. Knepley 1284a75e6a4aSMatthew G. Knepley Input Parameter: 1285a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1286a75e6a4aSMatthew G. Knepley 1287a75e6a4aSMatthew G. Knepley Level: intermediate 1288a75e6a4aSMatthew G. Knepley 1289a75e6a4aSMatthew G. Knepley Options Database Keys: 1290a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1291a75e6a4aSMatthew G. Knepley 1292a75e6a4aSMatthew G. Knepley Environmental variables: 1293a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1294a75e6a4aSMatthew G. Knepley 1295a75e6a4aSMatthew G. Knepley Notes: 1296a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1297a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1298a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1299a75e6a4aSMatthew G. Knepley 1300a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1301a75e6a4aSMatthew G. Knepley @*/ 1302a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1303a75e6a4aSMatthew G. Knepley { 1304a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1305a75e6a4aSMatthew G. Knepley PetscBool flg; 1306a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1307a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1308a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1309a75e6a4aSMatthew G. Knepley 1310a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1311a75e6a4aSMatthew G. Knepley ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1312a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 131312801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1314a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1315a75e6a4aSMatthew G. Knepley } 131647435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1317a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1318a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1319a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1320a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1321a75e6a4aSMatthew G. Knepley if (!flg) { 1322a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1323a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1324a75e6a4aSMatthew G. Knepley } 1325a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1326a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1327a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1328a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 132947435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1330a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1331a75e6a4aSMatthew G. Knepley } 1332a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1333a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1334a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1335a75e6a4aSMatthew G. Knepley } 1336