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 */ 235c6c1daeSBarry Smith } PetscViewer_HDF5; 245c6c1daeSBarry Smith 25eb5a92b4SVaclav Hapla struct _n_HDF5ReadCtx { 26eb5a92b4SVaclav Hapla hid_t file, group, dataset, dataspace, plist; 27eb5a92b4SVaclav Hapla PetscInt timestep; 2809dabeb0SVaclav Hapla PetscBool complexVal, dim2, horizontal; 29eb5a92b4SVaclav Hapla }; 30eb5a92b4SVaclav Hapla typedef struct _n_HDF5ReadCtx* HDF5ReadCtx; 31eb5a92b4SVaclav Hapla 326c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 336c132bc1SVaclav Hapla { 346c132bc1SVaclav Hapla const char *group; 356c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 366c132bc1SVaclav Hapla PetscErrorCode ierr; 376c132bc1SVaclav Hapla 386c132bc1SVaclav Hapla PetscFunctionBegin; 396c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 406c132bc1SVaclav Hapla ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 416c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 426c132bc1SVaclav Hapla ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 436c132bc1SVaclav Hapla ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 446c132bc1SVaclav Hapla PetscFunctionReturn(0); 456c132bc1SVaclav Hapla } 466c132bc1SVaclav Hapla 47577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 48577e0e04SVaclav Hapla { 49577e0e04SVaclav Hapla PetscBool has; 50577e0e04SVaclav Hapla const char *group; 51577e0e04SVaclav Hapla PetscErrorCode ierr; 52577e0e04SVaclav Hapla 53577e0e04SVaclav Hapla PetscFunctionBegin; 54577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 55577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 56577e0e04SVaclav Hapla if (!has) { 57577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 58577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 59577e0e04SVaclav Hapla } 60577e0e04SVaclav Hapla PetscFunctionReturn(0); 61577e0e04SVaclav Hapla } 62577e0e04SVaclav Hapla 634416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 6482ea9b62SBarry Smith { 6582ea9b62SBarry Smith PetscErrorCode ierr; 6682ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6782ea9b62SBarry Smith 6882ea9b62SBarry Smith PetscFunctionBegin; 6982ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 7082ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 719a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 7282ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7382ea9b62SBarry Smith PetscFunctionReturn(0); 7482ea9b62SBarry Smith } 7582ea9b62SBarry Smith 765c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 775c6c1daeSBarry Smith { 785c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 795c6c1daeSBarry Smith PetscErrorCode ierr; 805c6c1daeSBarry Smith 815c6c1daeSBarry Smith PetscFunctionBegin; 825c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 83729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 845c6c1daeSBarry Smith PetscFunctionReturn(0); 855c6c1daeSBarry Smith } 865c6c1daeSBarry Smith 879b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 885c6c1daeSBarry Smith { 895c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 905c6c1daeSBarry Smith PetscErrorCode ierr; 915c6c1daeSBarry Smith 925c6c1daeSBarry Smith PetscFunctionBegin; 935c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 945c6c1daeSBarry Smith while (hdf5->groups) { 955c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 965c6c1daeSBarry Smith 975c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 985c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 995c6c1daeSBarry Smith hdf5->groups = tmp; 1005c6c1daeSBarry Smith } 1015c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1020b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 103d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1040b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 105058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 106058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1075c6c1daeSBarry Smith PetscFunctionReturn(0); 1085c6c1daeSBarry Smith } 1095c6c1daeSBarry Smith 1109b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1115c6c1daeSBarry Smith { 1125c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1135c6c1daeSBarry Smith 1145c6c1daeSBarry Smith PetscFunctionBegin; 1155c6c1daeSBarry Smith hdf5->btype = type; 1165c6c1daeSBarry Smith PetscFunctionReturn(0); 1175c6c1daeSBarry Smith } 1185c6c1daeSBarry Smith 1190b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1200b62783dSVaclav Hapla { 1210b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1220b62783dSVaclav Hapla 1230b62783dSVaclav Hapla PetscFunctionBegin; 1240b62783dSVaclav Hapla *type = hdf5->btype; 1250b62783dSVaclav Hapla PetscFunctionReturn(0); 1260b62783dSVaclav Hapla } 1270b62783dSVaclav Hapla 1289b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 12982ea9b62SBarry Smith { 13082ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13182ea9b62SBarry Smith 13282ea9b62SBarry Smith PetscFunctionBegin; 13382ea9b62SBarry Smith hdf5->basedimension2 = flg; 13482ea9b62SBarry Smith PetscFunctionReturn(0); 13582ea9b62SBarry Smith } 13682ea9b62SBarry Smith 137df863907SAlex Fikl /*@ 13882ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 13982ea9b62SBarry Smith dimension of 2. 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith Logically Collective on PetscViewer 14282ea9b62SBarry Smith 14382ea9b62SBarry Smith Input Parameters: 14482ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 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 14782ea9b62SBarry Smith Options Database: 14882ea9b62SBarry 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 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith 15195452b02SPatrick Sanan Notes: 15295452b02SPatrick 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 15382ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15482ea9b62SBarry Smith 15582ea9b62SBarry Smith Level: intermediate 15682ea9b62SBarry Smith 15782ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith @*/ 16082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16182ea9b62SBarry Smith { 16282ea9b62SBarry Smith PetscErrorCode ierr; 16382ea9b62SBarry Smith 16482ea9b62SBarry Smith PetscFunctionBegin; 16582ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16682ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16782ea9b62SBarry Smith PetscFunctionReturn(0); 16882ea9b62SBarry Smith } 16982ea9b62SBarry Smith 170df863907SAlex Fikl /*@ 17182ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17282ea9b62SBarry Smith dimension of 2. 17382ea9b62SBarry Smith 17482ea9b62SBarry Smith Logically Collective on PetscViewer 17582ea9b62SBarry Smith 17682ea9b62SBarry Smith Input Parameter: 17782ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 17882ea9b62SBarry Smith 17982ea9b62SBarry Smith Output Parameter: 18082ea9b62SBarry 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 18182ea9b62SBarry Smith 18295452b02SPatrick Sanan Notes: 18395452b02SPatrick 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 18482ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18582ea9b62SBarry Smith 18682ea9b62SBarry Smith Level: intermediate 18782ea9b62SBarry Smith 18882ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 18982ea9b62SBarry Smith 19082ea9b62SBarry Smith @*/ 19182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19282ea9b62SBarry Smith { 19382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19482ea9b62SBarry Smith 19582ea9b62SBarry Smith PetscFunctionBegin; 19682ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19782ea9b62SBarry Smith *flg = hdf5->basedimension2; 19882ea9b62SBarry Smith PetscFunctionReturn(0); 19982ea9b62SBarry Smith } 20082ea9b62SBarry Smith 2019b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2029a0502c6SHåkon Strandenes { 2039a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2049a0502c6SHåkon Strandenes 2059a0502c6SHåkon Strandenes PetscFunctionBegin; 2069a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2079a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2089a0502c6SHåkon Strandenes } 2099a0502c6SHåkon Strandenes 210df863907SAlex Fikl /*@ 2119a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2129a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2139a0502c6SHåkon Strandenes 2149a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2159a0502c6SHåkon Strandenes 2169a0502c6SHåkon Strandenes Input Parameters: 2179a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2189a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2199a0502c6SHåkon Strandenes 2209a0502c6SHåkon Strandenes Options Database: 2219a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2229a0502c6SHåkon Strandenes 2239a0502c6SHåkon Strandenes 22495452b02SPatrick Sanan Notes: 22595452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2269a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2279a0502c6SHåkon Strandenes 2289a0502c6SHåkon Strandenes Level: intermediate 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2319a0502c6SHåkon Strandenes PetscReal 2329a0502c6SHåkon Strandenes 2339a0502c6SHåkon Strandenes @*/ 2349a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2359a0502c6SHåkon Strandenes { 2369a0502c6SHåkon Strandenes PetscErrorCode ierr; 2379a0502c6SHåkon Strandenes 2389a0502c6SHåkon Strandenes PetscFunctionBegin; 2399a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2409a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2419a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2429a0502c6SHåkon Strandenes } 2439a0502c6SHåkon Strandenes 244df863907SAlex Fikl /*@ 2459a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2469a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2479a0502c6SHåkon Strandenes 2489a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2499a0502c6SHåkon Strandenes 2509a0502c6SHåkon Strandenes Input Parameter: 2519a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2529a0502c6SHåkon Strandenes 2539a0502c6SHåkon Strandenes Output Parameter: 2549a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2559a0502c6SHåkon Strandenes 25695452b02SPatrick Sanan Notes: 25795452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2589a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2599a0502c6SHåkon Strandenes 2609a0502c6SHåkon Strandenes Level: intermediate 2619a0502c6SHåkon Strandenes 2629a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2639a0502c6SHåkon Strandenes PetscReal 2649a0502c6SHåkon Strandenes 2659a0502c6SHåkon Strandenes @*/ 2669a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2679a0502c6SHåkon Strandenes { 2689a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2699a0502c6SHåkon Strandenes 2709a0502c6SHåkon Strandenes PetscFunctionBegin; 2719a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2729a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2739a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2749a0502c6SHåkon Strandenes } 2759a0502c6SHåkon Strandenes 2769b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2775c6c1daeSBarry Smith { 2785c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2795c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2805c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2815c6c1daeSBarry Smith #endif 2825c6c1daeSBarry Smith hid_t plist_id; 2835c6c1daeSBarry Smith PetscErrorCode ierr; 2845c6c1daeSBarry Smith 2855c6c1daeSBarry Smith PetscFunctionBegin; 286f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 287f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2885c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2895c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 290729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2915c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 292729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2935c6c1daeSBarry Smith #endif 2945c6c1daeSBarry Smith /* Create or open the file collectively */ 2955c6c1daeSBarry Smith switch (hdf5->btype) { 2965c6c1daeSBarry Smith case FILE_MODE_READ: 297729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 2985c6c1daeSBarry Smith break; 2995c6c1daeSBarry Smith case FILE_MODE_APPEND: 300729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 3015c6c1daeSBarry Smith break; 3025c6c1daeSBarry Smith case FILE_MODE_WRITE: 303729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 3045c6c1daeSBarry Smith break; 3055c6c1daeSBarry Smith default: 3065c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3075c6c1daeSBarry Smith } 3085c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 309729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 3105c6c1daeSBarry Smith PetscFunctionReturn(0); 3115c6c1daeSBarry Smith } 3125c6c1daeSBarry Smith 313d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 314d1232d7fSToby Isaac { 315d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 316d1232d7fSToby Isaac 317d1232d7fSToby Isaac PetscFunctionBegin; 318d1232d7fSToby Isaac *name = vhdf5->filename; 319d1232d7fSToby Isaac PetscFunctionReturn(0); 320d1232d7fSToby Isaac } 321d1232d7fSToby Isaac 322b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 323b723ab35SVaclav Hapla { 324*5dc64a97SVaclav Hapla /* 325b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 326b723ab35SVaclav Hapla PetscErrorCode ierr; 327*5dc64a97SVaclav Hapla */ 328b723ab35SVaclav Hapla 329b723ab35SVaclav Hapla PetscFunctionBegin; 330b723ab35SVaclav Hapla PetscFunctionReturn(0); 331b723ab35SVaclav Hapla } 332b723ab35SVaclav Hapla 3338556b5ebSBarry Smith /*MC 3348556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 3358556b5ebSBarry Smith 3368556b5ebSBarry Smith 3378556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 3388556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 3398556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 3408556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 3418556b5ebSBarry Smith 3421b266c99SBarry Smith Level: beginner 3438556b5ebSBarry Smith M*/ 344d1232d7fSToby Isaac 3458cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 3465c6c1daeSBarry Smith { 3475c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 3485c6c1daeSBarry Smith PetscErrorCode ierr; 3495c6c1daeSBarry Smith 3505c6c1daeSBarry Smith PetscFunctionBegin; 351b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 3525c6c1daeSBarry Smith 3535c6c1daeSBarry Smith v->data = (void*) hdf5; 3545c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 35582ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 356b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 3575c6c1daeSBarry Smith v->ops->flush = 0; 3585c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 3595c6c1daeSBarry Smith hdf5->filename = 0; 3605c6c1daeSBarry Smith hdf5->timestep = -1; 3610298fd71SBarry Smith hdf5->groups = NULL; 3625c6c1daeSBarry Smith 3630b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 364d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 3650b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 3660b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 36782ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 3689a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 3695c6c1daeSBarry Smith PetscFunctionReturn(0); 3705c6c1daeSBarry Smith } 3715c6c1daeSBarry Smith 3725c6c1daeSBarry Smith /*@C 3735c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 3745c6c1daeSBarry Smith 3755c6c1daeSBarry Smith Collective on MPI_Comm 3765c6c1daeSBarry Smith 3775c6c1daeSBarry Smith Input Parameters: 3785c6c1daeSBarry Smith + comm - MPI communicator 3795c6c1daeSBarry Smith . name - name of file 3805c6c1daeSBarry Smith - type - type of file 3815c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 3825c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 3835c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 3845c6c1daeSBarry Smith 3855c6c1daeSBarry Smith Output Parameter: 3865c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 3875c6c1daeSBarry Smith 38882ea9b62SBarry Smith Options Database: 38982ea9b62SBarry 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 3909a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 39182ea9b62SBarry Smith 3925c6c1daeSBarry Smith Level: beginner 3935c6c1daeSBarry Smith 3945c6c1daeSBarry Smith Note: 3955c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 3965c6c1daeSBarry Smith 3975c6c1daeSBarry Smith Concepts: HDF5 files 3985c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 3995c6c1daeSBarry Smith 4006a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4019a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 402a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4035c6c1daeSBarry Smith @*/ 4045c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4055c6c1daeSBarry Smith { 4065c6c1daeSBarry Smith PetscErrorCode ierr; 4075c6c1daeSBarry Smith 4085c6c1daeSBarry Smith PetscFunctionBegin; 4095c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 4105c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 4115c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 4125c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 4135c6c1daeSBarry Smith PetscFunctionReturn(0); 4145c6c1daeSBarry Smith } 4155c6c1daeSBarry Smith 4165c6c1daeSBarry Smith /*@C 4175c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 4185c6c1daeSBarry Smith 4195c6c1daeSBarry Smith Not collective 4205c6c1daeSBarry Smith 4215c6c1daeSBarry Smith Input Parameter: 4225c6c1daeSBarry Smith . viewer - the PetscViewer 4235c6c1daeSBarry Smith 4245c6c1daeSBarry Smith Output Parameter: 4255c6c1daeSBarry Smith . file_id - The file id 4265c6c1daeSBarry Smith 4275c6c1daeSBarry Smith Level: intermediate 4285c6c1daeSBarry Smith 4295c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 4305c6c1daeSBarry Smith @*/ 4315c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 4325c6c1daeSBarry Smith { 4335c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4345c6c1daeSBarry Smith 4355c6c1daeSBarry Smith PetscFunctionBegin; 4365c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4375c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 4385c6c1daeSBarry Smith PetscFunctionReturn(0); 4395c6c1daeSBarry Smith } 4405c6c1daeSBarry Smith 4415c6c1daeSBarry Smith /*@C 4425c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 4435c6c1daeSBarry Smith 4445c6c1daeSBarry Smith Not collective 4455c6c1daeSBarry Smith 4465c6c1daeSBarry Smith Input Parameters: 4475c6c1daeSBarry Smith + viewer - the PetscViewer 4485c6c1daeSBarry Smith - name - The group name 4495c6c1daeSBarry Smith 4505c6c1daeSBarry Smith Level: intermediate 4515c6c1daeSBarry Smith 452874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 4535c6c1daeSBarry Smith @*/ 4545c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 4555c6c1daeSBarry Smith { 4565c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4575c6c1daeSBarry Smith GroupList *groupNode; 4585c6c1daeSBarry Smith PetscErrorCode ierr; 4595c6c1daeSBarry Smith 4605c6c1daeSBarry Smith PetscFunctionBegin; 4615c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 4625c6c1daeSBarry Smith PetscValidCharPointer(name,2); 46395dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 4645c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 465a297a907SKarl Rupp 4665c6c1daeSBarry Smith groupNode->next = hdf5->groups; 4675c6c1daeSBarry Smith hdf5->groups = groupNode; 4685c6c1daeSBarry Smith PetscFunctionReturn(0); 4695c6c1daeSBarry Smith } 4705c6c1daeSBarry Smith 4713ef9c667SSatish Balay /*@ 4725c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 4735c6c1daeSBarry Smith 4745c6c1daeSBarry Smith Not collective 4755c6c1daeSBarry Smith 4765c6c1daeSBarry Smith Input Parameter: 4775c6c1daeSBarry Smith . viewer - the PetscViewer 4785c6c1daeSBarry Smith 4795c6c1daeSBarry Smith Level: intermediate 4805c6c1daeSBarry Smith 481874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 4825c6c1daeSBarry Smith @*/ 4835c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 4845c6c1daeSBarry Smith { 4855c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4865c6c1daeSBarry Smith GroupList *groupNode; 4875c6c1daeSBarry Smith PetscErrorCode ierr; 4885c6c1daeSBarry Smith 4895c6c1daeSBarry Smith PetscFunctionBegin; 4905c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 49182f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 4925c6c1daeSBarry Smith groupNode = hdf5->groups; 4935c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 4945c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 4955c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 4965c6c1daeSBarry Smith PetscFunctionReturn(0); 4975c6c1daeSBarry Smith } 4985c6c1daeSBarry Smith 4995c6c1daeSBarry Smith /*@C 500874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 501874270d9SVaclav Hapla If none has been assigned, returns NULL. 5025c6c1daeSBarry Smith 5035c6c1daeSBarry Smith Not collective 5045c6c1daeSBarry Smith 5055c6c1daeSBarry Smith Input Parameter: 5065c6c1daeSBarry Smith . viewer - the PetscViewer 5075c6c1daeSBarry Smith 5085c6c1daeSBarry Smith Output Parameter: 5095c6c1daeSBarry Smith . name - The group name 5105c6c1daeSBarry Smith 5115c6c1daeSBarry Smith Level: intermediate 5125c6c1daeSBarry Smith 513874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 5145c6c1daeSBarry Smith @*/ 5155c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 5165c6c1daeSBarry Smith { 5175c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 5185c6c1daeSBarry Smith 5195c6c1daeSBarry Smith PetscFunctionBegin; 5205c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 521c959eef4SJed Brown PetscValidPointer(name,2); 522a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 5230298fd71SBarry Smith else *name = NULL; 5245c6c1daeSBarry Smith PetscFunctionReturn(0); 5255c6c1daeSBarry Smith } 5265c6c1daeSBarry Smith 5278c557b5aSVaclav Hapla /*@ 528874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 529874270d9SVaclav Hapla and return this group's ID and file ID. 530874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 531874270d9SVaclav Hapla 532874270d9SVaclav Hapla Not collective 533874270d9SVaclav Hapla 534874270d9SVaclav Hapla Input Parameter: 535874270d9SVaclav Hapla . viewer - the PetscViewer 536874270d9SVaclav Hapla 537874270d9SVaclav Hapla Output Parameter: 538874270d9SVaclav Hapla + fileId - The HDF5 file ID 539874270d9SVaclav Hapla - groupId - The HDF5 group ID 540874270d9SVaclav Hapla 541e74616beSVaclav Hapla Notes: 542e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 543e74616beSVaclav Hapla 544874270d9SVaclav Hapla Level: intermediate 545874270d9SVaclav Hapla 546874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 547874270d9SVaclav Hapla @*/ 54854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 54954dbf706SVaclav Hapla { 55090e03537SVaclav Hapla hid_t file_id; 55190e03537SVaclav Hapla H5O_type_t type; 55254dbf706SVaclav Hapla const char *groupName = NULL; 553e74616beSVaclav Hapla PetscBool create; 55454dbf706SVaclav Hapla PetscErrorCode ierr; 55554dbf706SVaclav Hapla 55654dbf706SVaclav Hapla PetscFunctionBegin; 557e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 55854dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 55954dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 560e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 56190e03537SVaclav Hapla if (type != H5O_TYPE_GROUP) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s resolves to something which is not a group", groupName); 56290e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 56354dbf706SVaclav Hapla *fileId = file_id; 56454dbf706SVaclav Hapla PetscFunctionReturn(0); 56554dbf706SVaclav Hapla } 56654dbf706SVaclav Hapla 5673ef9c667SSatish Balay /*@ 5685c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 5695c6c1daeSBarry Smith 5705c6c1daeSBarry Smith Not collective 5715c6c1daeSBarry Smith 5725c6c1daeSBarry Smith Input Parameter: 5735c6c1daeSBarry Smith . viewer - the PetscViewer 5745c6c1daeSBarry Smith 5755c6c1daeSBarry Smith Level: intermediate 5765c6c1daeSBarry Smith 5775c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 5785c6c1daeSBarry Smith @*/ 5795c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 5805c6c1daeSBarry Smith { 5815c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5825c6c1daeSBarry Smith 5835c6c1daeSBarry Smith PetscFunctionBegin; 5845c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5855c6c1daeSBarry Smith ++hdf5->timestep; 5865c6c1daeSBarry Smith PetscFunctionReturn(0); 5875c6c1daeSBarry Smith } 5885c6c1daeSBarry Smith 5893ef9c667SSatish Balay /*@ 5905c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 5915c6c1daeSBarry Smith of -1 disables blocking with timesteps. 5925c6c1daeSBarry Smith 5935c6c1daeSBarry Smith Not collective 5945c6c1daeSBarry Smith 5955c6c1daeSBarry Smith Input Parameters: 5965c6c1daeSBarry Smith + viewer - the PetscViewer 5975c6c1daeSBarry Smith - timestep - The timestep number 5985c6c1daeSBarry Smith 5995c6c1daeSBarry Smith Level: intermediate 6005c6c1daeSBarry Smith 6015c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6025c6c1daeSBarry Smith @*/ 6035c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6045c6c1daeSBarry Smith { 6055c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6065c6c1daeSBarry Smith 6075c6c1daeSBarry Smith PetscFunctionBegin; 6085c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6095c6c1daeSBarry Smith hdf5->timestep = timestep; 6105c6c1daeSBarry Smith PetscFunctionReturn(0); 6115c6c1daeSBarry Smith } 6125c6c1daeSBarry Smith 6133ef9c667SSatish Balay /*@ 6145c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 6155c6c1daeSBarry Smith 6165c6c1daeSBarry Smith Not collective 6175c6c1daeSBarry Smith 6185c6c1daeSBarry Smith Input Parameter: 6195c6c1daeSBarry Smith . viewer - the PetscViewer 6205c6c1daeSBarry Smith 6215c6c1daeSBarry Smith Output Parameter: 6225c6c1daeSBarry Smith . timestep - The timestep number 6235c6c1daeSBarry Smith 6245c6c1daeSBarry Smith Level: intermediate 6255c6c1daeSBarry Smith 6265c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 6275c6c1daeSBarry Smith @*/ 6285c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 6295c6c1daeSBarry Smith { 6305c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6315c6c1daeSBarry Smith 6325c6c1daeSBarry Smith PetscFunctionBegin; 6335c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6345c6c1daeSBarry Smith PetscValidPointer(timestep,2); 6355c6c1daeSBarry Smith *timestep = hdf5->timestep; 6365c6c1daeSBarry Smith PetscFunctionReturn(0); 6375c6c1daeSBarry Smith } 6385c6c1daeSBarry Smith 63936bce6e8SMatthew G. Knepley /*@C 64036bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 64136bce6e8SMatthew G. Knepley 64236bce6e8SMatthew G. Knepley Not collective 64336bce6e8SMatthew G. Knepley 64436bce6e8SMatthew G. Knepley Input Parameter: 64536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 64636bce6e8SMatthew G. Knepley 64736bce6e8SMatthew G. Knepley Output Parameter: 64836bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 64936bce6e8SMatthew G. Knepley 65036bce6e8SMatthew G. Knepley Level: advanced 65136bce6e8SMatthew G. Knepley 65236bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 65336bce6e8SMatthew G. Knepley @*/ 65436bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 65536bce6e8SMatthew G. Knepley { 65636bce6e8SMatthew G. Knepley PetscFunctionBegin; 65736bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 65836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 65936bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 66036bce6e8SMatthew G. Knepley #else 66136bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 66236bce6e8SMatthew G. Knepley #endif 66336bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 66436bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 66536bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 666c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 667de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 66836bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 66936bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 67036bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 6717e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 67236bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 67336bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 67436bce6e8SMatthew G. Knepley } 67536bce6e8SMatthew G. Knepley 67636bce6e8SMatthew G. Knepley /*@C 67736bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 67836bce6e8SMatthew G. Knepley 67936bce6e8SMatthew G. Knepley Not collective 68036bce6e8SMatthew G. Knepley 68136bce6e8SMatthew G. Knepley Input Parameter: 68236bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 68336bce6e8SMatthew G. Knepley 68436bce6e8SMatthew G. Knepley Output Parameter: 68536bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 68636bce6e8SMatthew G. Knepley 68736bce6e8SMatthew G. Knepley Level: advanced 68836bce6e8SMatthew G. Knepley 68936bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 69036bce6e8SMatthew G. Knepley @*/ 69136bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 69236bce6e8SMatthew G. Knepley { 69336bce6e8SMatthew G. Knepley PetscFunctionBegin; 69436bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 69536bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 69636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 69736bce6e8SMatthew G. Knepley #else 69836bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 69936bce6e8SMatthew G. Knepley #endif 70036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 70136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 70236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 70336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 70436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 70536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7067e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 70736bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 70836bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 70936bce6e8SMatthew G. Knepley } 71036bce6e8SMatthew G. Knepley 711df863907SAlex Fikl /*@C 712b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 71336bce6e8SMatthew G. Knepley 71436bce6e8SMatthew G. Knepley Input Parameters: 71536bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 71657d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 71736bce6e8SMatthew G. Knepley . name - The attribute name 71836bce6e8SMatthew G. Knepley . datatype - The attribute type 71936bce6e8SMatthew G. Knepley - value - The attribute value 72036bce6e8SMatthew G. Knepley 72136bce6e8SMatthew G. Knepley Level: advanced 72236bce6e8SMatthew G. Knepley 723577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 72436bce6e8SMatthew G. Knepley @*/ 72557d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 72636bce6e8SMatthew G. Knepley { 72757d22f7dSVaclav Hapla char *parent; 72860568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 729080f144cSVaclav Hapla PetscBool has; 73036bce6e8SMatthew G. Knepley PetscErrorCode ierr; 73136bce6e8SMatthew G. Knepley 73236bce6e8SMatthew G. Knepley PetscFunctionBegin; 7335cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 73457d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 735c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 736b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 73757d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 738bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 739b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 74036bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 7417e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 7427e97c476SMatthew G. Knepley size_t len; 7437e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 744729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 7457e97c476SMatthew G. Knepley } 74636bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 747729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 74860568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 749080f144cSVaclav Hapla if (has) { 750080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 751080f144cSVaclav Hapla } else { 75260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 753080f144cSVaclav Hapla } 754729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 755729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 756729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 75760568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 758729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 75957d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 76036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 76136bce6e8SMatthew G. Knepley } 76236bce6e8SMatthew G. Knepley 763df863907SAlex Fikl /*@C 764577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 765577e0e04SVaclav Hapla 766577e0e04SVaclav Hapla Input Parameters: 767577e0e04SVaclav Hapla + viewer - The HDF5 viewer 768577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 769577e0e04SVaclav Hapla . name - The attribute name 770577e0e04SVaclav Hapla . datatype - The attribute type 771577e0e04SVaclav Hapla - value - The attribute value 772577e0e04SVaclav Hapla 773577e0e04SVaclav Hapla Notes: 774577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 775577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 776577e0e04SVaclav Hapla 777577e0e04SVaclav Hapla Level: advanced 778577e0e04SVaclav Hapla 779577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 780577e0e04SVaclav Hapla @*/ 781577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 782577e0e04SVaclav Hapla { 783577e0e04SVaclav Hapla PetscErrorCode ierr; 784577e0e04SVaclav Hapla 785577e0e04SVaclav Hapla PetscFunctionBegin; 786577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 787577e0e04SVaclav Hapla PetscValidHeader(obj,2); 788577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 789b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 790577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 791577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 792577e0e04SVaclav Hapla PetscFunctionReturn(0); 793577e0e04SVaclav Hapla } 794577e0e04SVaclav Hapla 795577e0e04SVaclav Hapla /*@C 796b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 797ad0c61feSMatthew G. Knepley 798ad0c61feSMatthew G. Knepley Input Parameters: 799ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 80057d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 801ad0c61feSMatthew G. Knepley . name - The attribute name 802ad0c61feSMatthew G. Knepley - datatype - The attribute type 803ad0c61feSMatthew G. Knepley 804ad0c61feSMatthew G. Knepley Output Parameter: 805ad0c61feSMatthew G. Knepley . value - The attribute value 806ad0c61feSMatthew G. Knepley 807ad0c61feSMatthew G. Knepley Level: advanced 808ad0c61feSMatthew G. Knepley 809577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 810ad0c61feSMatthew G. Knepley @*/ 81157d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 812ad0c61feSMatthew G. Knepley { 81357d22f7dSVaclav Hapla char *parent; 814f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 815969748fdSVaclav Hapla PetscBool has; 816ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 817ad0c61feSMatthew G. Knepley 818ad0c61feSMatthew G. Knepley PetscFunctionBegin; 8195cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 82057d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 821c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 822b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 82357d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 824969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 825969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 826969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 827ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 828ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 82960568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 83060568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 831f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 832f0b58479SMatthew G. Knepley size_t len; 833e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 83415b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 835f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 836f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 837f0b58479SMatthew G. Knepley } 83870efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 839729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 840e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 841e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 84257d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 843ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 844ad0c61feSMatthew G. Knepley } 845ad0c61feSMatthew G. Knepley 846577e0e04SVaclav Hapla /*@C 847577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 848577e0e04SVaclav Hapla 849577e0e04SVaclav Hapla Input Parameters: 850577e0e04SVaclav Hapla + viewer - The HDF5 viewer 851577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 852577e0e04SVaclav Hapla . name - The attribute name 853577e0e04SVaclav Hapla - datatype - The attribute type 854577e0e04SVaclav Hapla 855577e0e04SVaclav Hapla Output Parameter: 856577e0e04SVaclav Hapla . value - The attribute value 857577e0e04SVaclav Hapla 858577e0e04SVaclav Hapla Notes: 859577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 860577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 861577e0e04SVaclav Hapla 862577e0e04SVaclav Hapla Level: advanced 863577e0e04SVaclav Hapla 864577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 865577e0e04SVaclav Hapla @*/ 866577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 867577e0e04SVaclav Hapla { 868577e0e04SVaclav Hapla PetscErrorCode ierr; 869577e0e04SVaclav Hapla 870577e0e04SVaclav Hapla PetscFunctionBegin; 871577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 872577e0e04SVaclav Hapla PetscValidHeader(obj,2); 873577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 874b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 875577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 876577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 877577e0e04SVaclav Hapla PetscFunctionReturn(0); 878577e0e04SVaclav Hapla } 879577e0e04SVaclav Hapla 880a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 881a75c3fd4SVaclav Hapla { 882a75c3fd4SVaclav Hapla htri_t exists; 883a75c3fd4SVaclav Hapla hid_t group; 884a75c3fd4SVaclav Hapla 885a75c3fd4SVaclav Hapla PetscFunctionBegin; 886a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 887a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 888a75c3fd4SVaclav Hapla if (!exists && createGroup) { 889a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 890a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 891a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 892a75c3fd4SVaclav Hapla } 893a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 894a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 895a75c3fd4SVaclav Hapla } 896a75c3fd4SVaclav Hapla 897bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 8985cdeb108SMatthew G. Knepley { 89990e03537SVaclav Hapla const char rootGroupName[] = "/"; 9005cdeb108SMatthew G. Knepley hid_t h5; 901e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 90215b861d2SVaclav Hapla PetscInt i; 90315b861d2SVaclav Hapla int n; 90485688ae2SVaclav Hapla char **hierarchy; 90585688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 9065cdeb108SMatthew G. Knepley PetscErrorCode ierr; 9075cdeb108SMatthew G. Knepley 9085cdeb108SMatthew G. Knepley PetscFunctionBegin; 9095cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 91090e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 91190e03537SVaclav Hapla else name = rootGroupName; 912ccd66a83SVaclav Hapla if (has) { 91356cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 9145cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 915ccd66a83SVaclav Hapla } 916ccd66a83SVaclav Hapla if (otype) { 917ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 91856cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 919ccd66a83SVaclav Hapla } 920c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 92185688ae2SVaclav Hapla 92285688ae2SVaclav Hapla /* 92385688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 92485688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 92585688ae2SVaclav Hapla 1) whether it's a valid link 92685688ae2SVaclav Hapla 2) whether this link resolves to an object 92785688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 92885688ae2SVaclav Hapla */ 92985688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 93085688ae2SVaclav Hapla if (!n) { 93185688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 932ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 933ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 93485688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 93585688ae2SVaclav Hapla PetscFunctionReturn(0); 93685688ae2SVaclav Hapla } 93785688ae2SVaclav Hapla for (i=0; i<n; i++) { 93885688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 93985688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 940a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 941a75c3fd4SVaclav Hapla if (!exists) break; 94285688ae2SVaclav Hapla } 94385688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 94485688ae2SVaclav Hapla 94585688ae2SVaclav Hapla /* If the object exists, get its type */ 946ccd66a83SVaclav Hapla if (exists && otype) { 9475cdeb108SMatthew G. Knepley H5O_info_t info; 9485cdeb108SMatthew G. Knepley 94976276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 95004633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 95156cc0592SVaclav Hapla *otype = info.type; 9525cdeb108SMatthew G. Knepley } 953ccd66a83SVaclav Hapla if (has) *has = exists; 9545cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 9555cdeb108SMatthew G. Knepley } 9565cdeb108SMatthew G. Knepley 957ecce1506SVaclav Hapla /*@ 9580a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 9590a9f38efSVaclav Hapla 9600a9f38efSVaclav Hapla Input Parameters: 9610a9f38efSVaclav Hapla . viewer - The HDF5 viewer 9620a9f38efSVaclav Hapla 9630a9f38efSVaclav Hapla Output Parameter: 9640a9f38efSVaclav Hapla . has - Flag for group existence 9650a9f38efSVaclav Hapla 9660a9f38efSVaclav Hapla Notes: 9670a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 9680a9f38efSVaclav Hapla 9690a9f38efSVaclav Hapla Level: advanced 9700a9f38efSVaclav Hapla 9710a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 9720a9f38efSVaclav Hapla @*/ 9730a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 9740a9f38efSVaclav Hapla { 9750a9f38efSVaclav Hapla H5O_type_t type; 9760a9f38efSVaclav Hapla const char *name; 9770a9f38efSVaclav Hapla PetscErrorCode ierr; 9780a9f38efSVaclav Hapla 9790a9f38efSVaclav Hapla PetscFunctionBegin; 9800a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9810a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 9820a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 9830a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 9840a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 9850a9f38efSVaclav Hapla PetscFunctionReturn(0); 9860a9f38efSVaclav Hapla } 9870a9f38efSVaclav Hapla 9880a9f38efSVaclav Hapla /*@ 989e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 990ecce1506SVaclav Hapla 991ecce1506SVaclav Hapla Input Parameters: 992ecce1506SVaclav Hapla + viewer - The HDF5 viewer 993ecce1506SVaclav Hapla - obj - The named object 994ecce1506SVaclav Hapla 995ecce1506SVaclav Hapla Output Parameter: 996ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 997ecce1506SVaclav Hapla 998e3f143f7SVaclav Hapla Notes: 999e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1000e3f143f7SVaclav Hapla 1001ecce1506SVaclav Hapla Level: advanced 1002ecce1506SVaclav Hapla 1003e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1004ecce1506SVaclav Hapla @*/ 1005ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1006ecce1506SVaclav Hapla { 100756cc0592SVaclav Hapla H5O_type_t type; 10086c132bc1SVaclav Hapla char *path; 1009ecce1506SVaclav Hapla PetscErrorCode ierr; 1010ecce1506SVaclav Hapla 1011ecce1506SVaclav Hapla PetscFunctionBegin; 1012c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1013c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1014c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1015ecce1506SVaclav Hapla *has = PETSC_FALSE; 1016e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 10176c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 10186c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 101956cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 10206c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1021ecce1506SVaclav Hapla PetscFunctionReturn(0); 1022ecce1506SVaclav Hapla } 1023ecce1506SVaclav Hapla 1024df863907SAlex Fikl /*@C 1025ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1026e7bdbf8eSMatthew G. Knepley 1027e7bdbf8eSMatthew G. Knepley Input Parameters: 1028e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 102910e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1030e7bdbf8eSMatthew G. Knepley - name - The attribute name 1031e7bdbf8eSMatthew G. Knepley 1032e7bdbf8eSMatthew G. Knepley Output Parameter: 1033e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1034e7bdbf8eSMatthew G. Knepley 1035e7bdbf8eSMatthew G. Knepley Level: advanced 1036e7bdbf8eSMatthew G. Knepley 1037577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1038e7bdbf8eSMatthew G. Knepley @*/ 103910e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1040e7bdbf8eSMatthew G. Knepley { 104110e69e8fSVaclav Hapla char *parent; 1042e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1043e7bdbf8eSMatthew G. Knepley 1044e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 10455cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 104610e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1047c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1048c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 104910e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1050bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 105110e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 105210e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 105306db490cSVaclav Hapla PetscFunctionReturn(0); 105406db490cSVaclav Hapla } 105506db490cSVaclav Hapla 1056577e0e04SVaclav Hapla /*@C 1057577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1058577e0e04SVaclav Hapla 1059577e0e04SVaclav Hapla Input Parameters: 1060577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1061577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1062577e0e04SVaclav Hapla - name - The attribute name 1063577e0e04SVaclav Hapla 1064577e0e04SVaclav Hapla Output Parameter: 1065577e0e04SVaclav Hapla . has - Flag for attribute existence 1066577e0e04SVaclav Hapla 1067577e0e04SVaclav Hapla Notes: 1068577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1069577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1070577e0e04SVaclav Hapla 1071577e0e04SVaclav Hapla Level: advanced 1072577e0e04SVaclav Hapla 1073577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1074577e0e04SVaclav Hapla @*/ 1075577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1076577e0e04SVaclav Hapla { 1077577e0e04SVaclav Hapla PetscErrorCode ierr; 1078577e0e04SVaclav Hapla 1079577e0e04SVaclav Hapla PetscFunctionBegin; 1080577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1081577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1082577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1083577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1084577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1085577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1086577e0e04SVaclav Hapla PetscFunctionReturn(0); 1087577e0e04SVaclav Hapla } 1088577e0e04SVaclav Hapla 108906db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 109006db490cSVaclav Hapla { 109106db490cSVaclav Hapla hid_t h5; 109206db490cSVaclav Hapla htri_t hhas; 109306db490cSVaclav Hapla PetscErrorCode ierr; 109406db490cSVaclav Hapla 109506db490cSVaclav Hapla PetscFunctionBegin; 109606db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 10972f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 10985cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1099e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1100e7bdbf8eSMatthew G. Knepley } 1101e7bdbf8eSMatthew G. Knepley 1102eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1103a5e1feadSVaclav Hapla { 1104fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1105a5e1feadSVaclav Hapla PetscErrorCode ierr; 1106a5e1feadSVaclav Hapla 1107a5e1feadSVaclav Hapla PetscFunctionBegin; 110869a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 110969a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 111069a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 111169a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 11129568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 1113adb363a2SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr); 1114aa54fc5aSVaclav Hapla if (h->complexVal) {ierr = PetscViewerHDF5ReadAttribute(viewer,name,"complex",PETSC_BOOL,&h->complexVal);CHKERRQ(ierr);} 111509dabeb0SVaclav Hapla /* MATLAB stores column vectors horizontally */ 1116adb363a2SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1117b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 1118b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1119b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1120b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1121b8ef406cSVaclav Hapla #endif 1122b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 112369a06e7bSVaclav Hapla *ctx = h; 112469a06e7bSVaclav Hapla PetscFunctionReturn(0); 112569a06e7bSVaclav Hapla } 112669a06e7bSVaclav Hapla 1127eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 112869a06e7bSVaclav Hapla { 112969a06e7bSVaclav Hapla HDF5ReadCtx h; 113069a06e7bSVaclav Hapla PetscErrorCode ierr; 113169a06e7bSVaclav Hapla 113269a06e7bSVaclav Hapla PetscFunctionBegin; 113369a06e7bSVaclav Hapla h = *ctx; 1134b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 113590e03537SVaclav Hapla PetscStackCallHDF5(H5Gclose,(h->group)); 113669a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 113769a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 113869a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 113969a06e7bSVaclav Hapla PetscFunctionReturn(0); 114069a06e7bSVaclav Hapla } 114169a06e7bSVaclav Hapla 1142eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 114369a06e7bSVaclav Hapla { 114469a06e7bSVaclav Hapla int rdim, dim; 114569a06e7bSVaclav Hapla hsize_t dims[4]; 114609dabeb0SVaclav Hapla PetscInt bsInd, lenInd, bs, len, N; 11478374c777SVaclav Hapla PetscLayout map; 11488374c777SVaclav Hapla PetscErrorCode ierr; 114969a06e7bSVaclav Hapla 115069a06e7bSVaclav Hapla PetscFunctionBegin; 11518374c777SVaclav Hapla if (!(*map_)) { 11528374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 11538374c777SVaclav Hapla } 11548374c777SVaclav Hapla map = *map_; 11559787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1156a5e1feadSVaclav Hapla dim = 0; 11579568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1158a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 11599568d583SVaclav Hapla if (ctx->complexVal) ++dim; 11609787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 116169a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 11629787f08bSVaclav Hapla /* calculate expected dimension indices */ 11639787f08bSVaclav Hapla lenInd = 0; 11649568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 11659787f08bSVaclav Hapla bsInd = lenInd + 1; 11669568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 11679787f08bSVaclav Hapla if (rdim == dim) { 116845e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 11699787f08bSVaclav Hapla } else if (rdim == dim+1) { 117045e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 11719568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1172a5e1feadSVaclav Hapla } else { 1173b54f1188SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim); 1174a5e1feadSVaclav Hapla } 117509dabeb0SVaclav Hapla len = dims[lenInd]; 117609dabeb0SVaclav Hapla if (ctx->horizontal) { 1177b54f1188SVaclav Hapla if (len != 1) SETERRQ(PetscObjectComm((PetscObject)viewer), 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."); 117809dabeb0SVaclav Hapla len = bs; 117909dabeb0SVaclav Hapla bs = 1; 118009dabeb0SVaclav Hapla } 118109dabeb0SVaclav Hapla N = (PetscInt) len*bs; 11828374c777SVaclav Hapla 11838374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 11848374c777SVaclav Hapla if (map->bs < 0) { 118545e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 118645e21b7fSVaclav 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); 11878374c777SVaclav Hapla if (map->N < 0) { 118845e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 1189b54f1188SVaclav Hapla } else if (map->N != N) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_FILE_UNEXPECTED, "Global size of array in file is %D, not %D as expected",N,map->N); 119069a06e7bSVaclav Hapla PetscFunctionReturn(0); 119169a06e7bSVaclav Hapla } 119269a06e7bSVaclav Hapla 1193eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 119416f6402dSVaclav Hapla { 119516f6402dSVaclav Hapla hsize_t count[4], offset[4]; 119616f6402dSVaclav Hapla int dim; 119716f6402dSVaclav Hapla PetscInt bs, n, low; 119816f6402dSVaclav Hapla PetscErrorCode ierr; 119916f6402dSVaclav Hapla 120016f6402dSVaclav Hapla PetscFunctionBegin; 120116f6402dSVaclav Hapla /* Compute local size and ownership range */ 120216f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 120316f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 120416f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 120516f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 120616f6402dSVaclav Hapla 120716f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 120816f6402dSVaclav Hapla dim = 0; 120916f6402dSVaclav Hapla if (ctx->timestep >= 0) { 121016f6402dSVaclav Hapla count[dim] = 1; 121116f6402dSVaclav Hapla offset[dim] = ctx->timestep; 121216f6402dSVaclav Hapla ++dim; 121316f6402dSVaclav Hapla } 121409dabeb0SVaclav Hapla if (ctx->horizontal) { 121509dabeb0SVaclav Hapla count[dim] = 1; 121609dabeb0SVaclav Hapla offset[dim] = 0; 121709dabeb0SVaclav Hapla ++dim; 121809dabeb0SVaclav Hapla } 121916f6402dSVaclav Hapla { 122016f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 122116f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 122216f6402dSVaclav Hapla ++dim; 122316f6402dSVaclav Hapla } 122416f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 122509dabeb0SVaclav Hapla if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 122616f6402dSVaclav Hapla count[dim] = bs; 122716f6402dSVaclav Hapla offset[dim] = 0; 122816f6402dSVaclav Hapla ++dim; 122916f6402dSVaclav Hapla } 123016f6402dSVaclav Hapla if (ctx->complexVal) { 123116f6402dSVaclav Hapla count[dim] = 2; 123216f6402dSVaclav Hapla offset[dim] = 0; 123316f6402dSVaclav Hapla ++dim; 123416f6402dSVaclav Hapla } 123516f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 123616f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 123716f6402dSVaclav Hapla PetscFunctionReturn(0); 123816f6402dSVaclav Hapla } 123916f6402dSVaclav Hapla 1240eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1241ef2d82ceSVaclav Hapla { 1242ef2d82ceSVaclav Hapla PetscFunctionBegin; 1243ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1244ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1245ef2d82ceSVaclav Hapla } 1246ef2d82ceSVaclav Hapla 1247dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1248a25c73c6SVaclav Hapla { 1249fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1250fbd37863SVaclav Hapla hid_t memspace=0; 1251a25c73c6SVaclav Hapla size_t unitsize; 1252a25c73c6SVaclav Hapla void *arr; 1253a25c73c6SVaclav Hapla PetscErrorCode ierr; 1254a25c73c6SVaclav Hapla 1255a25c73c6SVaclav Hapla PetscFunctionBegin; 1256eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 12575a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 12585a89fdf4SVaclav Hapla if (!h->complexVal) { 1259c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1260c76551afSVaclav 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."); 12615a89fdf4SVaclav Hapla } 12625a89fdf4SVaclav Hapla #else 12635a89fdf4SVaclav 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."); 12645a89fdf4SVaclav Hapla #endif 1265e9e90110SVaclav Hapla 1266e9e90110SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1267e9e90110SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1268e9e90110SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 1269e9e90110SVaclav Hapla 12704fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 12714fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 1272dff35581SVaclav 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); 12734fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 12744fc17bcdSVaclav Hapla 1275eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1276a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1277eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1278a25c73c6SVaclav Hapla *newarr = arr; 1279a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1280a25c73c6SVaclav Hapla } 1281a25c73c6SVaclav Hapla 1282c1aaad9cSVaclav Hapla /*@C 1283c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1284c1aaad9cSVaclav Hapla 1285c1aaad9cSVaclav Hapla Input Parameters: 1286c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1287c1aaad9cSVaclav Hapla - name - The vector name 1288c1aaad9cSVaclav Hapla 1289c1aaad9cSVaclav Hapla Output Parameter: 1290c1aaad9cSVaclav Hapla + bs - block size 1291c1aaad9cSVaclav Hapla - N - global size 1292c1aaad9cSVaclav Hapla 1293c1aaad9cSVaclav Hapla Note: 1294c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1295c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1296c1aaad9cSVaclav Hapla 1297c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1298c1aaad9cSVaclav Hapla 1299c1aaad9cSVaclav Hapla Level: advanced 1300c1aaad9cSVaclav Hapla 1301c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1302c1aaad9cSVaclav Hapla @*/ 130369a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 130469a06e7bSVaclav Hapla { 1305fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 13068374c777SVaclav Hapla PetscLayout map=NULL; 130769a06e7bSVaclav Hapla PetscErrorCode ierr; 130869a06e7bSVaclav Hapla 130969a06e7bSVaclav Hapla PetscFunctionBegin; 1310c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1311eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1312eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1313eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 13148374c777SVaclav Hapla if (bs) *bs = map->bs; 13158374c777SVaclav Hapla if (N) *N = map->N; 13168374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1317a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1318a5e1feadSVaclav Hapla } 1319a5e1feadSVaclav Hapla 1320a75e6a4aSMatthew G. Knepley /* 1321a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1322a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1323a75e6a4aSMatthew G. Knepley */ 1324d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1325a75e6a4aSMatthew G. Knepley 1326a75e6a4aSMatthew G. Knepley /*@C 1327a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1328a75e6a4aSMatthew G. Knepley 1329a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1330a75e6a4aSMatthew G. Knepley 1331a75e6a4aSMatthew G. Knepley Input Parameter: 1332a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1333a75e6a4aSMatthew G. Knepley 1334a75e6a4aSMatthew G. Knepley Level: intermediate 1335a75e6a4aSMatthew G. Knepley 1336a75e6a4aSMatthew G. Knepley Options Database Keys: 1337a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1338a75e6a4aSMatthew G. Knepley 1339a75e6a4aSMatthew G. Knepley Environmental variables: 1340a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1341a75e6a4aSMatthew G. Knepley 1342a75e6a4aSMatthew G. Knepley Notes: 1343a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1344a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1345a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1346a75e6a4aSMatthew G. Knepley 1347a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1348a75e6a4aSMatthew G. Knepley @*/ 1349a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1350a75e6a4aSMatthew G. Knepley { 1351a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1352a75e6a4aSMatthew G. Knepley PetscBool flg; 1353a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1354a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1355a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1356a75e6a4aSMatthew G. Knepley 1357a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1358a75e6a4aSMatthew 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);} 1359a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 136012801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1361a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1362a75e6a4aSMatthew G. Knepley } 136347435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1364a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1365a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1366a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1367a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1368a75e6a4aSMatthew G. Knepley if (!flg) { 1369a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1370a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1371a75e6a4aSMatthew G. Knepley } 1372a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1373a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1374a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1375a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 137647435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1377a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1378a75e6a4aSMatthew G. Knepley } 1379a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1380a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1381a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1382a75e6a4aSMatthew G. Knepley } 1383