1af0996ceSBarry Smith #include <petsc/private/viewerimpl.h> /*I "petscsys.h" I*/ 2d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 3ded06c3fSVaclav Hapla #if (H5_VERS_MAJOR * 10000 + H5_VERS_MINOR * 100 + H5_VERS_RELEASE < 10800) 4ded06c3fSVaclav Hapla #error "PETSc needs HDF5 version >= 1.8.0" 5ded06c3fSVaclav Hapla #endif 65c6c1daeSBarry Smith 7bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 906db490cSVaclav Hapla 105c6c1daeSBarry Smith typedef struct GroupList { 115c6c1daeSBarry Smith const char *name; 125c6c1daeSBarry Smith struct GroupList *next; 135c6c1daeSBarry Smith } GroupList; 145c6c1daeSBarry Smith 155c6c1daeSBarry Smith typedef struct { 165c6c1daeSBarry Smith char *filename; 175c6c1daeSBarry Smith PetscFileMode btype; 185c6c1daeSBarry Smith hid_t file_id; 195c6c1daeSBarry Smith PetscInt timestep; 205c6c1daeSBarry Smith GroupList *groups; 2182ea9b62SBarry Smith PetscBool basedimension2; /* save vectors and DMDA vectors with a dimension of at least 2 even if the bs/dof is 1 */ 229a0502c6SHåkon Strandenes PetscBool spoutput; /* write data in single precision even if PETSc is compiled with double precision PetscReal */ 23058bd781SVaclav Hapla char *mataij_iname; 24058bd781SVaclav Hapla char *mataij_jname; 25058bd781SVaclav Hapla char *mataij_aname; 26058bd781SVaclav Hapla char *mataij_cname; 275c6c1daeSBarry Smith } PetscViewer_HDF5; 285c6c1daeSBarry Smith 29eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx { 30eb5a92b4SVaclav Hapla hid_t file, group, dataset, dataspace, plist; 31eb5a92b4SVaclav Hapla PetscInt timestep; 3209dabeb0SVaclav Hapla PetscBool complexVal, dim2, horizontal; 33eb5a92b4SVaclav Hapla }; 34eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 35eb5a92b4SVaclav Hapla 366c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 376c132bc1SVaclav Hapla { 386c132bc1SVaclav Hapla const char *group; 396c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 406c132bc1SVaclav Hapla PetscErrorCode ierr; 416c132bc1SVaclav Hapla 426c132bc1SVaclav Hapla PetscFunctionBegin; 436c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 446c132bc1SVaclav Hapla ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 456c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 466c132bc1SVaclav Hapla ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 476c132bc1SVaclav Hapla ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 486c132bc1SVaclav Hapla PetscFunctionReturn(0); 496c132bc1SVaclav Hapla } 506c132bc1SVaclav Hapla 514416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 5282ea9b62SBarry Smith { 5382ea9b62SBarry Smith PetscErrorCode ierr; 5482ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 5582ea9b62SBarry Smith 5682ea9b62SBarry Smith PetscFunctionBegin; 5782ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 5882ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 599a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 6082ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 6182ea9b62SBarry Smith PetscFunctionReturn(0); 6282ea9b62SBarry Smith } 6382ea9b62SBarry Smith 645c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 655c6c1daeSBarry Smith { 665c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 675c6c1daeSBarry Smith PetscErrorCode ierr; 685c6c1daeSBarry Smith 695c6c1daeSBarry Smith PetscFunctionBegin; 705c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 71729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 725c6c1daeSBarry Smith PetscFunctionReturn(0); 735c6c1daeSBarry Smith } 745c6c1daeSBarry Smith 755c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 765c6c1daeSBarry Smith { 775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 785c6c1daeSBarry Smith PetscErrorCode ierr; 795c6c1daeSBarry Smith 805c6c1daeSBarry Smith PetscFunctionBegin; 815c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 825c6c1daeSBarry Smith while (hdf5->groups) { 835c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 845c6c1daeSBarry Smith 855c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 865c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 875c6c1daeSBarry Smith hdf5->groups = tmp; 885c6c1daeSBarry Smith } 8961912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 9061912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 9161912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 9261912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 935c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 940b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 95d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 960b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 97058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 98058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 99058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 100058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 1015c6c1daeSBarry Smith PetscFunctionReturn(0); 1025c6c1daeSBarry Smith } 1035c6c1daeSBarry Smith 1045c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1055c6c1daeSBarry Smith { 1065c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1075c6c1daeSBarry Smith 1085c6c1daeSBarry Smith PetscFunctionBegin; 1095c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1105c6c1daeSBarry Smith hdf5->btype = type; 1115c6c1daeSBarry Smith PetscFunctionReturn(0); 1125c6c1daeSBarry Smith } 1135c6c1daeSBarry Smith 11482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 11582ea9b62SBarry Smith { 11682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 11782ea9b62SBarry Smith 11882ea9b62SBarry Smith PetscFunctionBegin; 11982ea9b62SBarry Smith hdf5->basedimension2 = flg; 12082ea9b62SBarry Smith PetscFunctionReturn(0); 12182ea9b62SBarry Smith } 12282ea9b62SBarry Smith 123df863907SAlex Fikl /*@ 12482ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 12582ea9b62SBarry Smith dimension of 2. 12682ea9b62SBarry Smith 12782ea9b62SBarry Smith Logically Collective on PetscViewer 12882ea9b62SBarry Smith 12982ea9b62SBarry Smith Input Parameters: 13082ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 13182ea9b62SBarry 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 13282ea9b62SBarry Smith 13382ea9b62SBarry Smith Options Database: 13482ea9b62SBarry 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 13582ea9b62SBarry Smith 13682ea9b62SBarry Smith 13795452b02SPatrick Sanan Notes: 13895452b02SPatrick 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 13982ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith Level: intermediate 14282ea9b62SBarry Smith 14382ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 14482ea9b62SBarry Smith 14582ea9b62SBarry Smith @*/ 14682ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 14782ea9b62SBarry Smith { 14882ea9b62SBarry Smith PetscErrorCode ierr; 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith PetscFunctionBegin; 15182ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 15282ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 15382ea9b62SBarry Smith PetscFunctionReturn(0); 15482ea9b62SBarry Smith } 15582ea9b62SBarry Smith 156df863907SAlex Fikl /*@ 15782ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 15882ea9b62SBarry Smith dimension of 2. 15982ea9b62SBarry Smith 16082ea9b62SBarry Smith Logically Collective on PetscViewer 16182ea9b62SBarry Smith 16282ea9b62SBarry Smith Input Parameter: 16382ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 16482ea9b62SBarry Smith 16582ea9b62SBarry Smith Output Parameter: 16682ea9b62SBarry 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 16782ea9b62SBarry Smith 16895452b02SPatrick Sanan Notes: 16995452b02SPatrick 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 17082ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 17182ea9b62SBarry Smith 17282ea9b62SBarry Smith Level: intermediate 17382ea9b62SBarry Smith 17482ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 17582ea9b62SBarry Smith 17682ea9b62SBarry Smith @*/ 17782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 17882ea9b62SBarry Smith { 17982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 18082ea9b62SBarry Smith 18182ea9b62SBarry Smith PetscFunctionBegin; 18282ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 18382ea9b62SBarry Smith *flg = hdf5->basedimension2; 18482ea9b62SBarry Smith PetscFunctionReturn(0); 18582ea9b62SBarry Smith } 18682ea9b62SBarry Smith 1879a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1889a0502c6SHåkon Strandenes { 1899a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1909a0502c6SHåkon Strandenes 1919a0502c6SHåkon Strandenes PetscFunctionBegin; 1929a0502c6SHåkon Strandenes hdf5->spoutput = flg; 1939a0502c6SHåkon Strandenes PetscFunctionReturn(0); 1949a0502c6SHåkon Strandenes } 1959a0502c6SHåkon Strandenes 196df863907SAlex Fikl /*@ 1979a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 1989a0502c6SHåkon Strandenes compiled with double precision PetscReal. 1999a0502c6SHåkon Strandenes 2009a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2019a0502c6SHåkon Strandenes 2029a0502c6SHåkon Strandenes Input Parameters: 2039a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2049a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2059a0502c6SHåkon Strandenes 2069a0502c6SHåkon Strandenes Options Database: 2079a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2089a0502c6SHåkon Strandenes 2099a0502c6SHåkon Strandenes 21095452b02SPatrick Sanan Notes: 21195452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2129a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2139a0502c6SHåkon Strandenes 2149a0502c6SHåkon Strandenes Level: intermediate 2159a0502c6SHåkon Strandenes 2169a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2179a0502c6SHåkon Strandenes PetscReal 2189a0502c6SHåkon Strandenes 2199a0502c6SHåkon Strandenes @*/ 2209a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2219a0502c6SHåkon Strandenes { 2229a0502c6SHåkon Strandenes PetscErrorCode ierr; 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes PetscFunctionBegin; 2259a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2269a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2279a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2289a0502c6SHåkon Strandenes } 2299a0502c6SHåkon Strandenes 230df863907SAlex Fikl /*@ 2319a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2329a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2339a0502c6SHåkon Strandenes 2349a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2359a0502c6SHåkon Strandenes 2369a0502c6SHåkon Strandenes Input Parameter: 2379a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2389a0502c6SHåkon Strandenes 2399a0502c6SHåkon Strandenes Output Parameter: 2409a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2419a0502c6SHåkon Strandenes 24295452b02SPatrick Sanan Notes: 24395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2449a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2459a0502c6SHåkon Strandenes 2469a0502c6SHåkon Strandenes Level: intermediate 2479a0502c6SHåkon Strandenes 2489a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2499a0502c6SHåkon Strandenes PetscReal 2509a0502c6SHåkon Strandenes 2519a0502c6SHåkon Strandenes @*/ 2529a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2539a0502c6SHåkon Strandenes { 2549a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2559a0502c6SHåkon Strandenes 2569a0502c6SHåkon Strandenes PetscFunctionBegin; 2579a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2589a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2599a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2609a0502c6SHåkon Strandenes } 2619a0502c6SHåkon Strandenes 2625c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2635c6c1daeSBarry Smith { 2645c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2655c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2665c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2675c6c1daeSBarry Smith #endif 2685c6c1daeSBarry Smith hid_t plist_id; 2695c6c1daeSBarry Smith PetscErrorCode ierr; 2705c6c1daeSBarry Smith 2715c6c1daeSBarry Smith PetscFunctionBegin; 272f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 273f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2745c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2755c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 276729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2775c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 278729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2795c6c1daeSBarry Smith #endif 2805c6c1daeSBarry Smith /* Create or open the file collectively */ 2815c6c1daeSBarry Smith switch (hdf5->btype) { 2825c6c1daeSBarry Smith case FILE_MODE_READ: 283729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 2845c6c1daeSBarry Smith break; 2855c6c1daeSBarry Smith case FILE_MODE_APPEND: 286729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 2875c6c1daeSBarry Smith break; 2885c6c1daeSBarry Smith case FILE_MODE_WRITE: 289729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 2905c6c1daeSBarry Smith break; 2915c6c1daeSBarry Smith default: 2925c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 2935c6c1daeSBarry Smith } 2945c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 295729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 2965c6c1daeSBarry Smith PetscFunctionReturn(0); 2975c6c1daeSBarry Smith } 2985c6c1daeSBarry Smith 299d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 300d1232d7fSToby Isaac { 301d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 302d1232d7fSToby Isaac 303d1232d7fSToby Isaac PetscFunctionBegin; 304d1232d7fSToby Isaac *name = vhdf5->filename; 305d1232d7fSToby Isaac PetscFunctionReturn(0); 306d1232d7fSToby Isaac } 307d1232d7fSToby Isaac 308058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 309058bd781SVaclav Hapla { 310058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 311058bd781SVaclav Hapla PetscErrorCode ierr; 312058bd781SVaclav Hapla 313058bd781SVaclav Hapla PetscFunctionBegin; 314058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 315058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 316058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 317058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 318058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 319058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 320058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 321058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 322058bd781SVaclav Hapla PetscFunctionReturn(0); 323058bd781SVaclav Hapla } 324058bd781SVaclav Hapla 325058bd781SVaclav Hapla /*@C 326058bd781SVaclav 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. 327058bd781SVaclav Hapla 328058bd781SVaclav Hapla Collective on PetscViewer 329058bd781SVaclav Hapla 330058bd781SVaclav Hapla Input Parameters: 331058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 332058bd781SVaclav 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 333058bd781SVaclav Hapla . jname - name of dataset j representing column indices 334058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 335058bd781SVaclav Hapla - cname - name of attribute stoting column count 336058bd781SVaclav Hapla 337058bd781SVaclav Hapla Level: advanced 338058bd781SVaclav Hapla 339058bd781SVaclav Hapla Notes: 340058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 341058bd781SVaclav Hapla 342058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 343058bd781SVaclav Hapla @*/ 344058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 345058bd781SVaclav Hapla { 346058bd781SVaclav Hapla PetscErrorCode ierr; 347058bd781SVaclav Hapla 348058bd781SVaclav Hapla PetscFunctionBegin; 349058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 350058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 351058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 352058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 353058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 354058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 355058bd781SVaclav Hapla PetscFunctionReturn(0); 356058bd781SVaclav Hapla } 357058bd781SVaclav Hapla 358058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 359058bd781SVaclav Hapla { 360058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 361058bd781SVaclav Hapla 362058bd781SVaclav Hapla PetscFunctionBegin; 363058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 364058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 365058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 366058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 367058bd781SVaclav Hapla PetscFunctionReturn(0); 368058bd781SVaclav Hapla } 369058bd781SVaclav Hapla 370058bd781SVaclav Hapla /*@C 371058bd781SVaclav 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. 372058bd781SVaclav Hapla 373058bd781SVaclav Hapla Collective on PetscViewer 374058bd781SVaclav Hapla 375058bd781SVaclav Hapla Input Parameters: 376058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 377058bd781SVaclav Hapla 378058bd781SVaclav Hapla Output Parameters: 379058bd781SVaclav 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 380058bd781SVaclav Hapla . jname - name of dataset j representing column indices 381058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 382058bd781SVaclav Hapla - cname - name of attribute stoting column count 383058bd781SVaclav Hapla 384058bd781SVaclav Hapla Level: advanced 385058bd781SVaclav Hapla 386058bd781SVaclav Hapla Notes: 387058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 388058bd781SVaclav Hapla 389058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 390058bd781SVaclav Hapla @*/ 391058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 392058bd781SVaclav Hapla { 393058bd781SVaclav Hapla PetscErrorCode ierr; 394058bd781SVaclav Hapla 395058bd781SVaclav Hapla PetscFunctionBegin; 396058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 397058bd781SVaclav Hapla PetscValidPointer(iname,2); 398058bd781SVaclav Hapla PetscValidPointer(jname,3); 399058bd781SVaclav Hapla PetscValidPointer(aname,4); 400058bd781SVaclav Hapla PetscValidPointer(cname,5); 401058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 402058bd781SVaclav Hapla PetscFunctionReturn(0); 403058bd781SVaclav Hapla } 404058bd781SVaclav Hapla 4058556b5ebSBarry Smith /*MC 4068556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4078556b5ebSBarry Smith 4088556b5ebSBarry Smith 4098556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4108556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4118556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4128556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4138556b5ebSBarry Smith 4141b266c99SBarry Smith Level: beginner 4158556b5ebSBarry Smith M*/ 416d1232d7fSToby Isaac 4178cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4185c6c1daeSBarry Smith { 4195c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4205c6c1daeSBarry Smith PetscErrorCode ierr; 4215c6c1daeSBarry Smith 4225c6c1daeSBarry Smith PetscFunctionBegin; 423b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4245c6c1daeSBarry Smith 4255c6c1daeSBarry Smith v->data = (void*) hdf5; 4265c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 42782ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 4285c6c1daeSBarry Smith v->ops->flush = 0; 4295c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4305c6c1daeSBarry Smith hdf5->filename = 0; 4315c6c1daeSBarry Smith hdf5->timestep = -1; 4320298fd71SBarry Smith hdf5->groups = NULL; 4335c6c1daeSBarry Smith 434058bd781SVaclav Hapla /* ir and jc are deliberately swapped as MATLAB uses column-major format */ 435058bd781SVaclav Hapla ierr = PetscStrallocpy("jc", &hdf5->mataij_iname);CHKERRQ(ierr); 436058bd781SVaclav Hapla ierr = PetscStrallocpy("ir", &hdf5->mataij_jname);CHKERRQ(ierr); 437058bd781SVaclav Hapla ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr); 438058bd781SVaclav Hapla ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr); 439058bd781SVaclav Hapla 4400b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 441d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4420b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 44382ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4449a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 445058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 446058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4475c6c1daeSBarry Smith PetscFunctionReturn(0); 4485c6c1daeSBarry Smith } 4495c6c1daeSBarry Smith 4505c6c1daeSBarry Smith /*@C 4515c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4525c6c1daeSBarry Smith 4535c6c1daeSBarry Smith Collective on MPI_Comm 4545c6c1daeSBarry Smith 4555c6c1daeSBarry Smith Input Parameters: 4565c6c1daeSBarry Smith + comm - MPI communicator 4575c6c1daeSBarry Smith . name - name of file 4585c6c1daeSBarry Smith - type - type of file 4595c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4605c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4615c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4625c6c1daeSBarry Smith 4635c6c1daeSBarry Smith Output Parameter: 4645c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4655c6c1daeSBarry Smith 46682ea9b62SBarry Smith Options Database: 46782ea9b62SBarry 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 4689a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 46982ea9b62SBarry Smith 4705c6c1daeSBarry Smith Level: beginner 4715c6c1daeSBarry Smith 4725c6c1daeSBarry Smith Note: 4735c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4745c6c1daeSBarry Smith 4755c6c1daeSBarry Smith Concepts: HDF5 files 4765c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4775c6c1daeSBarry Smith 4786a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4799a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 480a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4815c6c1daeSBarry Smith @*/ 4825c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4835c6c1daeSBarry Smith { 4845c6c1daeSBarry Smith PetscErrorCode ierr; 4855c6c1daeSBarry Smith 4865c6c1daeSBarry Smith PetscFunctionBegin; 4875c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 4885c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 4895c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 4905c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 4915c6c1daeSBarry Smith PetscFunctionReturn(0); 4925c6c1daeSBarry Smith } 4935c6c1daeSBarry Smith 4945c6c1daeSBarry Smith /*@C 4955c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 4965c6c1daeSBarry Smith 4975c6c1daeSBarry Smith Not collective 4985c6c1daeSBarry Smith 4995c6c1daeSBarry Smith Input Parameter: 5005c6c1daeSBarry Smith . viewer - the PetscViewer 5015c6c1daeSBarry Smith 5025c6c1daeSBarry Smith Output Parameter: 5035c6c1daeSBarry Smith . file_id - The file id 5045c6c1daeSBarry Smith 5055c6c1daeSBarry Smith Level: intermediate 5065c6c1daeSBarry Smith 5075c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5085c6c1daeSBarry Smith @*/ 5095c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5105c6c1daeSBarry Smith { 5115c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5125c6c1daeSBarry Smith 5135c6c1daeSBarry Smith PetscFunctionBegin; 5145c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5155c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5165c6c1daeSBarry Smith PetscFunctionReturn(0); 5175c6c1daeSBarry Smith } 5185c6c1daeSBarry Smith 5195c6c1daeSBarry Smith /*@C 5205c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5215c6c1daeSBarry Smith 5225c6c1daeSBarry Smith Not collective 5235c6c1daeSBarry Smith 5245c6c1daeSBarry Smith Input Parameters: 5255c6c1daeSBarry Smith + viewer - the PetscViewer 5265c6c1daeSBarry Smith - name - The group name 5275c6c1daeSBarry Smith 5285c6c1daeSBarry Smith Level: intermediate 5295c6c1daeSBarry Smith 530874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5315c6c1daeSBarry Smith @*/ 5325c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5335c6c1daeSBarry Smith { 5345c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5355c6c1daeSBarry Smith GroupList *groupNode; 5365c6c1daeSBarry Smith PetscErrorCode ierr; 5375c6c1daeSBarry Smith 5385c6c1daeSBarry Smith PetscFunctionBegin; 5395c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5405c6c1daeSBarry Smith PetscValidCharPointer(name,2); 54195dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5425c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 543a297a907SKarl Rupp 5445c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5455c6c1daeSBarry Smith hdf5->groups = groupNode; 5465c6c1daeSBarry Smith PetscFunctionReturn(0); 5475c6c1daeSBarry Smith } 5485c6c1daeSBarry Smith 5493ef9c667SSatish Balay /*@ 5505c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5515c6c1daeSBarry Smith 5525c6c1daeSBarry Smith Not collective 5535c6c1daeSBarry Smith 5545c6c1daeSBarry Smith Input Parameter: 5555c6c1daeSBarry Smith . viewer - the PetscViewer 5565c6c1daeSBarry Smith 5575c6c1daeSBarry Smith Level: intermediate 5585c6c1daeSBarry Smith 559874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5605c6c1daeSBarry Smith @*/ 5615c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5625c6c1daeSBarry Smith { 5635c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5645c6c1daeSBarry Smith GroupList *groupNode; 5655c6c1daeSBarry Smith PetscErrorCode ierr; 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith PetscFunctionBegin; 5685c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 56982f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5705c6c1daeSBarry Smith groupNode = hdf5->groups; 5715c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5725c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5735c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5745c6c1daeSBarry Smith PetscFunctionReturn(0); 5755c6c1daeSBarry Smith } 5765c6c1daeSBarry Smith 5775c6c1daeSBarry Smith /*@C 578874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 579874270d9SVaclav Hapla If none has been assigned, returns NULL. 5805c6c1daeSBarry Smith 5815c6c1daeSBarry Smith Not collective 5825c6c1daeSBarry Smith 5835c6c1daeSBarry Smith Input Parameter: 5845c6c1daeSBarry Smith . viewer - the PetscViewer 5855c6c1daeSBarry Smith 5865c6c1daeSBarry Smith Output Parameter: 5875c6c1daeSBarry Smith . name - The group name 5885c6c1daeSBarry Smith 5895c6c1daeSBarry Smith Level: intermediate 5905c6c1daeSBarry Smith 591874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 5925c6c1daeSBarry Smith @*/ 5935c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 5945c6c1daeSBarry Smith { 5955c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 5965c6c1daeSBarry Smith 5975c6c1daeSBarry Smith PetscFunctionBegin; 5985c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 599c959eef4SJed Brown PetscValidPointer(name,2); 600a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6010298fd71SBarry Smith else *name = NULL; 6025c6c1daeSBarry Smith PetscFunctionReturn(0); 6035c6c1daeSBarry Smith } 6045c6c1daeSBarry Smith 6058c557b5aSVaclav Hapla /*@ 606874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 607874270d9SVaclav Hapla and return this group's ID and file ID. 608874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 609874270d9SVaclav Hapla 610874270d9SVaclav Hapla Not collective 611874270d9SVaclav Hapla 612874270d9SVaclav Hapla Input Parameter: 613874270d9SVaclav Hapla . viewer - the PetscViewer 614874270d9SVaclav Hapla 615874270d9SVaclav Hapla Output Parameter: 616874270d9SVaclav Hapla + fileId - The HDF5 file ID 617874270d9SVaclav Hapla - groupId - The HDF5 group ID 618874270d9SVaclav Hapla 619874270d9SVaclav Hapla Level: intermediate 620874270d9SVaclav Hapla 621874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 622874270d9SVaclav Hapla @*/ 62354dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 62454dbf706SVaclav Hapla { 62554dbf706SVaclav Hapla hid_t file_id, group; 62654dbf706SVaclav Hapla htri_t found; 62754dbf706SVaclav Hapla const char *groupName = NULL; 62854dbf706SVaclav Hapla PetscErrorCode ierr; 62954dbf706SVaclav Hapla 63054dbf706SVaclav Hapla PetscFunctionBegin; 63154dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 63254dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 63354dbf706SVaclav Hapla /* Open group */ 63454dbf706SVaclav Hapla if (groupName) { 63554dbf706SVaclav Hapla PetscBool root; 63654dbf706SVaclav Hapla 63754dbf706SVaclav Hapla ierr = PetscStrcmp(groupName, "/", &root);CHKERRQ(ierr); 63854dbf706SVaclav Hapla PetscStackCall("H5Lexists",found = H5Lexists(file_id, groupName, H5P_DEFAULT)); 63954dbf706SVaclav Hapla if (!root && (found <= 0)) { 64054dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(file_id, groupName, 0, H5P_DEFAULT, H5P_DEFAULT)); 64154dbf706SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 64254dbf706SVaclav Hapla } 64354dbf706SVaclav Hapla PetscStackCallHDF5Return(group,H5Gopen2,(file_id, groupName, H5P_DEFAULT)); 64454dbf706SVaclav Hapla } else group = file_id; 64554dbf706SVaclav Hapla 64654dbf706SVaclav Hapla *fileId = file_id; 64754dbf706SVaclav Hapla *groupId = group; 64854dbf706SVaclav Hapla PetscFunctionReturn(0); 64954dbf706SVaclav Hapla } 65054dbf706SVaclav Hapla 6513ef9c667SSatish Balay /*@ 6525c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6535c6c1daeSBarry Smith 6545c6c1daeSBarry Smith Not collective 6555c6c1daeSBarry Smith 6565c6c1daeSBarry Smith Input Parameter: 6575c6c1daeSBarry Smith . viewer - the PetscViewer 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith Level: intermediate 6605c6c1daeSBarry Smith 6615c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6625c6c1daeSBarry Smith @*/ 6635c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6645c6c1daeSBarry Smith { 6655c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6665c6c1daeSBarry Smith 6675c6c1daeSBarry Smith PetscFunctionBegin; 6685c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6695c6c1daeSBarry Smith ++hdf5->timestep; 6705c6c1daeSBarry Smith PetscFunctionReturn(0); 6715c6c1daeSBarry Smith } 6725c6c1daeSBarry Smith 6733ef9c667SSatish Balay /*@ 6745c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6755c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6765c6c1daeSBarry Smith 6775c6c1daeSBarry Smith Not collective 6785c6c1daeSBarry Smith 6795c6c1daeSBarry Smith Input Parameters: 6805c6c1daeSBarry Smith + viewer - the PetscViewer 6815c6c1daeSBarry Smith - timestep - The timestep number 6825c6c1daeSBarry Smith 6835c6c1daeSBarry Smith Level: intermediate 6845c6c1daeSBarry Smith 6855c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6865c6c1daeSBarry Smith @*/ 6875c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6885c6c1daeSBarry Smith { 6895c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6905c6c1daeSBarry Smith 6915c6c1daeSBarry Smith PetscFunctionBegin; 6925c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6935c6c1daeSBarry Smith hdf5->timestep = timestep; 6945c6c1daeSBarry Smith PetscFunctionReturn(0); 6955c6c1daeSBarry Smith } 6965c6c1daeSBarry Smith 6973ef9c667SSatish Balay /*@ 6985c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 6995c6c1daeSBarry Smith 7005c6c1daeSBarry Smith Not collective 7015c6c1daeSBarry Smith 7025c6c1daeSBarry Smith Input Parameter: 7035c6c1daeSBarry Smith . viewer - the PetscViewer 7045c6c1daeSBarry Smith 7055c6c1daeSBarry Smith Output Parameter: 7065c6c1daeSBarry Smith . timestep - The timestep number 7075c6c1daeSBarry Smith 7085c6c1daeSBarry Smith Level: intermediate 7095c6c1daeSBarry Smith 7105c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7115c6c1daeSBarry Smith @*/ 7125c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7135c6c1daeSBarry Smith { 7145c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7155c6c1daeSBarry Smith 7165c6c1daeSBarry Smith PetscFunctionBegin; 7175c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7185c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7195c6c1daeSBarry Smith *timestep = hdf5->timestep; 7205c6c1daeSBarry Smith PetscFunctionReturn(0); 7215c6c1daeSBarry Smith } 7225c6c1daeSBarry Smith 72336bce6e8SMatthew G. Knepley /*@C 72436bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 72536bce6e8SMatthew G. Knepley 72636bce6e8SMatthew G. Knepley Not collective 72736bce6e8SMatthew G. Knepley 72836bce6e8SMatthew G. Knepley Input Parameter: 72936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 73036bce6e8SMatthew G. Knepley 73136bce6e8SMatthew G. Knepley Output Parameter: 73236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 73336bce6e8SMatthew G. Knepley 73436bce6e8SMatthew G. Knepley Level: advanced 73536bce6e8SMatthew G. Knepley 73636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 73736bce6e8SMatthew G. Knepley @*/ 73836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 73936bce6e8SMatthew G. Knepley { 74036bce6e8SMatthew G. Knepley PetscFunctionBegin; 74136bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 74236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 74336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 74436bce6e8SMatthew G. Knepley #else 74536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 74636bce6e8SMatthew G. Knepley #endif 74736bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 74836bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 74936bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 75036bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 75136bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 75236bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 75336bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 75436bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7557e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 75636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 75736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 75836bce6e8SMatthew G. Knepley } 75936bce6e8SMatthew G. Knepley 76036bce6e8SMatthew G. Knepley /*@C 76136bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 76236bce6e8SMatthew G. Knepley 76336bce6e8SMatthew G. Knepley Not collective 76436bce6e8SMatthew G. Knepley 76536bce6e8SMatthew G. Knepley Input Parameter: 76636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 76736bce6e8SMatthew G. Knepley 76836bce6e8SMatthew G. Knepley Output Parameter: 76936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 77036bce6e8SMatthew G. Knepley 77136bce6e8SMatthew G. Knepley Level: advanced 77236bce6e8SMatthew G. Knepley 77336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 77436bce6e8SMatthew G. Knepley @*/ 77536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 77636bce6e8SMatthew G. Knepley { 77736bce6e8SMatthew G. Knepley PetscFunctionBegin; 77836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 77936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 78036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 78136bce6e8SMatthew G. Knepley #else 78236bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 78336bce6e8SMatthew G. Knepley #endif 78436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 78536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 78636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 78736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 78836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 78936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7907e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 79136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 79236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 79336bce6e8SMatthew G. Knepley } 79436bce6e8SMatthew G. Knepley 795df863907SAlex Fikl /*@C 796b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 79736bce6e8SMatthew G. Knepley 79836bce6e8SMatthew G. Knepley Input Parameters: 79936bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 800*57d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 80136bce6e8SMatthew G. Knepley . name - The attribute name 80236bce6e8SMatthew G. Knepley . datatype - The attribute type 80336bce6e8SMatthew G. Knepley - value - The attribute value 80436bce6e8SMatthew G. Knepley 80536bce6e8SMatthew G. Knepley Level: advanced 80636bce6e8SMatthew G. Knepley 807*57d22f7dSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 80836bce6e8SMatthew G. Knepley @*/ 809*57d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 81036bce6e8SMatthew G. Knepley { 811*57d22f7dSVaclav Hapla char *parent; 81260568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 813080f144cSVaclav Hapla PetscBool has; 81436bce6e8SMatthew G. Knepley PetscErrorCode ierr; 81536bce6e8SMatthew G. Knepley 81636bce6e8SMatthew G. Knepley PetscFunctionBegin; 8175cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 818*57d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 819c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 82036bce6e8SMatthew G. Knepley PetscValidPointer(value, 4); 821*57d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 822bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 823b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 82436bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8257e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8267e97c476SMatthew G. Knepley size_t len; 8277e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 828729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8297e97c476SMatthew G. Knepley } 83036bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 831729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 83260568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 833080f144cSVaclav Hapla if (has) { 834080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 835080f144cSVaclav Hapla } else { 83660568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 837080f144cSVaclav Hapla } 838729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 839729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 840729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 84160568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 842729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 843*57d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 84436bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 84536bce6e8SMatthew G. Knepley } 84636bce6e8SMatthew G. Knepley 847df863907SAlex Fikl /*@C 848b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 849ad0c61feSMatthew G. Knepley 850ad0c61feSMatthew G. Knepley Input Parameters: 851ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 852*57d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 853ad0c61feSMatthew G. Knepley . name - The attribute name 854ad0c61feSMatthew G. Knepley - datatype - The attribute type 855ad0c61feSMatthew G. Knepley 856ad0c61feSMatthew G. Knepley Output Parameter: 857ad0c61feSMatthew G. Knepley . value - The attribute value 858ad0c61feSMatthew G. Knepley 859ad0c61feSMatthew G. Knepley Level: advanced 860ad0c61feSMatthew G. Knepley 861*57d22f7dSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 862ad0c61feSMatthew G. Knepley @*/ 863*57d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 864ad0c61feSMatthew G. Knepley { 865*57d22f7dSVaclav Hapla char *parent; 866f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 867ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 868ad0c61feSMatthew G. Knepley 869ad0c61feSMatthew G. Knepley PetscFunctionBegin; 8705cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 871*57d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 872c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 873ad0c61feSMatthew G. Knepley PetscValidPointer(value, 4); 874*57d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 875ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 876ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 87760568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 87860568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 879f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 880f0b58479SMatthew G. Knepley size_t len; 881e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 882f0b58479SMatthew G. Knepley PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 883f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 884f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 885f0b58479SMatthew G. Knepley } 88670efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 887729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 888e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 889e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 890*57d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 891ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 892ad0c61feSMatthew G. Knepley } 893ad0c61feSMatthew G. Knepley 894a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 895a75c3fd4SVaclav Hapla { 896a75c3fd4SVaclav Hapla htri_t exists; 897a75c3fd4SVaclav Hapla hid_t group; 898a75c3fd4SVaclav Hapla 899a75c3fd4SVaclav Hapla PetscFunctionBegin; 900a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 901a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 902a75c3fd4SVaclav Hapla if (!exists && createGroup) { 903a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 904a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 905a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 906a75c3fd4SVaclav Hapla } 907a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 908a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 909a75c3fd4SVaclav Hapla } 910a75c3fd4SVaclav Hapla 911bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 9125cdeb108SMatthew G. Knepley { 9135cdeb108SMatthew G. Knepley hid_t h5; 914a75c3fd4SVaclav Hapla PetscBool exists; 91585688ae2SVaclav Hapla PetscInt i,n; 91685688ae2SVaclav Hapla char **hierarchy; 91785688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 9185cdeb108SMatthew G. Knepley PetscErrorCode ierr; 9195cdeb108SMatthew G. Knepley 9205cdeb108SMatthew G. Knepley PetscFunctionBegin; 9215cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 92200385fc5SVaclav Hapla PetscValidCharPointer(name, 2); 923ccd66a83SVaclav Hapla if (has) { 92456cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 9255cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 926ccd66a83SVaclav Hapla } 927ccd66a83SVaclav Hapla if (otype) { 928ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 92956cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 930ccd66a83SVaclav Hapla } 931c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 93285688ae2SVaclav Hapla 93385688ae2SVaclav Hapla /* 93485688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 93585688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 93685688ae2SVaclav Hapla 1) whether it's a valid link 93785688ae2SVaclav Hapla 2) whether this link resolves to an object 93885688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 93985688ae2SVaclav Hapla */ 94085688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 94185688ae2SVaclav Hapla if (!n) { 94285688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 943ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 944ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 94585688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 94685688ae2SVaclav Hapla PetscFunctionReturn(0); 94785688ae2SVaclav Hapla } 94885688ae2SVaclav Hapla for (i=0; i<n; i++) { 94985688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 95085688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 951a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 952a75c3fd4SVaclav Hapla if (!exists) break; 95385688ae2SVaclav Hapla } 95485688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 95585688ae2SVaclav Hapla 95685688ae2SVaclav Hapla /* If the object exists, get its type */ 957ccd66a83SVaclav Hapla if (exists && otype) { 9585cdeb108SMatthew G. Knepley H5O_info_t info; 9595cdeb108SMatthew G. Knepley 96004633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 96156cc0592SVaclav Hapla *otype = info.type; 9625cdeb108SMatthew G. Knepley } 963ccd66a83SVaclav Hapla if (has) *has = exists; 9645cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 9655cdeb108SMatthew G. Knepley } 9665cdeb108SMatthew G. Knepley 967ecce1506SVaclav Hapla /*@ 9680a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 9690a9f38efSVaclav Hapla 9700a9f38efSVaclav Hapla Input Parameters: 9710a9f38efSVaclav Hapla . viewer - The HDF5 viewer 9720a9f38efSVaclav Hapla 9730a9f38efSVaclav Hapla Output Parameter: 9740a9f38efSVaclav Hapla . has - Flag for group existence 9750a9f38efSVaclav Hapla 9760a9f38efSVaclav Hapla Notes: 9770a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 9780a9f38efSVaclav Hapla 9790a9f38efSVaclav Hapla Level: advanced 9800a9f38efSVaclav Hapla 9810a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 9820a9f38efSVaclav Hapla @*/ 9830a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 9840a9f38efSVaclav Hapla { 9850a9f38efSVaclav Hapla H5O_type_t type; 9860a9f38efSVaclav Hapla const char *name; 9870a9f38efSVaclav Hapla PetscErrorCode ierr; 9880a9f38efSVaclav Hapla 9890a9f38efSVaclav Hapla PetscFunctionBegin; 9900a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9910a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 9920a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 9930a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 9940a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 9950a9f38efSVaclav Hapla PetscFunctionReturn(0); 9960a9f38efSVaclav Hapla } 9970a9f38efSVaclav Hapla 9980a9f38efSVaclav Hapla /*@ 999e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1000ecce1506SVaclav Hapla 1001ecce1506SVaclav Hapla Input Parameters: 1002ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1003ecce1506SVaclav Hapla - obj - The named object 1004ecce1506SVaclav Hapla 1005ecce1506SVaclav Hapla Output Parameter: 1006ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1007ecce1506SVaclav Hapla 1008e3f143f7SVaclav Hapla Notes: 1009e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1010e3f143f7SVaclav Hapla 1011ecce1506SVaclav Hapla Level: advanced 1012ecce1506SVaclav Hapla 1013e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1014ecce1506SVaclav Hapla @*/ 1015ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1016ecce1506SVaclav Hapla { 101756cc0592SVaclav Hapla H5O_type_t type; 10186c132bc1SVaclav Hapla char *path; 1019ecce1506SVaclav Hapla PetscErrorCode ierr; 1020ecce1506SVaclav Hapla 1021ecce1506SVaclav Hapla PetscFunctionBegin; 1022c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1023c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1024c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1025ecce1506SVaclav Hapla *has = PETSC_FALSE; 1026e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 10276c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 10286c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 102956cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 10306c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1031ecce1506SVaclav Hapla PetscFunctionReturn(0); 1032ecce1506SVaclav Hapla } 1033ecce1506SVaclav Hapla 1034df863907SAlex Fikl /*@C 1035ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1036e7bdbf8eSMatthew G. Knepley 1037e7bdbf8eSMatthew G. Knepley Input Parameters: 1038e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 103910e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1040e7bdbf8eSMatthew G. Knepley - name - The attribute name 1041e7bdbf8eSMatthew G. Knepley 1042e7bdbf8eSMatthew G. Knepley Output Parameter: 1043e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1044e7bdbf8eSMatthew G. Knepley 1045e7bdbf8eSMatthew G. Knepley Level: advanced 1046e7bdbf8eSMatthew G. Knepley 104710e69e8fSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1048e7bdbf8eSMatthew G. Knepley @*/ 104910e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1050e7bdbf8eSMatthew G. Knepley { 105110e69e8fSVaclav Hapla char *parent; 1052e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1053e7bdbf8eSMatthew G. Knepley 1054e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 10555cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 105610e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1057c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1058c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 105910e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1060bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 106110e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 106210e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 106306db490cSVaclav Hapla PetscFunctionReturn(0); 106406db490cSVaclav Hapla } 106506db490cSVaclav Hapla 106606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 106706db490cSVaclav Hapla { 106806db490cSVaclav Hapla hid_t h5; 106906db490cSVaclav Hapla htri_t hhas; 107006db490cSVaclav Hapla PetscErrorCode ierr; 107106db490cSVaclav Hapla 107206db490cSVaclav Hapla PetscFunctionBegin; 107306db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 10742f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 10755cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1076e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1077e7bdbf8eSMatthew G. Knepley } 1078e7bdbf8eSMatthew G. Knepley 1079eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1080a5e1feadSVaclav Hapla { 1081fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1082fbd37863SVaclav Hapla const char *groupname=NULL; 108369a06e7bSVaclav Hapla char vecgroup[PETSC_MAX_PATH_LEN]; 1084a5e1feadSVaclav Hapla PetscErrorCode ierr; 1085a5e1feadSVaclav Hapla 1086a5e1feadSVaclav Hapla PetscFunctionBegin; 108769a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 108869a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 108969a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 109069a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 10919568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 109269a06e7bSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer,&groupname);CHKERRQ(ierr); 1093b44ade6dSVaclav Hapla ierr = PetscSNPrintf(vecgroup,PETSC_MAX_PATH_LEN,"%s/%s",groupname ? groupname : "",name);CHKERRQ(ierr); 10949568d583SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"complex",&h->complexVal);CHKERRQ(ierr); 109509dabeb0SVaclav Hapla /* MATLAB stores column vectors horizontally */ 109609dabeb0SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,vecgroup,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1097b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 1098b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1099b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1100b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1101b8ef406cSVaclav Hapla #endif 1102b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 110369a06e7bSVaclav Hapla *ctx = h; 110469a06e7bSVaclav Hapla PetscFunctionReturn(0); 110569a06e7bSVaclav Hapla } 110669a06e7bSVaclav Hapla 1107eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 110869a06e7bSVaclav Hapla { 110969a06e7bSVaclav Hapla HDF5ReadCtx h; 111069a06e7bSVaclav Hapla PetscErrorCode ierr; 111169a06e7bSVaclav Hapla 111269a06e7bSVaclav Hapla PetscFunctionBegin; 111369a06e7bSVaclav Hapla h = *ctx; 1114b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 111569a06e7bSVaclav Hapla if (h->group != h->file) PetscStackCallHDF5(H5Gclose,(h->group)); 111669a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 111769a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 111869a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 111969a06e7bSVaclav Hapla PetscFunctionReturn(0); 112069a06e7bSVaclav Hapla } 112169a06e7bSVaclav Hapla 1122eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 112369a06e7bSVaclav Hapla { 112469a06e7bSVaclav Hapla int rdim, dim; 112569a06e7bSVaclav Hapla hsize_t dims[4]; 112609dabeb0SVaclav Hapla PetscInt bsInd, lenInd, bs, len, N; 11278374c777SVaclav Hapla PetscLayout map; 11288374c777SVaclav Hapla PetscErrorCode ierr; 112969a06e7bSVaclav Hapla 113069a06e7bSVaclav Hapla PetscFunctionBegin; 11318374c777SVaclav Hapla if (!(*map_)) { 11328374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 11338374c777SVaclav Hapla } 11348374c777SVaclav Hapla map = *map_; 11359787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1136a5e1feadSVaclav Hapla dim = 0; 11379568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1138a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 11399568d583SVaclav Hapla if (ctx->complexVal) ++dim; 11409787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 114169a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 11429787f08bSVaclav Hapla /* calculate expected dimension indices */ 11439787f08bSVaclav Hapla lenInd = 0; 11449568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 11459787f08bSVaclav Hapla bsInd = lenInd + 1; 11469568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 11479787f08bSVaclav Hapla if (rdim == dim) { 114845e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 11499787f08bSVaclav Hapla } else if (rdim == dim+1) { 115045e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 11519568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1152a5e1feadSVaclav Hapla } else { 11539787f08bSVaclav Hapla SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_FILE_UNEXPECTED, "Dimension of array in file %d not %d as expected", rdim, dim); 1154a5e1feadSVaclav Hapla } 115509dabeb0SVaclav Hapla len = dims[lenInd]; 115609dabeb0SVaclav Hapla if (ctx->horizontal) { 115709dabeb0SVaclav 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."); 115809dabeb0SVaclav Hapla len = bs; 115909dabeb0SVaclav Hapla bs = 1; 116009dabeb0SVaclav Hapla } 116109dabeb0SVaclav Hapla N = (PetscInt) len*bs; 11628374c777SVaclav Hapla 11638374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 11648374c777SVaclav Hapla if (map->bs < 0) { 116545e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 116645e21b7fSVaclav 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); 11678374c777SVaclav Hapla if (map->N < 0) { 116845e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 116945e21b7fSVaclav 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); 117069a06e7bSVaclav Hapla PetscFunctionReturn(0); 117169a06e7bSVaclav Hapla } 117269a06e7bSVaclav Hapla 1173eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 117416f6402dSVaclav Hapla { 117516f6402dSVaclav Hapla hsize_t count[4], offset[4]; 117616f6402dSVaclav Hapla int dim; 117716f6402dSVaclav Hapla PetscInt bs, n, low; 117816f6402dSVaclav Hapla PetscErrorCode ierr; 117916f6402dSVaclav Hapla 118016f6402dSVaclav Hapla PetscFunctionBegin; 118116f6402dSVaclav Hapla /* Compute local size and ownership range */ 118216f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 118316f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 118416f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 118516f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 118616f6402dSVaclav Hapla 118716f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 118816f6402dSVaclav Hapla dim = 0; 118916f6402dSVaclav Hapla if (ctx->timestep >= 0) { 119016f6402dSVaclav Hapla count[dim] = 1; 119116f6402dSVaclav Hapla offset[dim] = ctx->timestep; 119216f6402dSVaclav Hapla ++dim; 119316f6402dSVaclav Hapla } 119409dabeb0SVaclav Hapla if (ctx->horizontal) { 119509dabeb0SVaclav Hapla count[dim] = 1; 119609dabeb0SVaclav Hapla offset[dim] = 0; 119709dabeb0SVaclav Hapla ++dim; 119809dabeb0SVaclav Hapla } 119916f6402dSVaclav Hapla { 120016f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 120116f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 120216f6402dSVaclav Hapla ++dim; 120316f6402dSVaclav Hapla } 120416f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 120509dabeb0SVaclav Hapla if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 120616f6402dSVaclav Hapla count[dim] = bs; 120716f6402dSVaclav Hapla offset[dim] = 0; 120816f6402dSVaclav Hapla ++dim; 120916f6402dSVaclav Hapla } 121016f6402dSVaclav Hapla if (ctx->complexVal) { 121116f6402dSVaclav Hapla count[dim] = 2; 121216f6402dSVaclav Hapla offset[dim] = 0; 121316f6402dSVaclav Hapla ++dim; 121416f6402dSVaclav Hapla } 121516f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 121616f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 121716f6402dSVaclav Hapla PetscFunctionReturn(0); 121816f6402dSVaclav Hapla } 121916f6402dSVaclav Hapla 1220eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1221ef2d82ceSVaclav Hapla { 1222ef2d82ceSVaclav Hapla PetscFunctionBegin; 1223ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1224ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1225ef2d82ceSVaclav Hapla } 1226ef2d82ceSVaclav Hapla 1227dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1228a25c73c6SVaclav Hapla { 1229fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1230fbd37863SVaclav Hapla hid_t memspace=0; 1231a25c73c6SVaclav Hapla size_t unitsize; 1232a25c73c6SVaclav Hapla void *arr; 1233a25c73c6SVaclav Hapla PetscErrorCode ierr; 1234a25c73c6SVaclav Hapla 1235a25c73c6SVaclav Hapla PetscFunctionBegin; 1236eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1237eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1238a25c73c6SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1239eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 12404fc17bcdSVaclav Hapla 12415a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 12425a89fdf4SVaclav Hapla if (!h->complexVal) { 1243c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1244c76551afSVaclav 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."); 12455a89fdf4SVaclav Hapla } 12465a89fdf4SVaclav Hapla #else 12475a89fdf4SVaclav 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."); 12485a89fdf4SVaclav Hapla #endif 12494fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 12504fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 1251dff35581SVaclav 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); 12524fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 12534fc17bcdSVaclav Hapla 1254eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1255a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1256eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1257a25c73c6SVaclav Hapla *newarr = arr; 1258a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1259a25c73c6SVaclav Hapla } 1260a25c73c6SVaclav Hapla 1261c1aaad9cSVaclav Hapla /*@C 1262c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1263c1aaad9cSVaclav Hapla 1264c1aaad9cSVaclav Hapla Input Parameters: 1265c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1266c1aaad9cSVaclav Hapla - name - The vector name 1267c1aaad9cSVaclav Hapla 1268c1aaad9cSVaclav Hapla Output Parameter: 1269c1aaad9cSVaclav Hapla + bs - block size 1270c1aaad9cSVaclav Hapla - N - global size 1271c1aaad9cSVaclav Hapla 1272c1aaad9cSVaclav Hapla Note: 1273c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1274c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1275c1aaad9cSVaclav Hapla 1276c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1277c1aaad9cSVaclav Hapla 1278c1aaad9cSVaclav Hapla Level: advanced 1279c1aaad9cSVaclav Hapla 1280c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1281c1aaad9cSVaclav Hapla @*/ 128269a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 128369a06e7bSVaclav Hapla { 1284fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 12858374c777SVaclav Hapla PetscLayout map=NULL; 128669a06e7bSVaclav Hapla PetscErrorCode ierr; 128769a06e7bSVaclav Hapla 128869a06e7bSVaclav Hapla PetscFunctionBegin; 1289c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1290eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1291eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1292eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 12938374c777SVaclav Hapla if (bs) *bs = map->bs; 12948374c777SVaclav Hapla if (N) *N = map->N; 12958374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1296a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1297a5e1feadSVaclav Hapla } 1298a5e1feadSVaclav Hapla 1299a75e6a4aSMatthew G. Knepley /* 1300a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1301a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1302a75e6a4aSMatthew G. Knepley */ 1303d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1304a75e6a4aSMatthew G. Knepley 1305a75e6a4aSMatthew G. Knepley /*@C 1306a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1307a75e6a4aSMatthew G. Knepley 1308a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1309a75e6a4aSMatthew G. Knepley 1310a75e6a4aSMatthew G. Knepley Input Parameter: 1311a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1312a75e6a4aSMatthew G. Knepley 1313a75e6a4aSMatthew G. Knepley Level: intermediate 1314a75e6a4aSMatthew G. Knepley 1315a75e6a4aSMatthew G. Knepley Options Database Keys: 1316a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1317a75e6a4aSMatthew G. Knepley 1318a75e6a4aSMatthew G. Knepley Environmental variables: 1319a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1320a75e6a4aSMatthew G. Knepley 1321a75e6a4aSMatthew G. Knepley Notes: 1322a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1323a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1324a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1325a75e6a4aSMatthew G. Knepley 1326a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1327a75e6a4aSMatthew G. Knepley @*/ 1328a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1329a75e6a4aSMatthew G. Knepley { 1330a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1331a75e6a4aSMatthew G. Knepley PetscBool flg; 1332a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1333a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1334a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1335a75e6a4aSMatthew G. Knepley 1336a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1337a75e6a4aSMatthew 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);} 1338a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 133912801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1340a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1341a75e6a4aSMatthew G. Knepley } 134247435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1343a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1344a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1345a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1346a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1347a75e6a4aSMatthew G. Knepley if (!flg) { 1348a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1349a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1350a75e6a4aSMatthew G. Knepley } 1351a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1352a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1353a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1354a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 135547435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1356a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1357a75e6a4aSMatthew G. Knepley } 1358a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1359a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1360a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1361a75e6a4aSMatthew G. Knepley } 1362