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 75c6c1daeSBarry Smith typedef struct GroupList { 85c6c1daeSBarry Smith const char *name; 95c6c1daeSBarry Smith struct GroupList *next; 105c6c1daeSBarry Smith } GroupList; 115c6c1daeSBarry Smith 125c6c1daeSBarry Smith typedef struct { 135c6c1daeSBarry Smith char *filename; 145c6c1daeSBarry Smith PetscFileMode btype; 155c6c1daeSBarry Smith hid_t file_id; 165c6c1daeSBarry Smith PetscInt timestep; 175c6c1daeSBarry Smith GroupList *groups; 1882ea9b62SBarry Smith PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 199a0502c6SHåkon Strandenes PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 20058bd781SVaclav Hapla char *mataij_iname; 21058bd781SVaclav Hapla char *mataij_jname; 22058bd781SVaclav Hapla char *mataij_aname; 23058bd781SVaclav Hapla char *mataij_cname; 245c6c1daeSBarry Smith } PetscViewer_HDF5; 255c6c1daeSBarry Smith 26eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx { 27eb5a92b4SVaclav Hapla hid_t file, group, dataset, dataspace, plist; 28eb5a92b4SVaclav Hapla PetscInt timestep; 2909dabeb0SVaclav Hapla PetscBool complexVal, dim2, horizontal; 30eb5a92b4SVaclav Hapla }; 31eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 32eb5a92b4SVaclav Hapla 334416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 3482ea9b62SBarry Smith { 3582ea9b62SBarry Smith PetscErrorCode ierr; 3682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 3782ea9b62SBarry Smith 3882ea9b62SBarry Smith PetscFunctionBegin; 3982ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 4082ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 419a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 4282ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4382ea9b62SBarry Smith PetscFunctionReturn(0); 4482ea9b62SBarry Smith } 4582ea9b62SBarry Smith 465c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 475c6c1daeSBarry Smith { 485c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 495c6c1daeSBarry Smith PetscErrorCode ierr; 505c6c1daeSBarry Smith 515c6c1daeSBarry Smith PetscFunctionBegin; 525c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 53729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 545c6c1daeSBarry Smith PetscFunctionReturn(0); 555c6c1daeSBarry Smith } 565c6c1daeSBarry Smith 575c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 585c6c1daeSBarry Smith { 595c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 605c6c1daeSBarry Smith PetscErrorCode ierr; 615c6c1daeSBarry Smith 625c6c1daeSBarry Smith PetscFunctionBegin; 635c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 645c6c1daeSBarry Smith while (hdf5->groups) { 655c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 665c6c1daeSBarry Smith 675c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 685c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 695c6c1daeSBarry Smith hdf5->groups = tmp; 705c6c1daeSBarry Smith } 7161912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 7261912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 7361912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 7461912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 755c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 760b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 77d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 780b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 79058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 80058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 81058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 82058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 835c6c1daeSBarry Smith PetscFunctionReturn(0); 845c6c1daeSBarry Smith } 855c6c1daeSBarry Smith 865c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 875c6c1daeSBarry Smith { 885c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 895c6c1daeSBarry Smith 905c6c1daeSBarry Smith PetscFunctionBegin; 915c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 925c6c1daeSBarry Smith hdf5->btype = type; 935c6c1daeSBarry Smith PetscFunctionReturn(0); 945c6c1daeSBarry Smith } 955c6c1daeSBarry Smith 9682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 9782ea9b62SBarry Smith { 9882ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9982ea9b62SBarry Smith 10082ea9b62SBarry Smith PetscFunctionBegin; 10182ea9b62SBarry Smith hdf5->basedimension2 = flg; 10282ea9b62SBarry Smith PetscFunctionReturn(0); 10382ea9b62SBarry Smith } 10482ea9b62SBarry Smith 105df863907SAlex Fikl /*@ 10682ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 10782ea9b62SBarry Smith dimension of 2. 10882ea9b62SBarry Smith 10982ea9b62SBarry Smith Logically Collective on PetscViewer 11082ea9b62SBarry Smith 11182ea9b62SBarry Smith Input Parameters: 11282ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 11382ea9b62SBarry 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 11482ea9b62SBarry Smith 11582ea9b62SBarry Smith Options Database: 11682ea9b62SBarry 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 11782ea9b62SBarry Smith 11882ea9b62SBarry Smith 11995452b02SPatrick Sanan Notes: 12095452b02SPatrick 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 12182ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 12282ea9b62SBarry Smith 12382ea9b62SBarry Smith Level: intermediate 12482ea9b62SBarry Smith 12582ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 12682ea9b62SBarry Smith 12782ea9b62SBarry Smith @*/ 12882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 12982ea9b62SBarry Smith { 13082ea9b62SBarry Smith PetscErrorCode ierr; 13182ea9b62SBarry Smith 13282ea9b62SBarry Smith PetscFunctionBegin; 13382ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 13482ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 13582ea9b62SBarry Smith PetscFunctionReturn(0); 13682ea9b62SBarry Smith } 13782ea9b62SBarry Smith 138df863907SAlex Fikl /*@ 13982ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14082ea9b62SBarry Smith dimension of 2. 14182ea9b62SBarry Smith 14282ea9b62SBarry Smith Logically Collective on PetscViewer 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith Input Parameter: 14582ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith Output Parameter: 14882ea9b62SBarry 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 14982ea9b62SBarry Smith 15095452b02SPatrick Sanan Notes: 15195452b02SPatrick 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 15282ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15382ea9b62SBarry Smith 15482ea9b62SBarry Smith Level: intermediate 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15782ea9b62SBarry Smith 15882ea9b62SBarry Smith @*/ 15982ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 16082ea9b62SBarry Smith { 16182ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 16282ea9b62SBarry Smith 16382ea9b62SBarry Smith PetscFunctionBegin; 16482ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16582ea9b62SBarry Smith *flg = hdf5->basedimension2; 16682ea9b62SBarry Smith PetscFunctionReturn(0); 16782ea9b62SBarry Smith } 16882ea9b62SBarry Smith 1699a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1709a0502c6SHåkon Strandenes { 1719a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1729a0502c6SHåkon Strandenes 1739a0502c6SHåkon Strandenes PetscFunctionBegin; 1749a0502c6SHåkon Strandenes hdf5->spoutput = flg; 1759a0502c6SHåkon Strandenes PetscFunctionReturn(0); 1769a0502c6SHåkon Strandenes } 1779a0502c6SHåkon Strandenes 178df863907SAlex Fikl /*@ 1799a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 1809a0502c6SHåkon Strandenes compiled with double precision PetscReal. 1819a0502c6SHåkon Strandenes 1829a0502c6SHåkon Strandenes Logically Collective on PetscViewer 1839a0502c6SHåkon Strandenes 1849a0502c6SHåkon Strandenes Input Parameters: 1859a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 1869a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 1879a0502c6SHåkon Strandenes 1889a0502c6SHåkon Strandenes Options Database: 1899a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 1909a0502c6SHåkon Strandenes 1919a0502c6SHåkon Strandenes 19295452b02SPatrick Sanan Notes: 19395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 1949a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 1959a0502c6SHåkon Strandenes 1969a0502c6SHåkon Strandenes Level: intermediate 1979a0502c6SHåkon Strandenes 1989a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 1999a0502c6SHåkon Strandenes PetscReal 2009a0502c6SHåkon Strandenes 2019a0502c6SHåkon Strandenes @*/ 2029a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2039a0502c6SHåkon Strandenes { 2049a0502c6SHåkon Strandenes PetscErrorCode ierr; 2059a0502c6SHåkon Strandenes 2069a0502c6SHåkon Strandenes PetscFunctionBegin; 2079a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2089a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2099a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2109a0502c6SHåkon Strandenes } 2119a0502c6SHåkon Strandenes 212df863907SAlex Fikl /*@ 2139a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2149a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2159a0502c6SHåkon Strandenes 2169a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2179a0502c6SHåkon Strandenes 2189a0502c6SHåkon Strandenes Input Parameter: 2199a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2209a0502c6SHåkon Strandenes 2219a0502c6SHåkon Strandenes Output Parameter: 2229a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2239a0502c6SHåkon Strandenes 22495452b02SPatrick Sanan Notes: 22595452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2269a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2279a0502c6SHåkon Strandenes 2289a0502c6SHåkon Strandenes Level: intermediate 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2319a0502c6SHåkon Strandenes PetscReal 2329a0502c6SHåkon Strandenes 2339a0502c6SHåkon Strandenes @*/ 2349a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2359a0502c6SHåkon Strandenes { 2369a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2379a0502c6SHåkon Strandenes 2389a0502c6SHåkon Strandenes PetscFunctionBegin; 2399a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2409a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2419a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2429a0502c6SHåkon Strandenes } 2439a0502c6SHåkon Strandenes 2445c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2455c6c1daeSBarry Smith { 2465c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2475c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2485c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2495c6c1daeSBarry Smith #endif 2505c6c1daeSBarry Smith hid_t plist_id; 2515c6c1daeSBarry Smith PetscErrorCode ierr; 2525c6c1daeSBarry Smith 2535c6c1daeSBarry Smith PetscFunctionBegin; 254f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 255f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2565c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2575c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 258729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2595c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 260729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2615c6c1daeSBarry Smith #endif 2625c6c1daeSBarry Smith /* Create or open the file collectively */ 2635c6c1daeSBarry Smith switch (hdf5->btype) { 2645c6c1daeSBarry Smith case FILE_MODE_READ: 265729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 2665c6c1daeSBarry Smith break; 2675c6c1daeSBarry Smith case FILE_MODE_APPEND: 268729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 2695c6c1daeSBarry Smith break; 2705c6c1daeSBarry Smith case FILE_MODE_WRITE: 271729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 2725c6c1daeSBarry Smith break; 2735c6c1daeSBarry Smith default: 2745c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 2755c6c1daeSBarry Smith } 2765c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 277729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 2785c6c1daeSBarry Smith PetscFunctionReturn(0); 2795c6c1daeSBarry Smith } 2805c6c1daeSBarry Smith 281d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 282d1232d7fSToby Isaac { 283d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 284d1232d7fSToby Isaac 285d1232d7fSToby Isaac PetscFunctionBegin; 286d1232d7fSToby Isaac *name = vhdf5->filename; 287d1232d7fSToby Isaac PetscFunctionReturn(0); 288d1232d7fSToby Isaac } 289d1232d7fSToby Isaac 290058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 291058bd781SVaclav Hapla { 292058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 293058bd781SVaclav Hapla PetscErrorCode ierr; 294058bd781SVaclav Hapla 295058bd781SVaclav Hapla PetscFunctionBegin; 296058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 297058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 298058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 299058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 300058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 301058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 302058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 303058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 304058bd781SVaclav Hapla PetscFunctionReturn(0); 305058bd781SVaclav Hapla } 306058bd781SVaclav Hapla 307058bd781SVaclav Hapla /*@C 308058bd781SVaclav 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. 309058bd781SVaclav Hapla 310058bd781SVaclav Hapla Collective on PetscViewer 311058bd781SVaclav Hapla 312058bd781SVaclav Hapla Input Parameters: 313058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 314058bd781SVaclav 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 315058bd781SVaclav Hapla . jname - name of dataset j representing column indices 316058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 317058bd781SVaclav Hapla - cname - name of attribute stoting column count 318058bd781SVaclav Hapla 319058bd781SVaclav Hapla Level: advanced 320058bd781SVaclav Hapla 321058bd781SVaclav Hapla Notes: 322058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 323058bd781SVaclav Hapla 324058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 325058bd781SVaclav Hapla @*/ 326058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 327058bd781SVaclav Hapla { 328058bd781SVaclav Hapla PetscErrorCode ierr; 329058bd781SVaclav Hapla 330058bd781SVaclav Hapla PetscFunctionBegin; 331058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 332058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 333058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 334058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 335058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 336058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 337058bd781SVaclav Hapla PetscFunctionReturn(0); 338058bd781SVaclav Hapla } 339058bd781SVaclav Hapla 340058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 341058bd781SVaclav Hapla { 342058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 343058bd781SVaclav Hapla 344058bd781SVaclav Hapla PetscFunctionBegin; 345058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 346058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 347058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 348058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 349058bd781SVaclav Hapla PetscFunctionReturn(0); 350058bd781SVaclav Hapla } 351058bd781SVaclav Hapla 352058bd781SVaclav Hapla /*@C 353058bd781SVaclav 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. 354058bd781SVaclav Hapla 355058bd781SVaclav Hapla Collective on PetscViewer 356058bd781SVaclav Hapla 357058bd781SVaclav Hapla Input Parameters: 358058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 359058bd781SVaclav Hapla 360058bd781SVaclav Hapla Output Parameters: 361058bd781SVaclav 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 362058bd781SVaclav Hapla . jname - name of dataset j representing column indices 363058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 364058bd781SVaclav Hapla - cname - name of attribute stoting column count 365058bd781SVaclav Hapla 366058bd781SVaclav Hapla Level: advanced 367058bd781SVaclav Hapla 368058bd781SVaclav Hapla Notes: 369058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 370058bd781SVaclav Hapla 371058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 372058bd781SVaclav Hapla @*/ 373058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 374058bd781SVaclav Hapla { 375058bd781SVaclav Hapla PetscErrorCode ierr; 376058bd781SVaclav Hapla 377058bd781SVaclav Hapla PetscFunctionBegin; 378058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 379058bd781SVaclav Hapla PetscValidPointer(iname,2); 380058bd781SVaclav Hapla PetscValidPointer(jname,3); 381058bd781SVaclav Hapla PetscValidPointer(aname,4); 382058bd781SVaclav Hapla PetscValidPointer(cname,5); 383058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 384058bd781SVaclav Hapla PetscFunctionReturn(0); 385058bd781SVaclav Hapla } 386058bd781SVaclav Hapla 3878556b5ebSBarry Smith /*MC 3888556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 3898556b5ebSBarry Smith 3908556b5ebSBarry Smith 3918556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 3928556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 3938556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 3948556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 3958556b5ebSBarry Smith 3961b266c99SBarry Smith Level: beginner 3978556b5ebSBarry Smith M*/ 398d1232d7fSToby Isaac 3998cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4005c6c1daeSBarry Smith { 4015c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4025c6c1daeSBarry Smith PetscErrorCode ierr; 4035c6c1daeSBarry Smith 4045c6c1daeSBarry Smith PetscFunctionBegin; 405b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4065c6c1daeSBarry Smith 4075c6c1daeSBarry Smith v->data = (void*) hdf5; 4085c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 40982ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 4105c6c1daeSBarry Smith v->ops->flush = 0; 4115c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4125c6c1daeSBarry Smith hdf5->filename = 0; 4135c6c1daeSBarry Smith hdf5->timestep = -1; 4140298fd71SBarry Smith hdf5->groups = NULL; 4155c6c1daeSBarry Smith 416058bd781SVaclav Hapla /* ir and jc are deliberately swapped as MATLAB uses column-major format */ 417058bd781SVaclav Hapla ierr = PetscStrallocpy("jc", &hdf5->mataij_iname);CHKERRQ(ierr); 418058bd781SVaclav Hapla ierr = PetscStrallocpy("ir", &hdf5->mataij_jname);CHKERRQ(ierr); 419058bd781SVaclav Hapla ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr); 420058bd781SVaclav Hapla ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr); 421058bd781SVaclav Hapla 4220b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 423d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4240b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 42582ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4269a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 427058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 428058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4295c6c1daeSBarry Smith PetscFunctionReturn(0); 4305c6c1daeSBarry Smith } 4315c6c1daeSBarry Smith 4325c6c1daeSBarry Smith /*@C 4335c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4345c6c1daeSBarry Smith 4355c6c1daeSBarry Smith Collective on MPI_Comm 4365c6c1daeSBarry Smith 4375c6c1daeSBarry Smith Input Parameters: 4385c6c1daeSBarry Smith + comm - MPI communicator 4395c6c1daeSBarry Smith . name - name of file 4405c6c1daeSBarry Smith - type - type of file 4415c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4425c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4435c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4445c6c1daeSBarry Smith 4455c6c1daeSBarry Smith Output Parameter: 4465c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4475c6c1daeSBarry Smith 44882ea9b62SBarry Smith Options Database: 44982ea9b62SBarry 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 4509a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 45182ea9b62SBarry Smith 4525c6c1daeSBarry Smith Level: beginner 4535c6c1daeSBarry Smith 4545c6c1daeSBarry Smith Note: 4555c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4565c6c1daeSBarry Smith 4575c6c1daeSBarry Smith Concepts: HDF5 files 4585c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4595c6c1daeSBarry Smith 4606a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4619a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 462a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4635c6c1daeSBarry Smith @*/ 4645c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4655c6c1daeSBarry Smith { 4665c6c1daeSBarry Smith PetscErrorCode ierr; 4675c6c1daeSBarry Smith 4685c6c1daeSBarry Smith PetscFunctionBegin; 4695c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 4705c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 4715c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 4725c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 4735c6c1daeSBarry Smith PetscFunctionReturn(0); 4745c6c1daeSBarry Smith } 4755c6c1daeSBarry Smith 4765c6c1daeSBarry Smith /*@C 4775c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 4785c6c1daeSBarry Smith 4795c6c1daeSBarry Smith Not collective 4805c6c1daeSBarry Smith 4815c6c1daeSBarry Smith Input Parameter: 4825c6c1daeSBarry Smith . viewer - the PetscViewer 4835c6c1daeSBarry Smith 4845c6c1daeSBarry Smith Output Parameter: 4855c6c1daeSBarry Smith . file_id - The file id 4865c6c1daeSBarry Smith 4875c6c1daeSBarry Smith Level: intermediate 4885c6c1daeSBarry Smith 4895c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 4905c6c1daeSBarry Smith @*/ 4915c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 4925c6c1daeSBarry Smith { 4935c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4945c6c1daeSBarry Smith 4955c6c1daeSBarry Smith PetscFunctionBegin; 4965c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4975c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 4985c6c1daeSBarry Smith PetscFunctionReturn(0); 4995c6c1daeSBarry Smith } 5005c6c1daeSBarry Smith 5015c6c1daeSBarry Smith /*@C 5025c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5035c6c1daeSBarry Smith 5045c6c1daeSBarry Smith Not collective 5055c6c1daeSBarry Smith 5065c6c1daeSBarry Smith Input Parameters: 5075c6c1daeSBarry Smith + viewer - the PetscViewer 5085c6c1daeSBarry Smith - name - The group name 5095c6c1daeSBarry Smith 5105c6c1daeSBarry Smith Level: intermediate 5115c6c1daeSBarry Smith 512874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5135c6c1daeSBarry Smith @*/ 5145c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5155c6c1daeSBarry Smith { 5165c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5175c6c1daeSBarry Smith GroupList *groupNode; 5185c6c1daeSBarry Smith PetscErrorCode ierr; 5195c6c1daeSBarry Smith 5205c6c1daeSBarry Smith PetscFunctionBegin; 5215c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5225c6c1daeSBarry Smith PetscValidCharPointer(name,2); 52395dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5245c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 525a297a907SKarl Rupp 5265c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5275c6c1daeSBarry Smith hdf5->groups = groupNode; 5285c6c1daeSBarry Smith PetscFunctionReturn(0); 5295c6c1daeSBarry Smith } 5305c6c1daeSBarry Smith 5313ef9c667SSatish Balay /*@ 5325c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5335c6c1daeSBarry Smith 5345c6c1daeSBarry Smith Not collective 5355c6c1daeSBarry Smith 5365c6c1daeSBarry Smith Input Parameter: 5375c6c1daeSBarry Smith . viewer - the PetscViewer 5385c6c1daeSBarry Smith 5395c6c1daeSBarry Smith Level: intermediate 5405c6c1daeSBarry Smith 541874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5425c6c1daeSBarry Smith @*/ 5435c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5445c6c1daeSBarry Smith { 5455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5465c6c1daeSBarry Smith GroupList *groupNode; 5475c6c1daeSBarry Smith PetscErrorCode ierr; 5485c6c1daeSBarry Smith 5495c6c1daeSBarry Smith PetscFunctionBegin; 5505c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 55182f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5525c6c1daeSBarry Smith groupNode = hdf5->groups; 5535c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5545c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5555c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5565c6c1daeSBarry Smith PetscFunctionReturn(0); 5575c6c1daeSBarry Smith } 5585c6c1daeSBarry Smith 5595c6c1daeSBarry Smith /*@C 560874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 561874270d9SVaclav Hapla If none has been assigned, returns NULL. 5625c6c1daeSBarry Smith 5635c6c1daeSBarry Smith Not collective 5645c6c1daeSBarry Smith 5655c6c1daeSBarry Smith Input Parameter: 5665c6c1daeSBarry Smith . viewer - the PetscViewer 5675c6c1daeSBarry Smith 5685c6c1daeSBarry Smith Output Parameter: 5695c6c1daeSBarry Smith . name - The group name 5705c6c1daeSBarry Smith 5715c6c1daeSBarry Smith Level: intermediate 5725c6c1daeSBarry Smith 573874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 5745c6c1daeSBarry Smith @*/ 5755c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 5765c6c1daeSBarry Smith { 5775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 5785c6c1daeSBarry Smith 5795c6c1daeSBarry Smith PetscFunctionBegin; 5805c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 581c959eef4SJed Brown PetscValidPointer(name,2); 582a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 5830298fd71SBarry Smith else *name = NULL; 5845c6c1daeSBarry Smith PetscFunctionReturn(0); 5855c6c1daeSBarry Smith } 5865c6c1daeSBarry Smith 5878c557b5aSVaclav Hapla /*@ 588874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 589874270d9SVaclav Hapla and return this group's ID and file ID. 590874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 591874270d9SVaclav Hapla 592874270d9SVaclav Hapla Not collective 593874270d9SVaclav Hapla 594874270d9SVaclav Hapla Input Parameter: 595874270d9SVaclav Hapla . viewer - the PetscViewer 596874270d9SVaclav Hapla 597874270d9SVaclav Hapla Output Parameter: 598874270d9SVaclav Hapla + fileId - The HDF5 file ID 599874270d9SVaclav Hapla - groupId - The HDF5 group ID 600874270d9SVaclav Hapla 601874270d9SVaclav Hapla Level: intermediate 602874270d9SVaclav Hapla 603874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 604874270d9SVaclav Hapla @*/ 60554dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 60654dbf706SVaclav Hapla { 60754dbf706SVaclav Hapla hid_t file_id, group; 60854dbf706SVaclav Hapla htri_t found; 60954dbf706SVaclav Hapla const char *groupName = NULL; 61054dbf706SVaclav Hapla PetscErrorCode ierr; 61154dbf706SVaclav Hapla 61254dbf706SVaclav Hapla PetscFunctionBegin; 61354dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 61454dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 61554dbf706SVaclav Hapla /* Open group */ 61654dbf706SVaclav Hapla if (groupName) { 61754dbf706SVaclav Hapla PetscBool root; 61854dbf706SVaclav Hapla 61954dbf706SVaclav Hapla ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr); 62054dbf706SVaclav Hapla PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT)); 62154dbf706SVaclav Hapla if (!root && (found <= 0)) { 62254dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT)); 62354dbf706SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 62454dbf706SVaclav Hapla } 62554dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT)); 62654dbf706SVaclav Hapla } else group = file_id; 62754dbf706SVaclav Hapla 62854dbf706SVaclav Hapla *fileId = file_id; 62954dbf706SVaclav Hapla *groupId = group; 63054dbf706SVaclav Hapla PetscFunctionReturn(0); 63154dbf706SVaclav Hapla } 63254dbf706SVaclav Hapla 6333ef9c667SSatish Balay /*@ 6345c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6355c6c1daeSBarry Smith 6365c6c1daeSBarry Smith Not collective 6375c6c1daeSBarry Smith 6385c6c1daeSBarry Smith Input Parameter: 6395c6c1daeSBarry Smith . viewer - the PetscViewer 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith Level: intermediate 6425c6c1daeSBarry Smith 6435c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6445c6c1daeSBarry Smith @*/ 6455c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6465c6c1daeSBarry Smith { 6475c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6485c6c1daeSBarry Smith 6495c6c1daeSBarry Smith PetscFunctionBegin; 6505c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6515c6c1daeSBarry Smith ++hdf5->timestep; 6525c6c1daeSBarry Smith PetscFunctionReturn(0); 6535c6c1daeSBarry Smith } 6545c6c1daeSBarry Smith 6553ef9c667SSatish Balay /*@ 6565c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6575c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith Not collective 6605c6c1daeSBarry Smith 6615c6c1daeSBarry Smith Input Parameters: 6625c6c1daeSBarry Smith + viewer - the PetscViewer 6635c6c1daeSBarry Smith - timestep - The timestep number 6645c6c1daeSBarry Smith 6655c6c1daeSBarry Smith Level: intermediate 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6685c6c1daeSBarry Smith @*/ 6695c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6705c6c1daeSBarry Smith { 6715c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6725c6c1daeSBarry Smith 6735c6c1daeSBarry Smith PetscFunctionBegin; 6745c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6755c6c1daeSBarry Smith hdf5->timestep = timestep; 6765c6c1daeSBarry Smith PetscFunctionReturn(0); 6775c6c1daeSBarry Smith } 6785c6c1daeSBarry Smith 6793ef9c667SSatish Balay /*@ 6805c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith Not collective 6835c6c1daeSBarry Smith 6845c6c1daeSBarry Smith Input Parameter: 6855c6c1daeSBarry Smith . viewer - the PetscViewer 6865c6c1daeSBarry Smith 6875c6c1daeSBarry Smith Output Parameter: 6885c6c1daeSBarry Smith . timestep - The timestep number 6895c6c1daeSBarry Smith 6905c6c1daeSBarry Smith Level: intermediate 6915c6c1daeSBarry Smith 6925c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 6935c6c1daeSBarry Smith @*/ 6945c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 6955c6c1daeSBarry Smith { 6965c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6975c6c1daeSBarry Smith 6985c6c1daeSBarry Smith PetscFunctionBegin; 6995c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7005c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7015c6c1daeSBarry Smith *timestep = hdf5->timestep; 7025c6c1daeSBarry Smith PetscFunctionReturn(0); 7035c6c1daeSBarry Smith } 7045c6c1daeSBarry Smith 70536bce6e8SMatthew G. Knepley /*@C 70636bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 70736bce6e8SMatthew G. Knepley 70836bce6e8SMatthew G. Knepley Not collective 70936bce6e8SMatthew G. Knepley 71036bce6e8SMatthew G. Knepley Input Parameter: 71136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 71236bce6e8SMatthew G. Knepley 71336bce6e8SMatthew G. Knepley Output Parameter: 71436bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 71536bce6e8SMatthew G. Knepley 71636bce6e8SMatthew G. Knepley Level: advanced 71736bce6e8SMatthew G. Knepley 71836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 71936bce6e8SMatthew G. Knepley @*/ 72036bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 72136bce6e8SMatthew G. Knepley { 72236bce6e8SMatthew G. Knepley PetscFunctionBegin; 72336bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 72436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 72536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 72636bce6e8SMatthew G. Knepley #else 72736bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 72836bce6e8SMatthew G. Knepley #endif 72936bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 73036bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 73136bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 73236bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 73336bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 73436bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 73536bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 73636bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7377e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 73836bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 73936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 74036bce6e8SMatthew G. Knepley } 74136bce6e8SMatthew G. Knepley 74236bce6e8SMatthew G. Knepley /*@C 74336bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 74436bce6e8SMatthew G. Knepley 74536bce6e8SMatthew G. Knepley Not collective 74636bce6e8SMatthew G. Knepley 74736bce6e8SMatthew G. Knepley Input Parameter: 74836bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 74936bce6e8SMatthew G. Knepley 75036bce6e8SMatthew G. Knepley Output Parameter: 75136bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 75236bce6e8SMatthew G. Knepley 75336bce6e8SMatthew G. Knepley Level: advanced 75436bce6e8SMatthew G. Knepley 75536bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 75636bce6e8SMatthew G. Knepley @*/ 75736bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 75836bce6e8SMatthew G. Knepley { 75936bce6e8SMatthew G. Knepley PetscFunctionBegin; 76036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 76136bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 76236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 76336bce6e8SMatthew G. Knepley #else 76436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 76536bce6e8SMatthew G. Knepley #endif 76636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 76736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 76836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 76936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 77036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 77136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7727e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 77336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 77436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 77536bce6e8SMatthew G. Knepley } 77636bce6e8SMatthew G. Knepley 777df863907SAlex Fikl /*@C 778b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 77936bce6e8SMatthew G. Knepley 78036bce6e8SMatthew G. Knepley Input Parameters: 78136bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 78236bce6e8SMatthew G. Knepley . parent - The parent name 78336bce6e8SMatthew G. Knepley . name - The attribute name 78436bce6e8SMatthew G. Knepley . datatype - The attribute type 78536bce6e8SMatthew G. Knepley - value - The attribute value 78636bce6e8SMatthew G. Knepley 78736bce6e8SMatthew G. Knepley Level: advanced 78836bce6e8SMatthew G. Knepley 789e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 79036bce6e8SMatthew G. Knepley @*/ 79136bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 79236bce6e8SMatthew G. Knepley { 79360568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 794080f144cSVaclav Hapla PetscBool has; 79536bce6e8SMatthew G. Knepley PetscErrorCode ierr; 79636bce6e8SMatthew G. Knepley 79736bce6e8SMatthew G. Knepley PetscFunctionBegin; 7985cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 79936bce6e8SMatthew G. Knepley PetscValidPointer(parent, 2); 80036bce6e8SMatthew G. Knepley PetscValidPointer(name, 3); 80136bce6e8SMatthew G. Knepley PetscValidPointer(value, 4); 802080f144cSVaclav Hapla 803080f144cSVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, parent, name, &has);CHKERRQ(ierr); 80436bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8057e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8067e97c476SMatthew G. Knepley size_t len; 8077e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 808729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8097e97c476SMatthew G. Knepley } 81036bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 811729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 81260568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 813080f144cSVaclav Hapla if (has) { 814080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 815080f144cSVaclav Hapla } else { 81660568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 817080f144cSVaclav Hapla } 818729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 819729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 820729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 82160568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 822729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 82336bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 82436bce6e8SMatthew G. Knepley } 82536bce6e8SMatthew G. Knepley 826df863907SAlex Fikl /*@C 827b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 828ad0c61feSMatthew G. Knepley 829ad0c61feSMatthew G. Knepley Input Parameters: 830ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 831ad0c61feSMatthew G. Knepley . parent - The parent name 832ad0c61feSMatthew G. Knepley . name - The attribute name 833ad0c61feSMatthew G. Knepley - datatype - The attribute type 834ad0c61feSMatthew G. Knepley 835ad0c61feSMatthew G. Knepley Output Parameter: 836ad0c61feSMatthew G. Knepley . value - The attribute value 837ad0c61feSMatthew G. Knepley 838ad0c61feSMatthew G. Knepley Level: advanced 839ad0c61feSMatthew G. Knepley 840e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 841ad0c61feSMatthew G. Knepley @*/ 842ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 843ad0c61feSMatthew G. Knepley { 844f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 845ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 846ad0c61feSMatthew G. Knepley 847ad0c61feSMatthew G. Knepley PetscFunctionBegin; 8485cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 849ad0c61feSMatthew G. Knepley PetscValidPointer(parent, 2); 850ad0c61feSMatthew G. Knepley PetscValidPointer(name, 3); 851ad0c61feSMatthew G. Knepley PetscValidPointer(value, 4); 852ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 853ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 85460568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 85560568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 856f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 857f0b58479SMatthew G. Knepley size_t len; 858e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 859f0b58479SMatthew G. Knepley PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 860f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 861f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 862f0b58479SMatthew G. Knepley } 86370efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 864729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 865e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 866e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 867ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 868ad0c61feSMatthew G. Knepley } 869ad0c61feSMatthew G. Knepley 870*7c528401SVaclav Hapla static PetscErrorCode PetscViewerHDF5HasObject_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 8715cdeb108SMatthew G. Knepley { 8725cdeb108SMatthew G. Knepley hid_t h5; 8731edf8cb3SVaclav Hapla htri_t exists; 87485688ae2SVaclav Hapla PetscInt i,n; 87585688ae2SVaclav Hapla char **hierarchy; 87685688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 8775cdeb108SMatthew G. Knepley PetscErrorCode ierr; 8785cdeb108SMatthew G. Knepley 8795cdeb108SMatthew G. Knepley PetscFunctionBegin; 8805cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 88100385fc5SVaclav Hapla PetscValidCharPointer(name, 2); 882ccd66a83SVaclav Hapla if (has) { 88356cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 8845cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 885ccd66a83SVaclav Hapla } 886ccd66a83SVaclav Hapla if (otype) { 887ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 88856cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 889ccd66a83SVaclav Hapla } 890c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 89185688ae2SVaclav Hapla 89285688ae2SVaclav Hapla /* 89385688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 89485688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 89585688ae2SVaclav Hapla 1) whether it's a valid link 89685688ae2SVaclav Hapla 2) whether this link resolves to an object 89785688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 89885688ae2SVaclav Hapla */ 89985688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 90085688ae2SVaclav Hapla if (!n) { 90185688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 902ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 903ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 90485688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 90585688ae2SVaclav Hapla PetscFunctionReturn(0); 90685688ae2SVaclav Hapla } 90785688ae2SVaclav Hapla for (i=0; i<n; i++) { 90885688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 90985688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 91085688ae2SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, buf, H5P_DEFAULT)); 91185688ae2SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, buf, H5P_DEFAULT)); 912*7c528401SVaclav Hapla if (!exists) { 913*7c528401SVaclav Hapla if (createGroup) { 914*7c528401SVaclav Hapla hid_t group; 915*7c528401SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, buf, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 916*7c528401SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 917*7c528401SVaclav Hapla exists = PETSC_TRUE; 918*7c528401SVaclav Hapla } else break; 919*7c528401SVaclav Hapla } 92085688ae2SVaclav Hapla } 92185688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 92285688ae2SVaclav Hapla 92385688ae2SVaclav Hapla /* If the object exists, get its type */ 924ccd66a83SVaclav Hapla if (exists && otype) { 9255cdeb108SMatthew G. Knepley H5O_info_t info; 9265cdeb108SMatthew G. Knepley 92704633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 92856cc0592SVaclav Hapla *otype = info.type; 9295cdeb108SMatthew G. Knepley } 930ccd66a83SVaclav Hapla if (has) *has = exists; 9315cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 9325cdeb108SMatthew G. Knepley } 9335cdeb108SMatthew G. Knepley 934ecce1506SVaclav Hapla /*@ 935ecce1506SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file 936ecce1506SVaclav Hapla 937ecce1506SVaclav Hapla Input Parameters: 938ecce1506SVaclav Hapla + viewer - The HDF5 viewer 939ecce1506SVaclav Hapla - obj - The named object 940ecce1506SVaclav Hapla 941ecce1506SVaclav Hapla Output Parameter: 942ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 943ecce1506SVaclav Hapla 944ecce1506SVaclav Hapla Level: advanced 945ecce1506SVaclav Hapla 946ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute() 947ecce1506SVaclav Hapla @*/ 948ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 949ecce1506SVaclav Hapla { 95056cc0592SVaclav Hapla H5O_type_t type; 951ecce1506SVaclav Hapla PetscErrorCode ierr; 952ecce1506SVaclav Hapla 953ecce1506SVaclav Hapla PetscFunctionBegin; 954ecce1506SVaclav Hapla *has = PETSC_FALSE; 95556cc0592SVaclav Hapla if (obj->name) { 956*7c528401SVaclav Hapla ierr = PetscViewerHDF5HasObject_Internal(viewer, obj->name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 95756cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 95856cc0592SVaclav Hapla } 959ecce1506SVaclav Hapla PetscFunctionReturn(0); 960ecce1506SVaclav Hapla } 961ecce1506SVaclav Hapla 962df863907SAlex Fikl /*@C 963ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 964e7bdbf8eSMatthew G. Knepley 965e7bdbf8eSMatthew G. Knepley Input Parameters: 966e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 967e7bdbf8eSMatthew G. Knepley . parent - The parent name 968e7bdbf8eSMatthew G. Knepley - name - The attribute name 969e7bdbf8eSMatthew G. Knepley 970e7bdbf8eSMatthew G. Knepley Output Parameter: 971e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 972e7bdbf8eSMatthew G. Knepley 973e7bdbf8eSMatthew G. Knepley Level: advanced 974e7bdbf8eSMatthew G. Knepley 975ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasObject() 976e7bdbf8eSMatthew G. Knepley @*/ 977e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 978e7bdbf8eSMatthew G. Knepley { 979729ab6d9SBarry Smith hid_t h5, dataset; 9805cdeb108SMatthew G. Knepley htri_t hhas; 9815cdeb108SMatthew G. Knepley PetscBool exists; 98256cc0592SVaclav Hapla H5O_type_t type; 983e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 984e7bdbf8eSMatthew G. Knepley 985e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 9865cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 987e7bdbf8eSMatthew G. Knepley PetscValidPointer(parent, 2); 988e7bdbf8eSMatthew G. Knepley PetscValidPointer(name, 3); 989e7bdbf8eSMatthew G. Knepley PetscValidPointer(has, 4); 990e7bdbf8eSMatthew G. Knepley *has = PETSC_FALSE; 991e7bdbf8eSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 992*7c528401SVaclav Hapla ierr = PetscViewerHDF5HasObject_Internal(viewer, parent, PETSC_FALSE, &exists, &type);CHKERRQ(ierr); 99356cc0592SVaclav Hapla if (!exists || type != H5O_TYPE_DATASET) PetscFunctionReturn(0); 9945eab803fSVaclav Hapla PetscStackCallHDF5Return(dataset, H5Dopen2, (h5, parent, H5P_DEFAULT)); 995a56d7b29SVaclav Hapla PetscStackCallHDF5Return(hhas, H5Aexists, (dataset, name)); 996729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 9975cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 998e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 999e7bdbf8eSMatthew G. Knepley } 1000e7bdbf8eSMatthew G. Knepley 1001eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1002a5e1feadSVaclav Hapla { 1003fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1004fbd37863SVaclav Hapla const char *groupname=NULL; 100569a06e7bSVaclav Hapla char vecgroup[PETSC_MAX_PATH_LEN]; 1006a5e1feadSVaclav Hapla PetscErrorCode ierr; 1007a5e1feadSVaclav Hapla 1008a5e1feadSVaclav Hapla PetscFunctionBegin; 100969a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 101069a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 101169a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 101269a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 10139568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 101469a06e7bSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr); 1015b44ade6dSVaclav Hapla ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr); 10169568d583SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr); 101709dabeb0SVaclav Hapla /* MATLAB stores column vectors horizontally */ 101809dabeb0SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1019b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 1020b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1021b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1022b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1023b8ef406cSVaclav Hapla #endif 1024b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 102569a06e7bSVaclav Hapla *ctx = h; 102669a06e7bSVaclav Hapla PetscFunctionReturn(0); 102769a06e7bSVaclav Hapla } 102869a06e7bSVaclav Hapla 1029eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 103069a06e7bSVaclav Hapla { 103169a06e7bSVaclav Hapla HDF5ReadCtx h; 103269a06e7bSVaclav Hapla PetscErrorCode ierr; 103369a06e7bSVaclav Hapla 103469a06e7bSVaclav Hapla PetscFunctionBegin; 103569a06e7bSVaclav Hapla h = *ctx; 1036b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 103769a06e7bSVaclav Hapla if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group)); 103869a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 103969a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 104069a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 104169a06e7bSVaclav Hapla PetscFunctionReturn(0); 104269a06e7bSVaclav Hapla } 104369a06e7bSVaclav Hapla 1044eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 104569a06e7bSVaclav Hapla { 104669a06e7bSVaclav Hapla int rdim, dim; 104769a06e7bSVaclav Hapla hsize_t dims[4]; 104809dabeb0SVaclav Hapla PetscInt bsInd, lenInd, bs, len, N; 10498374c777SVaclav Hapla PetscLayout map; 10508374c777SVaclav Hapla PetscErrorCode ierr; 105169a06e7bSVaclav Hapla 105269a06e7bSVaclav Hapla PetscFunctionBegin; 10538374c777SVaclav Hapla if (!(*map_)) { 10548374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 10558374c777SVaclav Hapla } 10568374c777SVaclav Hapla map = *map_; 10579787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1058a5e1feadSVaclav Hapla dim = 0; 10599568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1060a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 10619568d583SVaclav Hapla if (ctx->complexVal) ++dim; 10629787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 106369a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 10649787f08bSVaclav Hapla /* calculate expected dimension indices */ 10659787f08bSVaclav Hapla lenInd = 0; 10669568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 10679787f08bSVaclav Hapla bsInd = lenInd + 1; 10689568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 10699787f08bSVaclav Hapla if (rdim == dim) { 107045e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 10719787f08bSVaclav Hapla } else if (rdim == dim+1) { 107245e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 10739568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1074a5e1feadSVaclav Hapla } else { 10759787f08bSVaclav Hapla SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim); 1076a5e1feadSVaclav Hapla } 107709dabeb0SVaclav Hapla len = dims[lenInd]; 107809dabeb0SVaclav Hapla if (ctx->horizontal) { 107909dabeb0SVaclav 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."); 108009dabeb0SVaclav Hapla len = bs; 108109dabeb0SVaclav Hapla bs = 1; 108209dabeb0SVaclav Hapla } 108309dabeb0SVaclav Hapla N = (PetscInt) len*bs; 10848374c777SVaclav Hapla 10858374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 10868374c777SVaclav Hapla if (map->bs < 0) { 108745e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 108845e21b7fSVaclav 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); 10898374c777SVaclav Hapla if (map->N < 0) { 109045e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 109145e21b7fSVaclav 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); 109269a06e7bSVaclav Hapla PetscFunctionReturn(0); 109369a06e7bSVaclav Hapla } 109469a06e7bSVaclav Hapla 1095eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 109616f6402dSVaclav Hapla { 109716f6402dSVaclav Hapla hsize_t count[4], offset[4]; 109816f6402dSVaclav Hapla int dim; 109916f6402dSVaclav Hapla PetscInt bs, n, low; 110016f6402dSVaclav Hapla PetscErrorCode ierr; 110116f6402dSVaclav Hapla 110216f6402dSVaclav Hapla PetscFunctionBegin; 110316f6402dSVaclav Hapla /* Compute local size and ownership range */ 110416f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 110516f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 110616f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 110716f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 110816f6402dSVaclav Hapla 110916f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 111016f6402dSVaclav Hapla dim = 0; 111116f6402dSVaclav Hapla if (ctx->timestep >= 0) { 111216f6402dSVaclav Hapla count[dim] = 1; 111316f6402dSVaclav Hapla offset[dim] = ctx->timestep; 111416f6402dSVaclav Hapla ++dim; 111516f6402dSVaclav Hapla } 111609dabeb0SVaclav Hapla if (ctx->horizontal) { 111709dabeb0SVaclav Hapla count[dim] = 1; 111809dabeb0SVaclav Hapla offset[dim] = 0; 111909dabeb0SVaclav Hapla ++dim; 112009dabeb0SVaclav Hapla } 112116f6402dSVaclav Hapla { 112216f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 112316f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 112416f6402dSVaclav Hapla ++dim; 112516f6402dSVaclav Hapla } 112616f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 112709dabeb0SVaclav Hapla if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 112816f6402dSVaclav Hapla count[dim] = bs; 112916f6402dSVaclav Hapla offset[dim] = 0; 113016f6402dSVaclav Hapla ++dim; 113116f6402dSVaclav Hapla } 113216f6402dSVaclav Hapla if (ctx->complexVal) { 113316f6402dSVaclav Hapla count[dim] = 2; 113416f6402dSVaclav Hapla offset[dim] = 0; 113516f6402dSVaclav Hapla ++dim; 113616f6402dSVaclav Hapla } 113716f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 113816f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 113916f6402dSVaclav Hapla PetscFunctionReturn(0); 114016f6402dSVaclav Hapla } 114116f6402dSVaclav Hapla 1142eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1143ef2d82ceSVaclav Hapla { 1144ef2d82ceSVaclav Hapla PetscFunctionBegin; 1145ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1146ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1147ef2d82ceSVaclav Hapla } 1148ef2d82ceSVaclav Hapla 1149dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1150a25c73c6SVaclav Hapla { 1151fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1152fbd37863SVaclav Hapla hid_t memspace=0; 1153a25c73c6SVaclav Hapla size_t unitsize; 1154a25c73c6SVaclav Hapla void *arr; 1155a25c73c6SVaclav Hapla PetscErrorCode ierr; 1156a25c73c6SVaclav Hapla 1157a25c73c6SVaclav Hapla PetscFunctionBegin; 1158eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1159eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1160a25c73c6SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1161eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 11624fc17bcdSVaclav Hapla 11635a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 11645a89fdf4SVaclav Hapla if (!h->complexVal) { 1165c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1166c76551afSVaclav 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."); 11675a89fdf4SVaclav Hapla } 11685a89fdf4SVaclav Hapla #else 11695a89fdf4SVaclav 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."); 11705a89fdf4SVaclav Hapla #endif 11714fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 11724fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 1173dff35581SVaclav 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); 11744fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 11754fc17bcdSVaclav Hapla 1176eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1177a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1178eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1179a25c73c6SVaclav Hapla *newarr = arr; 1180a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1181a25c73c6SVaclav Hapla } 1182a25c73c6SVaclav Hapla 1183c1aaad9cSVaclav Hapla /*@C 1184c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1185c1aaad9cSVaclav Hapla 1186c1aaad9cSVaclav Hapla Input Parameters: 1187c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1188c1aaad9cSVaclav Hapla - name - The vector name 1189c1aaad9cSVaclav Hapla 1190c1aaad9cSVaclav Hapla Output Parameter: 1191c1aaad9cSVaclav Hapla + bs - block size 1192c1aaad9cSVaclav Hapla - N - global size 1193c1aaad9cSVaclav Hapla 1194c1aaad9cSVaclav Hapla Note: 1195c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1196c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1197c1aaad9cSVaclav Hapla 1198c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1199c1aaad9cSVaclav Hapla 1200c1aaad9cSVaclav Hapla Level: advanced 1201c1aaad9cSVaclav Hapla 1202c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1203c1aaad9cSVaclav Hapla @*/ 120469a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 120569a06e7bSVaclav Hapla { 1206fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 12078374c777SVaclav Hapla PetscLayout map=NULL; 120869a06e7bSVaclav Hapla PetscErrorCode ierr; 120969a06e7bSVaclav Hapla 121069a06e7bSVaclav Hapla PetscFunctionBegin; 1211c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1212eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1213eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1214eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 12158374c777SVaclav Hapla if (bs) *bs = map->bs; 12168374c777SVaclav Hapla if (N) *N = map->N; 12178374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1218a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1219a5e1feadSVaclav Hapla } 1220a5e1feadSVaclav Hapla 1221a75e6a4aSMatthew G. Knepley /* 1222a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1223a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1224a75e6a4aSMatthew G. Knepley */ 1225d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1226a75e6a4aSMatthew G. Knepley 1227a75e6a4aSMatthew G. Knepley /*@C 1228a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1229a75e6a4aSMatthew G. Knepley 1230a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1231a75e6a4aSMatthew G. Knepley 1232a75e6a4aSMatthew G. Knepley Input Parameter: 1233a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1234a75e6a4aSMatthew G. Knepley 1235a75e6a4aSMatthew G. Knepley Level: intermediate 1236a75e6a4aSMatthew G. Knepley 1237a75e6a4aSMatthew G. Knepley Options Database Keys: 1238a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1239a75e6a4aSMatthew G. Knepley 1240a75e6a4aSMatthew G. Knepley Environmental variables: 1241a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1242a75e6a4aSMatthew G. Knepley 1243a75e6a4aSMatthew G. Knepley Notes: 1244a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1245a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1246a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1247a75e6a4aSMatthew G. Knepley 1248a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1249a75e6a4aSMatthew G. Knepley @*/ 1250a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1251a75e6a4aSMatthew G. Knepley { 1252a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1253a75e6a4aSMatthew G. Knepley PetscBool flg; 1254a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1255a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1256a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1257a75e6a4aSMatthew G. Knepley 1258a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1259a75e6a4aSMatthew 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);} 1260a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 126112801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1262a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1263a75e6a4aSMatthew G. Knepley } 126447435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1265a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1266a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1267a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1268a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1269a75e6a4aSMatthew G. Knepley if (!flg) { 1270a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1271a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1272a75e6a4aSMatthew G. Knepley } 1273a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1274a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1275a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1276a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 127747435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1278a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1279a75e6a4aSMatthew G. Knepley } 1280a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1281a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1282a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1283a75e6a4aSMatthew G. Knepley } 1284