1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I "petscsys.h" I*/ 2d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 35c6c1daeSBarry Smith 45c6c1daeSBarry Smith typedef struct GroupList { 55c6c1daeSBarry Smith const char *name; 65c6c1daeSBarry Smith struct GroupList *next; 75c6c1daeSBarry Smith } GroupList; 85c6c1daeSBarry Smith 95c6c1daeSBarry Smith typedef struct { 105c6c1daeSBarry Smith char *filename; 115c6c1daeSBarry Smith PetscFileMode btype; 125c6c1daeSBarry Smith hid_t file_id; 135c6c1daeSBarry Smith PetscInt timestep; 145c6c1daeSBarry Smith GroupList *groups; 1582ea9b62SBarry Smith PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 169a0502c6SHåkon Strandenes PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 17058bd781SVaclav Hapla char *mataij_iname; 18058bd781SVaclav Hapla char *mataij_jname; 19058bd781SVaclav Hapla char *mataij_aname; 20058bd781SVaclav Hapla char *mataij_cname; 215c6c1daeSBarry Smith } PetscViewer_HDF5; 225c6c1daeSBarry Smith 23eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx { 24eb5a92b4SVaclav Hapla hid_t file, group, dataset, dataspace, plist; 25eb5a92b4SVaclav Hapla PetscInt timestep; 26eb5a92b4SVaclav Hapla PetscBool complexVal, dim2; 27eb5a92b4SVaclav Hapla }; 28eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 29eb5a92b4SVaclav Hapla 304416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 3182ea9b62SBarry Smith { 3282ea9b62SBarry Smith PetscErrorCode ierr; 3382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 3482ea9b62SBarry Smith 3582ea9b62SBarry Smith PetscFunctionBegin; 3682ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 3782ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 389a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 3982ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 4082ea9b62SBarry Smith PetscFunctionReturn(0); 4182ea9b62SBarry Smith } 4282ea9b62SBarry Smith 435c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 445c6c1daeSBarry Smith { 455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 465c6c1daeSBarry Smith PetscErrorCode ierr; 475c6c1daeSBarry Smith 485c6c1daeSBarry Smith PetscFunctionBegin; 495c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 50729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 515c6c1daeSBarry Smith PetscFunctionReturn(0); 525c6c1daeSBarry Smith } 535c6c1daeSBarry Smith 545c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 555c6c1daeSBarry Smith { 565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 575c6c1daeSBarry Smith PetscErrorCode ierr; 585c6c1daeSBarry Smith 595c6c1daeSBarry Smith PetscFunctionBegin; 605c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 615c6c1daeSBarry Smith while (hdf5->groups) { 625c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 635c6c1daeSBarry Smith 645c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 655c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 665c6c1daeSBarry Smith hdf5->groups = tmp; 675c6c1daeSBarry Smith } 6861912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 6961912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 7061912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 7161912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 725c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 730b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 74d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 750b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 76058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 77058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 78058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 79058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 805c6c1daeSBarry Smith PetscFunctionReturn(0); 815c6c1daeSBarry Smith } 825c6c1daeSBarry Smith 835c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 845c6c1daeSBarry Smith { 855c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 865c6c1daeSBarry Smith 875c6c1daeSBarry Smith PetscFunctionBegin; 885c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 895c6c1daeSBarry Smith hdf5->btype = type; 905c6c1daeSBarry Smith PetscFunctionReturn(0); 915c6c1daeSBarry Smith } 925c6c1daeSBarry Smith 9382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 9482ea9b62SBarry Smith { 9582ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9682ea9b62SBarry Smith 9782ea9b62SBarry Smith PetscFunctionBegin; 9882ea9b62SBarry Smith hdf5->basedimension2 = flg; 9982ea9b62SBarry Smith PetscFunctionReturn(0); 10082ea9b62SBarry Smith } 10182ea9b62SBarry Smith 102df863907SAlex Fikl /*@ 10382ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 10482ea9b62SBarry Smith dimension of 2. 10582ea9b62SBarry Smith 10682ea9b62SBarry Smith Logically Collective on PetscViewer 10782ea9b62SBarry Smith 10882ea9b62SBarry Smith Input Parameters: 10982ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 11082ea9b62SBarry 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 11182ea9b62SBarry Smith 11282ea9b62SBarry Smith Options Database: 11382ea9b62SBarry 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 11482ea9b62SBarry Smith 11582ea9b62SBarry Smith 11695452b02SPatrick Sanan Notes: 11795452b02SPatrick 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 11882ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 11982ea9b62SBarry Smith 12082ea9b62SBarry Smith Level: intermediate 12182ea9b62SBarry Smith 12282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 12382ea9b62SBarry Smith 12482ea9b62SBarry Smith @*/ 12582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 12682ea9b62SBarry Smith { 12782ea9b62SBarry Smith PetscErrorCode ierr; 12882ea9b62SBarry Smith 12982ea9b62SBarry Smith PetscFunctionBegin; 13082ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 13182ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 13282ea9b62SBarry Smith PetscFunctionReturn(0); 13382ea9b62SBarry Smith } 13482ea9b62SBarry Smith 135df863907SAlex Fikl /*@ 13682ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 13782ea9b62SBarry Smith dimension of 2. 13882ea9b62SBarry Smith 13982ea9b62SBarry Smith Logically Collective on PetscViewer 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith Input Parameter: 14282ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith Output Parameter: 14582ea9b62SBarry 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 14682ea9b62SBarry Smith 14795452b02SPatrick Sanan Notes: 14895452b02SPatrick 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 14982ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15082ea9b62SBarry Smith 15182ea9b62SBarry Smith Level: intermediate 15282ea9b62SBarry Smith 15382ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15482ea9b62SBarry Smith 15582ea9b62SBarry Smith @*/ 15682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 15782ea9b62SBarry Smith { 15882ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 15982ea9b62SBarry Smith 16082ea9b62SBarry Smith PetscFunctionBegin; 16182ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16282ea9b62SBarry Smith *flg = hdf5->basedimension2; 16382ea9b62SBarry Smith PetscFunctionReturn(0); 16482ea9b62SBarry Smith } 16582ea9b62SBarry Smith 1669a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1679a0502c6SHåkon Strandenes { 1689a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1699a0502c6SHåkon Strandenes 1709a0502c6SHåkon Strandenes PetscFunctionBegin; 1719a0502c6SHåkon Strandenes hdf5->spoutput = flg; 1729a0502c6SHåkon Strandenes PetscFunctionReturn(0); 1739a0502c6SHåkon Strandenes } 1749a0502c6SHåkon Strandenes 175df863907SAlex Fikl /*@ 1769a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 1779a0502c6SHåkon Strandenes compiled with double precision PetscReal. 1789a0502c6SHåkon Strandenes 1799a0502c6SHåkon Strandenes Logically Collective on PetscViewer 1809a0502c6SHåkon Strandenes 1819a0502c6SHåkon Strandenes Input Parameters: 1829a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 1839a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 1849a0502c6SHåkon Strandenes 1859a0502c6SHåkon Strandenes Options Database: 1869a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 1879a0502c6SHåkon Strandenes 1889a0502c6SHåkon Strandenes 18995452b02SPatrick Sanan Notes: 19095452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 1919a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 1929a0502c6SHåkon Strandenes 1939a0502c6SHåkon Strandenes Level: intermediate 1949a0502c6SHåkon Strandenes 1959a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 1969a0502c6SHåkon Strandenes PetscReal 1979a0502c6SHåkon Strandenes 1989a0502c6SHåkon Strandenes @*/ 1999a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2009a0502c6SHåkon Strandenes { 2019a0502c6SHåkon Strandenes PetscErrorCode ierr; 2029a0502c6SHåkon Strandenes 2039a0502c6SHåkon Strandenes PetscFunctionBegin; 2049a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2059a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2069a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2079a0502c6SHåkon Strandenes } 2089a0502c6SHåkon Strandenes 209df863907SAlex Fikl /*@ 2109a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2119a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2129a0502c6SHåkon Strandenes 2139a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2149a0502c6SHåkon Strandenes 2159a0502c6SHåkon Strandenes Input Parameter: 2169a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2179a0502c6SHåkon Strandenes 2189a0502c6SHåkon Strandenes Output Parameter: 2199a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2209a0502c6SHåkon Strandenes 22195452b02SPatrick Sanan Notes: 22295452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2239a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2249a0502c6SHåkon Strandenes 2259a0502c6SHåkon Strandenes Level: intermediate 2269a0502c6SHåkon Strandenes 2279a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2289a0502c6SHåkon Strandenes PetscReal 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes @*/ 2319a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2329a0502c6SHåkon Strandenes { 2339a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes PetscFunctionBegin; 2369a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2379a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2389a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2399a0502c6SHåkon Strandenes } 2409a0502c6SHåkon Strandenes 2415c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2425c6c1daeSBarry Smith { 2435c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2445c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2455c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2465c6c1daeSBarry Smith #endif 2475c6c1daeSBarry Smith hid_t plist_id; 2485c6c1daeSBarry Smith PetscErrorCode ierr; 2495c6c1daeSBarry Smith 2505c6c1daeSBarry Smith PetscFunctionBegin; 251f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 252f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2535c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2545c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 255729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2565c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 257729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2585c6c1daeSBarry Smith #endif 2595c6c1daeSBarry Smith /* Create or open the file collectively */ 2605c6c1daeSBarry Smith switch (hdf5->btype) { 2615c6c1daeSBarry Smith case FILE_MODE_READ: 262729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 2635c6c1daeSBarry Smith break; 2645c6c1daeSBarry Smith case FILE_MODE_APPEND: 265729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 2665c6c1daeSBarry Smith break; 2675c6c1daeSBarry Smith case FILE_MODE_WRITE: 268729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 2695c6c1daeSBarry Smith break; 2705c6c1daeSBarry Smith default: 2715c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 2725c6c1daeSBarry Smith } 2735c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 274729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 2755c6c1daeSBarry Smith PetscFunctionReturn(0); 2765c6c1daeSBarry Smith } 2775c6c1daeSBarry Smith 278d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 279d1232d7fSToby Isaac { 280d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 281d1232d7fSToby Isaac 282d1232d7fSToby Isaac PetscFunctionBegin; 283d1232d7fSToby Isaac *name = vhdf5->filename; 284d1232d7fSToby Isaac PetscFunctionReturn(0); 285d1232d7fSToby Isaac } 286d1232d7fSToby Isaac 287058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 288058bd781SVaclav Hapla { 289058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 290058bd781SVaclav Hapla PetscErrorCode ierr; 291058bd781SVaclav Hapla 292058bd781SVaclav Hapla PetscFunctionBegin; 293058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 294058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 295058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 296058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 297058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 298058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 299058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 300058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 301058bd781SVaclav Hapla PetscFunctionReturn(0); 302058bd781SVaclav Hapla } 303058bd781SVaclav Hapla 304058bd781SVaclav Hapla /*@C 305058bd781SVaclav 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. 306058bd781SVaclav Hapla 307058bd781SVaclav Hapla Collective on PetscViewer 308058bd781SVaclav Hapla 309058bd781SVaclav Hapla Input Parameters: 310058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 311058bd781SVaclav 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 312058bd781SVaclav Hapla . jname - name of dataset j representing column indices 313058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 314058bd781SVaclav Hapla - cname - name of attribute stoting column count 315058bd781SVaclav Hapla 316058bd781SVaclav Hapla Level: advanced 317058bd781SVaclav Hapla 318058bd781SVaclav Hapla Notes: 319058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 320058bd781SVaclav Hapla 321058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 322058bd781SVaclav Hapla @*/ 323058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 324058bd781SVaclav Hapla { 325058bd781SVaclav Hapla PetscErrorCode ierr; 326058bd781SVaclav Hapla 327058bd781SVaclav Hapla PetscFunctionBegin; 328058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 329058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 330058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 331058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 332058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 333058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 334058bd781SVaclav Hapla PetscFunctionReturn(0); 335058bd781SVaclav Hapla } 336058bd781SVaclav Hapla 337058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 338058bd781SVaclav Hapla { 339058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 340058bd781SVaclav Hapla 341058bd781SVaclav Hapla PetscFunctionBegin; 342058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 343058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 344058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 345058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 346058bd781SVaclav Hapla PetscFunctionReturn(0); 347058bd781SVaclav Hapla } 348058bd781SVaclav Hapla 349058bd781SVaclav Hapla /*@C 350058bd781SVaclav 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. 351058bd781SVaclav Hapla 352058bd781SVaclav Hapla Collective on PetscViewer 353058bd781SVaclav Hapla 354058bd781SVaclav Hapla Input Parameters: 355058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 356058bd781SVaclav Hapla 357058bd781SVaclav Hapla Output Parameters: 358058bd781SVaclav 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 359058bd781SVaclav Hapla . jname - name of dataset j representing column indices 360058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 361058bd781SVaclav Hapla - cname - name of attribute stoting column count 362058bd781SVaclav Hapla 363058bd781SVaclav Hapla Level: advanced 364058bd781SVaclav Hapla 365058bd781SVaclav Hapla Notes: 366058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 367058bd781SVaclav Hapla 368058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 369058bd781SVaclav Hapla @*/ 370058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 371058bd781SVaclav Hapla { 372058bd781SVaclav Hapla PetscErrorCode ierr; 373058bd781SVaclav Hapla 374058bd781SVaclav Hapla PetscFunctionBegin; 375058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 376058bd781SVaclav Hapla PetscValidPointer(iname,2); 377058bd781SVaclav Hapla PetscValidPointer(jname,3); 378058bd781SVaclav Hapla PetscValidPointer(aname,4); 379058bd781SVaclav Hapla PetscValidPointer(cname,5); 380058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 381058bd781SVaclav Hapla PetscFunctionReturn(0); 382058bd781SVaclav Hapla } 383058bd781SVaclav Hapla 3848556b5ebSBarry Smith /*MC 3858556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 3868556b5ebSBarry Smith 3878556b5ebSBarry Smith 3888556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 3898556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 3908556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 3918556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 3928556b5ebSBarry Smith 3931b266c99SBarry Smith Level: beginner 3948556b5ebSBarry Smith M*/ 395d1232d7fSToby Isaac 3968cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 3975c6c1daeSBarry Smith { 3985c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 3995c6c1daeSBarry Smith PetscErrorCode ierr; 4005c6c1daeSBarry Smith 4015c6c1daeSBarry Smith PetscFunctionBegin; 402b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4035c6c1daeSBarry Smith 4045c6c1daeSBarry Smith v->data = (void*) hdf5; 4055c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 40682ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 4075c6c1daeSBarry Smith v->ops->flush = 0; 4085c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4095c6c1daeSBarry Smith hdf5->filename = 0; 4105c6c1daeSBarry Smith hdf5->timestep = -1; 4110298fd71SBarry Smith hdf5->groups = NULL; 4125c6c1daeSBarry Smith 413058bd781SVaclav Hapla /* ir and jc are deliberately swapped as MATLAB uses column-major format */ 414058bd781SVaclav Hapla ierr = PetscStrallocpy("jc", &hdf5->mataij_iname);CHKERRQ(ierr); 415058bd781SVaclav Hapla ierr = PetscStrallocpy("ir", &hdf5->mataij_jname);CHKERRQ(ierr); 416058bd781SVaclav Hapla ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr); 417058bd781SVaclav Hapla ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr); 418058bd781SVaclav Hapla 4190b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 420d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4210b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 42282ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4239a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 424058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 425058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4265c6c1daeSBarry Smith PetscFunctionReturn(0); 4275c6c1daeSBarry Smith } 4285c6c1daeSBarry Smith 4295c6c1daeSBarry Smith /*@C 4305c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4315c6c1daeSBarry Smith 4325c6c1daeSBarry Smith Collective on MPI_Comm 4335c6c1daeSBarry Smith 4345c6c1daeSBarry Smith Input Parameters: 4355c6c1daeSBarry Smith + comm - MPI communicator 4365c6c1daeSBarry Smith . name - name of file 4375c6c1daeSBarry Smith - type - type of file 4385c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4395c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4405c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4415c6c1daeSBarry Smith 4425c6c1daeSBarry Smith Output Parameter: 4435c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4445c6c1daeSBarry Smith 44582ea9b62SBarry Smith Options Database: 44682ea9b62SBarry 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 4479a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 44882ea9b62SBarry Smith 4495c6c1daeSBarry Smith Level: beginner 4505c6c1daeSBarry Smith 4515c6c1daeSBarry Smith Note: 4525c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4535c6c1daeSBarry Smith 4545c6c1daeSBarry Smith Concepts: HDF5 files 4555c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4565c6c1daeSBarry Smith 4576a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4589a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 459a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4605c6c1daeSBarry Smith @*/ 4615c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4625c6c1daeSBarry Smith { 4635c6c1daeSBarry Smith PetscErrorCode ierr; 4645c6c1daeSBarry Smith 4655c6c1daeSBarry Smith PetscFunctionBegin; 4665c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 4675c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 4685c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 4695c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 4705c6c1daeSBarry Smith PetscFunctionReturn(0); 4715c6c1daeSBarry Smith } 4725c6c1daeSBarry Smith 4735c6c1daeSBarry Smith /*@C 4745c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 4755c6c1daeSBarry Smith 4765c6c1daeSBarry Smith Not collective 4775c6c1daeSBarry Smith 4785c6c1daeSBarry Smith Input Parameter: 4795c6c1daeSBarry Smith . viewer - the PetscViewer 4805c6c1daeSBarry Smith 4815c6c1daeSBarry Smith Output Parameter: 4825c6c1daeSBarry Smith . file_id - The file id 4835c6c1daeSBarry Smith 4845c6c1daeSBarry Smith Level: intermediate 4855c6c1daeSBarry Smith 4865c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 4875c6c1daeSBarry Smith @*/ 4885c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 4895c6c1daeSBarry Smith { 4905c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4915c6c1daeSBarry Smith 4925c6c1daeSBarry Smith PetscFunctionBegin; 4935c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4945c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 4955c6c1daeSBarry Smith PetscFunctionReturn(0); 4965c6c1daeSBarry Smith } 4975c6c1daeSBarry Smith 4985c6c1daeSBarry Smith /*@C 4995c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5005c6c1daeSBarry Smith 5015c6c1daeSBarry Smith Not collective 5025c6c1daeSBarry Smith 5035c6c1daeSBarry Smith Input Parameters: 5045c6c1daeSBarry Smith + viewer - the PetscViewer 5055c6c1daeSBarry Smith - name - The group name 5065c6c1daeSBarry Smith 5075c6c1daeSBarry Smith Level: intermediate 5085c6c1daeSBarry Smith 509874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5105c6c1daeSBarry Smith @*/ 5115c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5125c6c1daeSBarry Smith { 5135c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5145c6c1daeSBarry Smith GroupList *groupNode; 5155c6c1daeSBarry Smith PetscErrorCode ierr; 5165c6c1daeSBarry Smith 5175c6c1daeSBarry Smith PetscFunctionBegin; 5185c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5195c6c1daeSBarry Smith PetscValidCharPointer(name,2); 52095dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5215c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 522a297a907SKarl Rupp 5235c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5245c6c1daeSBarry Smith hdf5->groups = groupNode; 5255c6c1daeSBarry Smith PetscFunctionReturn(0); 5265c6c1daeSBarry Smith } 5275c6c1daeSBarry Smith 5283ef9c667SSatish Balay /*@ 5295c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5305c6c1daeSBarry Smith 5315c6c1daeSBarry Smith Not collective 5325c6c1daeSBarry Smith 5335c6c1daeSBarry Smith Input Parameter: 5345c6c1daeSBarry Smith . viewer - the PetscViewer 5355c6c1daeSBarry Smith 5365c6c1daeSBarry Smith Level: intermediate 5375c6c1daeSBarry Smith 538874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5395c6c1daeSBarry Smith @*/ 5405c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5415c6c1daeSBarry Smith { 5425c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5435c6c1daeSBarry Smith GroupList *groupNode; 5445c6c1daeSBarry Smith PetscErrorCode ierr; 5455c6c1daeSBarry Smith 5465c6c1daeSBarry Smith PetscFunctionBegin; 5475c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 54882f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5495c6c1daeSBarry Smith groupNode = hdf5->groups; 5505c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5515c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5525c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5535c6c1daeSBarry Smith PetscFunctionReturn(0); 5545c6c1daeSBarry Smith } 5555c6c1daeSBarry Smith 5565c6c1daeSBarry Smith /*@C 557874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 558874270d9SVaclav Hapla If none has been assigned, returns NULL. 5595c6c1daeSBarry Smith 5605c6c1daeSBarry Smith Not collective 5615c6c1daeSBarry Smith 5625c6c1daeSBarry Smith Input Parameter: 5635c6c1daeSBarry Smith . viewer - the PetscViewer 5645c6c1daeSBarry Smith 5655c6c1daeSBarry Smith Output Parameter: 5665c6c1daeSBarry Smith . name - The group name 5675c6c1daeSBarry Smith 5685c6c1daeSBarry Smith Level: intermediate 5695c6c1daeSBarry Smith 570874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 5715c6c1daeSBarry Smith @*/ 5725c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 5735c6c1daeSBarry Smith { 5745c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 5755c6c1daeSBarry Smith 5765c6c1daeSBarry Smith PetscFunctionBegin; 5775c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 578c959eef4SJed Brown PetscValidPointer(name,2); 579a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 5800298fd71SBarry Smith else *name = NULL; 5815c6c1daeSBarry Smith PetscFunctionReturn(0); 5825c6c1daeSBarry Smith } 5835c6c1daeSBarry Smith 5848c557b5aSVaclav Hapla /*@ 585874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 586874270d9SVaclav Hapla and return this group's ID and file ID. 587874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 588874270d9SVaclav Hapla 589874270d9SVaclav Hapla Not collective 590874270d9SVaclav Hapla 591874270d9SVaclav Hapla Input Parameter: 592874270d9SVaclav Hapla . viewer - the PetscViewer 593874270d9SVaclav Hapla 594874270d9SVaclav Hapla Output Parameter: 595874270d9SVaclav Hapla + fileId - The HDF5 file ID 596874270d9SVaclav Hapla - groupId - The HDF5 group ID 597874270d9SVaclav Hapla 598874270d9SVaclav Hapla Level: intermediate 599874270d9SVaclav Hapla 600874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 601874270d9SVaclav Hapla @*/ 60254dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 60354dbf706SVaclav Hapla { 60454dbf706SVaclav Hapla hid_t file_id, group; 60554dbf706SVaclav Hapla htri_t found; 60654dbf706SVaclav Hapla const char *groupName = NULL; 60754dbf706SVaclav Hapla PetscErrorCode ierr; 60854dbf706SVaclav Hapla 60954dbf706SVaclav Hapla PetscFunctionBegin; 61054dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 61154dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 61254dbf706SVaclav Hapla /* Open group */ 61354dbf706SVaclav Hapla if (groupName) { 61454dbf706SVaclav Hapla PetscBool root; 61554dbf706SVaclav Hapla 61654dbf706SVaclav Hapla ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr); 61754dbf706SVaclav Hapla PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT)); 61854dbf706SVaclav Hapla if (!root && (found <= 0)) { 61954dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 62054dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT)); 62154dbf706SVaclav Hapla #else /* deprecated HDF5 1.6 API */ 62254dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate,(file_id, groupName, 0)); 62354dbf706SVaclav Hapla #endif 62454dbf706SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 62554dbf706SVaclav Hapla } 62654dbf706SVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 62754dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT)); 62854dbf706SVaclav Hapla #else 62954dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gopen,(file_id, groupName)); 63054dbf706SVaclav Hapla #endif 63154dbf706SVaclav Hapla } else group = file_id; 63254dbf706SVaclav Hapla 63354dbf706SVaclav Hapla *fileId = file_id; 63454dbf706SVaclav Hapla *groupId = group; 63554dbf706SVaclav Hapla PetscFunctionReturn(0); 63654dbf706SVaclav Hapla } 63754dbf706SVaclav Hapla 6383ef9c667SSatish Balay /*@ 6395c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith Not collective 6425c6c1daeSBarry Smith 6435c6c1daeSBarry Smith Input Parameter: 6445c6c1daeSBarry Smith . viewer - the PetscViewer 6455c6c1daeSBarry Smith 6465c6c1daeSBarry Smith Level: intermediate 6475c6c1daeSBarry Smith 6485c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6495c6c1daeSBarry Smith @*/ 6505c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6515c6c1daeSBarry Smith { 6525c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6535c6c1daeSBarry Smith 6545c6c1daeSBarry Smith PetscFunctionBegin; 6555c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6565c6c1daeSBarry Smith ++hdf5->timestep; 6575c6c1daeSBarry Smith PetscFunctionReturn(0); 6585c6c1daeSBarry Smith } 6595c6c1daeSBarry Smith 6603ef9c667SSatish Balay /*@ 6615c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6625c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith Not collective 6655c6c1daeSBarry Smith 6665c6c1daeSBarry Smith Input Parameters: 6675c6c1daeSBarry Smith + viewer - the PetscViewer 6685c6c1daeSBarry Smith - timestep - The timestep number 6695c6c1daeSBarry Smith 6705c6c1daeSBarry Smith Level: intermediate 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6735c6c1daeSBarry Smith @*/ 6745c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6755c6c1daeSBarry Smith { 6765c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6775c6c1daeSBarry Smith 6785c6c1daeSBarry Smith PetscFunctionBegin; 6795c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6805c6c1daeSBarry Smith hdf5->timestep = timestep; 6815c6c1daeSBarry Smith PetscFunctionReturn(0); 6825c6c1daeSBarry Smith } 6835c6c1daeSBarry Smith 6843ef9c667SSatish Balay /*@ 6855c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 6865c6c1daeSBarry Smith 6875c6c1daeSBarry Smith Not collective 6885c6c1daeSBarry Smith 6895c6c1daeSBarry Smith Input Parameter: 6905c6c1daeSBarry Smith . viewer - the PetscViewer 6915c6c1daeSBarry Smith 6925c6c1daeSBarry Smith Output Parameter: 6935c6c1daeSBarry Smith . timestep - The timestep number 6945c6c1daeSBarry Smith 6955c6c1daeSBarry Smith Level: intermediate 6965c6c1daeSBarry Smith 6975c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 6985c6c1daeSBarry Smith @*/ 6995c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7005c6c1daeSBarry Smith { 7015c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7025c6c1daeSBarry Smith 7035c6c1daeSBarry Smith PetscFunctionBegin; 7045c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7055c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7065c6c1daeSBarry Smith *timestep = hdf5->timestep; 7075c6c1daeSBarry Smith PetscFunctionReturn(0); 7085c6c1daeSBarry Smith } 7095c6c1daeSBarry Smith 71036bce6e8SMatthew G. Knepley /*@C 71136bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 71236bce6e8SMatthew G. Knepley 71336bce6e8SMatthew G. Knepley Not collective 71436bce6e8SMatthew G. Knepley 71536bce6e8SMatthew G. Knepley Input Parameter: 71636bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 71736bce6e8SMatthew G. Knepley 71836bce6e8SMatthew G. Knepley Output Parameter: 71936bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 72036bce6e8SMatthew G. Knepley 72136bce6e8SMatthew G. Knepley Level: advanced 72236bce6e8SMatthew G. Knepley 72336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 72436bce6e8SMatthew G. Knepley @*/ 72536bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 72636bce6e8SMatthew G. Knepley { 72736bce6e8SMatthew G. Knepley PetscFunctionBegin; 72836bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 72936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 73036bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 73136bce6e8SMatthew G. Knepley #else 73236bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 73336bce6e8SMatthew G. Knepley #endif 73436bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 73536bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 73636bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 73736bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 73836bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 73936bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 74036bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 74136bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7427e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 74336bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 74436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 74536bce6e8SMatthew G. Knepley } 74636bce6e8SMatthew G. Knepley 74736bce6e8SMatthew G. Knepley /*@C 74836bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 74936bce6e8SMatthew G. Knepley 75036bce6e8SMatthew G. Knepley Not collective 75136bce6e8SMatthew G. Knepley 75236bce6e8SMatthew G. Knepley Input Parameter: 75336bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 75436bce6e8SMatthew G. Knepley 75536bce6e8SMatthew G. Knepley Output Parameter: 75636bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 75736bce6e8SMatthew G. Knepley 75836bce6e8SMatthew G. Knepley Level: advanced 75936bce6e8SMatthew G. Knepley 76036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 76136bce6e8SMatthew G. Knepley @*/ 76236bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 76336bce6e8SMatthew G. Knepley { 76436bce6e8SMatthew G. Knepley PetscFunctionBegin; 76536bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 76636bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 76736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 76836bce6e8SMatthew G. Knepley #else 76936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 77036bce6e8SMatthew G. Knepley #endif 77136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 77236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 77336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 77436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 77536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 77636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7777e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 77836bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 77936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 78036bce6e8SMatthew G. Knepley } 78136bce6e8SMatthew G. Knepley 782df863907SAlex Fikl /*@C 78336bce6e8SMatthew G. Knepley PetscViewerHDF5WriteAttribute - Write a scalar attribute 78436bce6e8SMatthew G. Knepley 78536bce6e8SMatthew G. Knepley Input Parameters: 78636bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 78736bce6e8SMatthew G. Knepley . parent - The parent name 78836bce6e8SMatthew G. Knepley . name - The attribute name 78936bce6e8SMatthew G. Knepley . datatype - The attribute type 79036bce6e8SMatthew G. Knepley - value - The attribute value 79136bce6e8SMatthew G. Knepley 79236bce6e8SMatthew G. Knepley Level: advanced 79336bce6e8SMatthew G. Knepley 794e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute() 79536bce6e8SMatthew G. Knepley @*/ 79636bce6e8SMatthew G. Knepley PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 79736bce6e8SMatthew G. Knepley { 79860568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 79936bce6e8SMatthew G. Knepley PetscErrorCode ierr; 80036bce6e8SMatthew G. Knepley 80136bce6e8SMatthew G. Knepley PetscFunctionBegin; 8025cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 80336bce6e8SMatthew G. Knepley PetscValidPointer(parent, 2); 80436bce6e8SMatthew G. Knepley PetscValidPointer(name, 3); 80536bce6e8SMatthew G. Knepley PetscValidPointer(value, 4); 80636bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8077e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8087e97c476SMatthew G. Knepley size_t len; 8097e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 810729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8117e97c476SMatthew G. Knepley } 81236bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 813729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 81460568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 81536bce6e8SMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 81660568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 81736bce6e8SMatthew G. Knepley #else 81860568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate,(obj, name, dtype, dataspace, H5P_DEFAULT)); 81936bce6e8SMatthew G. Knepley #endif 820729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 821729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 822729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 82360568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 824729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 82536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 82636bce6e8SMatthew G. Knepley } 82736bce6e8SMatthew G. Knepley 828df863907SAlex Fikl /*@C 829ad0c61feSMatthew G. Knepley PetscViewerHDF5ReadAttribute - Read a scalar attribute 830ad0c61feSMatthew G. Knepley 831ad0c61feSMatthew G. Knepley Input Parameters: 832ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 833ad0c61feSMatthew G. Knepley . parent - The parent name 834ad0c61feSMatthew G. Knepley . name - The attribute name 835ad0c61feSMatthew G. Knepley - datatype - The attribute type 836ad0c61feSMatthew G. Knepley 837ad0c61feSMatthew G. Knepley Output Parameter: 838ad0c61feSMatthew G. Knepley . value - The attribute value 839ad0c61feSMatthew G. Knepley 840ad0c61feSMatthew G. Knepley Level: advanced 841ad0c61feSMatthew G. Knepley 842e7bdbf8eSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute() 843ad0c61feSMatthew G. Knepley @*/ 844ad0c61feSMatthew G. Knepley PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, void *value) 845ad0c61feSMatthew G. Knepley { 846f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 847ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 848ad0c61feSMatthew G. Knepley 849ad0c61feSMatthew G. Knepley PetscFunctionBegin; 8505cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 851ad0c61feSMatthew G. Knepley PetscValidPointer(parent, 2); 852ad0c61feSMatthew G. Knepley PetscValidPointer(name, 3); 853ad0c61feSMatthew G. Knepley PetscValidPointer(value, 4); 854ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 855ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 85660568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 85760568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 858f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 859f0b58479SMatthew G. Knepley size_t len; 860e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 861f0b58479SMatthew G. Knepley PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 862f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 863f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 864f0b58479SMatthew G. Knepley } 86570efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 866729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 867e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 868e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 869ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 870ad0c61feSMatthew G. Knepley } 871ad0c61feSMatthew G. Knepley 872*ecce1506SVaclav Hapla static PetscErrorCode PetscViewerHDF5HasObject_Internal(PetscViewer viewer, const char name[], H5O_type_t otype, PetscBool *has) 8735cdeb108SMatthew G. Knepley { 8745cdeb108SMatthew G. Knepley hid_t h5; 8755cdeb108SMatthew G. Knepley PetscErrorCode ierr; 8765cdeb108SMatthew G. Knepley 8775cdeb108SMatthew G. Knepley PetscFunctionBegin; 8785cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 8795cdeb108SMatthew G. Knepley PetscValidPointer(name, 2); 8805cdeb108SMatthew G. Knepley PetscValidPointer(has, 3); 8815cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 882c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 8835cdeb108SMatthew G. Knepley if (H5Lexists(h5, name, H5P_DEFAULT)) { 8845cdeb108SMatthew G. Knepley H5O_info_t info; 8855cdeb108SMatthew G. Knepley hid_t obj; 8865cdeb108SMatthew G. Knepley 887729ab6d9SBarry Smith PetscStackCallHDF5Return(obj,H5Oopen,(h5, name, H5P_DEFAULT)); 888729ab6d9SBarry Smith PetscStackCallHDF5(H5Oget_info,(obj, &info)); 8895cdeb108SMatthew G. Knepley if (otype == info.type) *has = PETSC_TRUE; 890729ab6d9SBarry Smith PetscStackCallHDF5(H5Oclose,(obj)); 8915cdeb108SMatthew G. Knepley } 8925cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 8935cdeb108SMatthew G. Knepley } 8945cdeb108SMatthew G. Knepley 895*ecce1506SVaclav Hapla /*@ 896*ecce1506SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file 897*ecce1506SVaclav Hapla 898*ecce1506SVaclav Hapla Input Parameters: 899*ecce1506SVaclav Hapla + viewer - The HDF5 viewer 900*ecce1506SVaclav Hapla - obj - The named object 901*ecce1506SVaclav Hapla 902*ecce1506SVaclav Hapla Output Parameter: 903*ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 904*ecce1506SVaclav Hapla 905*ecce1506SVaclav Hapla Level: advanced 906*ecce1506SVaclav Hapla 907*ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute() 908*ecce1506SVaclav Hapla @*/ 909*ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 910*ecce1506SVaclav Hapla { 911*ecce1506SVaclav Hapla PetscErrorCode ierr; 912*ecce1506SVaclav Hapla 913*ecce1506SVaclav Hapla PetscFunctionBegin; 914*ecce1506SVaclav Hapla *has = PETSC_FALSE; 915*ecce1506SVaclav Hapla if (obj->name) {ierr = PetscViewerHDF5HasObject_Internal(viewer, obj->name, H5O_TYPE_DATASET, has);CHKERRQ(ierr);} 916*ecce1506SVaclav Hapla PetscFunctionReturn(0); 917*ecce1506SVaclav Hapla } 918*ecce1506SVaclav Hapla 919df863907SAlex Fikl /*@C 920e7bdbf8eSMatthew G. Knepley PetscViewerHDF5HasAttribute - Check whether a scalar attribute exists 921e7bdbf8eSMatthew G. Knepley 922e7bdbf8eSMatthew G. Knepley Input Parameters: 923e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 924e7bdbf8eSMatthew G. Knepley . parent - The parent name 925e7bdbf8eSMatthew G. Knepley - name - The attribute name 926e7bdbf8eSMatthew G. Knepley 927e7bdbf8eSMatthew G. Knepley Output Parameter: 928e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 929e7bdbf8eSMatthew G. Knepley 930e7bdbf8eSMatthew G. Knepley Level: advanced 931e7bdbf8eSMatthew G. Knepley 932*ecce1506SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasObject() 933e7bdbf8eSMatthew G. Knepley @*/ 934e7bdbf8eSMatthew G. Knepley PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 935e7bdbf8eSMatthew G. Knepley { 936729ab6d9SBarry Smith hid_t h5, dataset; 9375cdeb108SMatthew G. Knepley htri_t hhas; 9385cdeb108SMatthew G. Knepley PetscBool exists; 939e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 940e7bdbf8eSMatthew G. Knepley 941e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 9425cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 943e7bdbf8eSMatthew G. Knepley PetscValidPointer(parent, 2); 944e7bdbf8eSMatthew G. Knepley PetscValidPointer(name, 3); 945e7bdbf8eSMatthew G. Knepley PetscValidPointer(has, 4); 946e7bdbf8eSMatthew G. Knepley *has = PETSC_FALSE; 947e7bdbf8eSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 948*ecce1506SVaclav Hapla ierr = PetscViewerHDF5HasObject_Internal(viewer, parent, H5O_TYPE_DATASET, &exists);CHKERRQ(ierr); 9495cdeb108SMatthew G. Knepley if (exists) { 950e7bdbf8eSMatthew G. Knepley #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 951729ab6d9SBarry Smith PetscStackCall("H5Dopen2",dataset = H5Dopen2(h5, parent, H5P_DEFAULT)); 952e7bdbf8eSMatthew G. Knepley #else 953729ab6d9SBarry Smith PetscStackCall("H5Dopen",dataset = H5Dopen(h5, parent)); 954e7bdbf8eSMatthew G. Knepley #endif 955729ab6d9SBarry Smith if (dataset < 0) PetscFunctionReturn(0); 956729ab6d9SBarry Smith PetscStackCall("H5Aexists",hhas = H5Aexists(dataset, name)); 957729ab6d9SBarry Smith if (hhas < 0) { 958729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 959729ab6d9SBarry Smith PetscFunctionReturn(0); 960729ab6d9SBarry Smith } 961729ab6d9SBarry Smith PetscStackCallHDF5(H5Dclose,(dataset)); 9625cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 9635cdeb108SMatthew G. Knepley } 964e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 965e7bdbf8eSMatthew G. Knepley } 966e7bdbf8eSMatthew G. Knepley 967eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 968a5e1feadSVaclav Hapla { 969fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 970fbd37863SVaclav Hapla const char *groupname=NULL; 97169a06e7bSVaclav Hapla char vecgroup[PETSC_MAX_PATH_LEN]; 972a5e1feadSVaclav Hapla PetscErrorCode ierr; 973a5e1feadSVaclav Hapla 974a5e1feadSVaclav Hapla PetscFunctionBegin; 97569a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 97669a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 977a5e1feadSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE >= 10800) 97869a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 979a5e1feadSVaclav Hapla #else 98069a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen,(h->group, name)); 981a5e1feadSVaclav Hapla #endif 98269a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 9839568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 98469a06e7bSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr); 985b44ade6dSVaclav Hapla ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr); 9869568d583SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr); 987b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 988b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 989b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 990b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 991b8ef406cSVaclav Hapla #endif 992b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 99369a06e7bSVaclav Hapla *ctx = h; 99469a06e7bSVaclav Hapla PetscFunctionReturn(0); 99569a06e7bSVaclav Hapla } 99669a06e7bSVaclav Hapla 997eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 99869a06e7bSVaclav Hapla { 99969a06e7bSVaclav Hapla HDF5ReadCtx h; 100069a06e7bSVaclav Hapla PetscErrorCode ierr; 100169a06e7bSVaclav Hapla 100269a06e7bSVaclav Hapla PetscFunctionBegin; 100369a06e7bSVaclav Hapla h = *ctx; 1004b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 100569a06e7bSVaclav Hapla if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group)); 100669a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 100769a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 100869a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 100969a06e7bSVaclav Hapla PetscFunctionReturn(0); 101069a06e7bSVaclav Hapla } 101169a06e7bSVaclav Hapla 1012eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 101369a06e7bSVaclav Hapla { 101469a06e7bSVaclav Hapla int rdim, dim; 101569a06e7bSVaclav Hapla hsize_t dims[4]; 101645e21b7fSVaclav Hapla PetscInt bsInd, lenInd, bs, N; 10178374c777SVaclav Hapla PetscLayout map; 10188374c777SVaclav Hapla PetscErrorCode ierr; 101969a06e7bSVaclav Hapla 102069a06e7bSVaclav Hapla PetscFunctionBegin; 10218374c777SVaclav Hapla if (!(*map_)) { 10228374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 10238374c777SVaclav Hapla } 10248374c777SVaclav Hapla map = *map_; 10259787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1026a5e1feadSVaclav Hapla dim = 0; 10279568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1028a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 10299568d583SVaclav Hapla if (ctx->complexVal) ++dim; 10309787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 103169a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 10329787f08bSVaclav Hapla /* calculate expected dimension indices */ 10339787f08bSVaclav Hapla lenInd = 0; 10349568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 10359787f08bSVaclav Hapla bsInd = lenInd + 1; 10369568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 10379787f08bSVaclav Hapla if (rdim == dim) { 103845e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 10399787f08bSVaclav Hapla } else if (rdim == dim+1) { 104045e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 10419568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1042a5e1feadSVaclav Hapla } else { 10439787f08bSVaclav Hapla SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim); 1044a5e1feadSVaclav Hapla } 104545e21b7fSVaclav Hapla N = (PetscInt) dims[lenInd]*bs; 10468374c777SVaclav Hapla 10478374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 10488374c777SVaclav Hapla if (map->bs < 0) { 104945e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 105045e21b7fSVaclav 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); 10518374c777SVaclav Hapla if (map->N < 0) { 105245e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 105345e21b7fSVaclav 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); 105469a06e7bSVaclav Hapla PetscFunctionReturn(0); 105569a06e7bSVaclav Hapla } 105669a06e7bSVaclav Hapla 1057eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 105816f6402dSVaclav Hapla { 105916f6402dSVaclav Hapla hsize_t count[4], offset[4]; 106016f6402dSVaclav Hapla int dim; 106116f6402dSVaclav Hapla PetscInt bs, n, low; 106216f6402dSVaclav Hapla PetscErrorCode ierr; 106316f6402dSVaclav Hapla 106416f6402dSVaclav Hapla PetscFunctionBegin; 106516f6402dSVaclav Hapla /* Compute local size and ownership range */ 106616f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 106716f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 106816f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 106916f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 107016f6402dSVaclav Hapla 107116f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 107216f6402dSVaclav Hapla dim = 0; 107316f6402dSVaclav Hapla if (ctx->timestep >= 0) { 107416f6402dSVaclav Hapla count[dim] = 1; 107516f6402dSVaclav Hapla offset[dim] = ctx->timestep; 107616f6402dSVaclav Hapla ++dim; 107716f6402dSVaclav Hapla } 107816f6402dSVaclav Hapla { 107916f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 108016f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 108116f6402dSVaclav Hapla ++dim; 108216f6402dSVaclav Hapla } 108316f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 108416f6402dSVaclav Hapla count[dim] = bs; 108516f6402dSVaclav Hapla offset[dim] = 0; 108616f6402dSVaclav Hapla ++dim; 108716f6402dSVaclav Hapla } 108816f6402dSVaclav Hapla if (ctx->complexVal) { 108916f6402dSVaclav Hapla count[dim] = 2; 109016f6402dSVaclav Hapla offset[dim] = 0; 109116f6402dSVaclav Hapla ++dim; 109216f6402dSVaclav Hapla } 109316f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 109416f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 109516f6402dSVaclav Hapla PetscFunctionReturn(0); 109616f6402dSVaclav Hapla } 109716f6402dSVaclav Hapla 1098eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1099ef2d82ceSVaclav Hapla { 1100ef2d82ceSVaclav Hapla PetscFunctionBegin; 1101ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1102ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1103ef2d82ceSVaclav Hapla } 1104ef2d82ceSVaclav Hapla 1105dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1106a25c73c6SVaclav Hapla { 1107fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1108fbd37863SVaclav Hapla hid_t memspace=0; 1109a25c73c6SVaclav Hapla size_t unitsize; 1110a25c73c6SVaclav Hapla void *arr; 1111a25c73c6SVaclav Hapla PetscErrorCode ierr; 1112a25c73c6SVaclav Hapla 1113a25c73c6SVaclav Hapla PetscFunctionBegin; 1114eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1115eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1116a25c73c6SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1117eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 11184fc17bcdSVaclav Hapla 11195a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 11205a89fdf4SVaclav Hapla if (!h->complexVal) { 1121c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1122c76551afSVaclav 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."); 11235a89fdf4SVaclav Hapla } 11245a89fdf4SVaclav Hapla #else 11255a89fdf4SVaclav 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."); 11265a89fdf4SVaclav Hapla #endif 11274fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 11284fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 11294fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 11304fc17bcdSVaclav Hapla 1131eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1132a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1133eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1134a25c73c6SVaclav Hapla *newarr = arr; 1135a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1136a25c73c6SVaclav Hapla } 1137a25c73c6SVaclav Hapla 1138c1aaad9cSVaclav Hapla /*@C 1139c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1140c1aaad9cSVaclav Hapla 1141c1aaad9cSVaclav Hapla Input Parameters: 1142c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1143c1aaad9cSVaclav Hapla - name - The vector name 1144c1aaad9cSVaclav Hapla 1145c1aaad9cSVaclav Hapla Output Parameter: 1146c1aaad9cSVaclav Hapla + bs - block size 1147c1aaad9cSVaclav Hapla - N - global size 1148c1aaad9cSVaclav Hapla 1149c1aaad9cSVaclav Hapla Note: 1150c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1151c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1152c1aaad9cSVaclav Hapla 1153c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1154c1aaad9cSVaclav Hapla 1155c1aaad9cSVaclav Hapla Level: advanced 1156c1aaad9cSVaclav Hapla 1157c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1158c1aaad9cSVaclav Hapla @*/ 115969a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 116069a06e7bSVaclav Hapla { 1161fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 11628374c777SVaclav Hapla PetscLayout map=NULL; 116369a06e7bSVaclav Hapla PetscErrorCode ierr; 116469a06e7bSVaclav Hapla 116569a06e7bSVaclav Hapla PetscFunctionBegin; 1166c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1167eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1168eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1169eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 11708374c777SVaclav Hapla if (bs) *bs = map->bs; 11718374c777SVaclav Hapla if (N) *N = map->N; 11728374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1173a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1174a5e1feadSVaclav Hapla } 1175a5e1feadSVaclav Hapla 1176a75e6a4aSMatthew G. Knepley /* 1177a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1178a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1179a75e6a4aSMatthew G. Knepley */ 1180d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1181a75e6a4aSMatthew G. Knepley 1182a75e6a4aSMatthew G. Knepley /*@C 1183a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1184a75e6a4aSMatthew G. Knepley 1185a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1186a75e6a4aSMatthew G. Knepley 1187a75e6a4aSMatthew G. Knepley Input Parameter: 1188a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1189a75e6a4aSMatthew G. Knepley 1190a75e6a4aSMatthew G. Knepley Level: intermediate 1191a75e6a4aSMatthew G. Knepley 1192a75e6a4aSMatthew G. Knepley Options Database Keys: 1193a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1194a75e6a4aSMatthew G. Knepley 1195a75e6a4aSMatthew G. Knepley Environmental variables: 1196a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1197a75e6a4aSMatthew G. Knepley 1198a75e6a4aSMatthew G. Knepley Notes: 1199a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1200a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1201a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1202a75e6a4aSMatthew G. Knepley 1203a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1204a75e6a4aSMatthew G. Knepley @*/ 1205a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1206a75e6a4aSMatthew G. Knepley { 1207a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1208a75e6a4aSMatthew G. Knepley PetscBool flg; 1209a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1210a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1211a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1212a75e6a4aSMatthew G. Knepley 1213a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1214a75e6a4aSMatthew 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);} 1215a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 121612801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1217a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1218a75e6a4aSMatthew G. Knepley } 121947435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1220a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1221a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1222a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1223a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1224a75e6a4aSMatthew G. Knepley if (!flg) { 1225a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1226a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1227a75e6a4aSMatthew G. Knepley } 1228a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1229a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1230a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1231a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 123247435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1233a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1234a75e6a4aSMatthew G. Knepley } 1235a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1236a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1237a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1238a75e6a4aSMatthew G. Knepley } 1239