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 51577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 52577e0e04SVaclav Hapla { 53577e0e04SVaclav Hapla PetscBool has; 54577e0e04SVaclav Hapla const char *group; 55577e0e04SVaclav Hapla PetscErrorCode ierr; 56577e0e04SVaclav Hapla 57577e0e04SVaclav Hapla PetscFunctionBegin; 58577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 59577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 60577e0e04SVaclav Hapla if (!has) { 61577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 62577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 63577e0e04SVaclav Hapla } 64577e0e04SVaclav Hapla PetscFunctionReturn(0); 65577e0e04SVaclav Hapla } 66577e0e04SVaclav Hapla 674416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 6882ea9b62SBarry Smith { 6982ea9b62SBarry Smith PetscErrorCode ierr; 7082ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 7182ea9b62SBarry Smith 7282ea9b62SBarry Smith PetscFunctionBegin; 7382ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 7482ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 759a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 7682ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7782ea9b62SBarry Smith PetscFunctionReturn(0); 7882ea9b62SBarry Smith } 7982ea9b62SBarry Smith 805c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 815c6c1daeSBarry Smith { 825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 835c6c1daeSBarry Smith PetscErrorCode ierr; 845c6c1daeSBarry Smith 855c6c1daeSBarry Smith PetscFunctionBegin; 865c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 87729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 885c6c1daeSBarry Smith PetscFunctionReturn(0); 895c6c1daeSBarry Smith } 905c6c1daeSBarry Smith 915c6c1daeSBarry Smith PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 925c6c1daeSBarry Smith { 935c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 945c6c1daeSBarry Smith PetscErrorCode ierr; 955c6c1daeSBarry Smith 965c6c1daeSBarry Smith PetscFunctionBegin; 975c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 985c6c1daeSBarry Smith while (hdf5->groups) { 995c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 1005c6c1daeSBarry Smith 1015c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 1025c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 1035c6c1daeSBarry Smith hdf5->groups = tmp; 1045c6c1daeSBarry Smith } 10561912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 10661912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 10761912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 10861912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 1095c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1100b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 111d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1120b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 113058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 114058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 115058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 116058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 1175c6c1daeSBarry Smith PetscFunctionReturn(0); 1185c6c1daeSBarry Smith } 1195c6c1daeSBarry Smith 1205c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1215c6c1daeSBarry Smith { 1225c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1235c6c1daeSBarry Smith 1245c6c1daeSBarry Smith PetscFunctionBegin; 1255c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 1265c6c1daeSBarry Smith hdf5->btype = type; 1275c6c1daeSBarry Smith PetscFunctionReturn(0); 1285c6c1daeSBarry Smith } 1295c6c1daeSBarry Smith 13082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13182ea9b62SBarry Smith { 13282ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13382ea9b62SBarry Smith 13482ea9b62SBarry Smith PetscFunctionBegin; 13582ea9b62SBarry Smith hdf5->basedimension2 = flg; 13682ea9b62SBarry Smith PetscFunctionReturn(0); 13782ea9b62SBarry Smith } 13882ea9b62SBarry Smith 139df863907SAlex Fikl /*@ 14082ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14182ea9b62SBarry Smith dimension of 2. 14282ea9b62SBarry Smith 14382ea9b62SBarry Smith Logically Collective on PetscViewer 14482ea9b62SBarry Smith 14582ea9b62SBarry Smith Input Parameters: 14682ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14782ea9b62SBarry 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 14882ea9b62SBarry Smith 14982ea9b62SBarry Smith Options Database: 15082ea9b62SBarry 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 15182ea9b62SBarry Smith 15282ea9b62SBarry Smith 15395452b02SPatrick Sanan Notes: 15495452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 15582ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15682ea9b62SBarry Smith 15782ea9b62SBarry Smith Level: intermediate 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16082ea9b62SBarry Smith 16182ea9b62SBarry Smith @*/ 16282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16382ea9b62SBarry Smith { 16482ea9b62SBarry Smith PetscErrorCode ierr; 16582ea9b62SBarry Smith 16682ea9b62SBarry Smith PetscFunctionBegin; 16782ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16882ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16982ea9b62SBarry Smith PetscFunctionReturn(0); 17082ea9b62SBarry Smith } 17182ea9b62SBarry Smith 172df863907SAlex Fikl /*@ 17382ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17482ea9b62SBarry Smith dimension of 2. 17582ea9b62SBarry Smith 17682ea9b62SBarry Smith Logically Collective on PetscViewer 17782ea9b62SBarry Smith 17882ea9b62SBarry Smith Input Parameter: 17982ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18082ea9b62SBarry Smith 18182ea9b62SBarry Smith Output Parameter: 18282ea9b62SBarry 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 18382ea9b62SBarry Smith 18495452b02SPatrick Sanan Notes: 18595452b02SPatrick 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 18682ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18782ea9b62SBarry Smith 18882ea9b62SBarry Smith Level: intermediate 18982ea9b62SBarry Smith 19082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 19182ea9b62SBarry Smith 19282ea9b62SBarry Smith @*/ 19382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19482ea9b62SBarry Smith { 19582ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19682ea9b62SBarry Smith 19782ea9b62SBarry Smith PetscFunctionBegin; 19882ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19982ea9b62SBarry Smith *flg = hdf5->basedimension2; 20082ea9b62SBarry Smith PetscFunctionReturn(0); 20182ea9b62SBarry Smith } 20282ea9b62SBarry Smith 2039a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2049a0502c6SHåkon Strandenes { 2059a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2069a0502c6SHåkon Strandenes 2079a0502c6SHåkon Strandenes PetscFunctionBegin; 2089a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2099a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2109a0502c6SHåkon Strandenes } 2119a0502c6SHåkon Strandenes 212df863907SAlex Fikl /*@ 2139a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2149a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2159a0502c6SHåkon Strandenes 2169a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2179a0502c6SHåkon Strandenes 2189a0502c6SHåkon Strandenes Input Parameters: 2199a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2209a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2219a0502c6SHåkon Strandenes 2229a0502c6SHåkon Strandenes Options Database: 2239a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2249a0502c6SHåkon Strandenes 2259a0502c6SHåkon Strandenes 22695452b02SPatrick Sanan Notes: 22795452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2289a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes Level: intermediate 2319a0502c6SHåkon Strandenes 2329a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2339a0502c6SHåkon Strandenes PetscReal 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes @*/ 2369a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2379a0502c6SHåkon Strandenes { 2389a0502c6SHåkon Strandenes PetscErrorCode ierr; 2399a0502c6SHåkon Strandenes 2409a0502c6SHåkon Strandenes PetscFunctionBegin; 2419a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2429a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2439a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2449a0502c6SHåkon Strandenes } 2459a0502c6SHåkon Strandenes 246df863907SAlex Fikl /*@ 2479a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2489a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2499a0502c6SHåkon Strandenes 2509a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2519a0502c6SHåkon Strandenes 2529a0502c6SHåkon Strandenes Input Parameter: 2539a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2549a0502c6SHåkon Strandenes 2559a0502c6SHåkon Strandenes Output Parameter: 2569a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2579a0502c6SHåkon Strandenes 25895452b02SPatrick Sanan Notes: 25995452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2609a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2619a0502c6SHåkon Strandenes 2629a0502c6SHåkon Strandenes Level: intermediate 2639a0502c6SHåkon Strandenes 2649a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2659a0502c6SHåkon Strandenes PetscReal 2669a0502c6SHåkon Strandenes 2679a0502c6SHåkon Strandenes @*/ 2689a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2699a0502c6SHåkon Strandenes { 2709a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2719a0502c6SHåkon Strandenes 2729a0502c6SHåkon Strandenes PetscFunctionBegin; 2739a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2749a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2759a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2769a0502c6SHåkon Strandenes } 2779a0502c6SHåkon Strandenes 2785c6c1daeSBarry Smith PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2795c6c1daeSBarry Smith { 2805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2815c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2825c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2835c6c1daeSBarry Smith #endif 2845c6c1daeSBarry Smith hid_t plist_id; 2855c6c1daeSBarry Smith PetscErrorCode ierr; 2865c6c1daeSBarry Smith 2875c6c1daeSBarry Smith PetscFunctionBegin; 288f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 289f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2905c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2915c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 292729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2935c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 294729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2955c6c1daeSBarry Smith #endif 2965c6c1daeSBarry Smith /* Create or open the file collectively */ 2975c6c1daeSBarry Smith switch (hdf5->btype) { 2985c6c1daeSBarry Smith case FILE_MODE_READ: 299729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3005c6c1daeSBarry Smith break; 3015c6c1daeSBarry Smith case FILE_MODE_APPEND: 302729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 3035c6c1daeSBarry Smith break; 3045c6c1daeSBarry Smith case FILE_MODE_WRITE: 305729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 3065c6c1daeSBarry Smith break; 3075c6c1daeSBarry Smith default: 3085c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3095c6c1daeSBarry Smith } 3105c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 311729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 3125c6c1daeSBarry Smith PetscFunctionReturn(0); 3135c6c1daeSBarry Smith } 3145c6c1daeSBarry Smith 315d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 316d1232d7fSToby Isaac { 317d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 318d1232d7fSToby Isaac 319d1232d7fSToby Isaac PetscFunctionBegin; 320d1232d7fSToby Isaac *name = vhdf5->filename; 321d1232d7fSToby Isaac PetscFunctionReturn(0); 322d1232d7fSToby Isaac } 323d1232d7fSToby Isaac 324058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 325058bd781SVaclav Hapla { 326058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 327058bd781SVaclav Hapla PetscErrorCode ierr; 328058bd781SVaclav Hapla 329058bd781SVaclav Hapla PetscFunctionBegin; 330058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 331058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 332058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 333058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 334058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 335058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 336058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 337058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 338058bd781SVaclav Hapla PetscFunctionReturn(0); 339058bd781SVaclav Hapla } 340058bd781SVaclav Hapla 341058bd781SVaclav Hapla /*@C 342058bd781SVaclav 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. 343058bd781SVaclav Hapla 344058bd781SVaclav Hapla Collective on PetscViewer 345058bd781SVaclav Hapla 346058bd781SVaclav Hapla Input Parameters: 347058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 348058bd781SVaclav 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 349058bd781SVaclav Hapla . jname - name of dataset j representing column indices 350058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 351058bd781SVaclav Hapla - cname - name of attribute stoting column count 352058bd781SVaclav Hapla 353058bd781SVaclav Hapla Level: advanced 354058bd781SVaclav Hapla 355058bd781SVaclav Hapla Notes: 356058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 357058bd781SVaclav Hapla 358058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 359058bd781SVaclav Hapla @*/ 360058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 361058bd781SVaclav Hapla { 362058bd781SVaclav Hapla PetscErrorCode ierr; 363058bd781SVaclav Hapla 364058bd781SVaclav Hapla PetscFunctionBegin; 365058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 366058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 367058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 368058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 369058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 370058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 371058bd781SVaclav Hapla PetscFunctionReturn(0); 372058bd781SVaclav Hapla } 373058bd781SVaclav Hapla 374058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 375058bd781SVaclav Hapla { 376058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 377058bd781SVaclav Hapla 378058bd781SVaclav Hapla PetscFunctionBegin; 379058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 380058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 381058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 382058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 383058bd781SVaclav Hapla PetscFunctionReturn(0); 384058bd781SVaclav Hapla } 385058bd781SVaclav Hapla 386058bd781SVaclav Hapla /*@C 387058bd781SVaclav 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. 388058bd781SVaclav Hapla 389058bd781SVaclav Hapla Collective on PetscViewer 390058bd781SVaclav Hapla 391058bd781SVaclav Hapla Input Parameters: 392058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 393058bd781SVaclav Hapla 394058bd781SVaclav Hapla Output Parameters: 395058bd781SVaclav 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 396058bd781SVaclav Hapla . jname - name of dataset j representing column indices 397058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 398058bd781SVaclav Hapla - cname - name of attribute stoting column count 399058bd781SVaclav Hapla 400058bd781SVaclav Hapla Level: advanced 401058bd781SVaclav Hapla 402058bd781SVaclav Hapla Notes: 403058bd781SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 404058bd781SVaclav Hapla 405058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 406058bd781SVaclav Hapla @*/ 407058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 408058bd781SVaclav Hapla { 409058bd781SVaclav Hapla PetscErrorCode ierr; 410058bd781SVaclav Hapla 411058bd781SVaclav Hapla PetscFunctionBegin; 412058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 413058bd781SVaclav Hapla PetscValidPointer(iname,2); 414058bd781SVaclav Hapla PetscValidPointer(jname,3); 415058bd781SVaclav Hapla PetscValidPointer(aname,4); 416058bd781SVaclav Hapla PetscValidPointer(cname,5); 417058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 418058bd781SVaclav Hapla PetscFunctionReturn(0); 419058bd781SVaclav Hapla } 420058bd781SVaclav Hapla 4218556b5ebSBarry Smith /*MC 4228556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4238556b5ebSBarry Smith 4248556b5ebSBarry Smith 4258556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4268556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4278556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4288556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4298556b5ebSBarry Smith 4301b266c99SBarry Smith Level: beginner 4318556b5ebSBarry Smith M*/ 432d1232d7fSToby Isaac 4338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4345c6c1daeSBarry Smith { 4355c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4365c6c1daeSBarry Smith PetscErrorCode ierr; 4375c6c1daeSBarry Smith 4385c6c1daeSBarry Smith PetscFunctionBegin; 439b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4405c6c1daeSBarry Smith 4415c6c1daeSBarry Smith v->data = (void*) hdf5; 4425c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 44382ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 4445c6c1daeSBarry Smith v->ops->flush = 0; 4455c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4465c6c1daeSBarry Smith hdf5->filename = 0; 4475c6c1daeSBarry Smith hdf5->timestep = -1; 4480298fd71SBarry Smith hdf5->groups = NULL; 4495c6c1daeSBarry Smith 450058bd781SVaclav Hapla /* ir and jc are deliberately swapped as MATLAB uses column-major format */ 451058bd781SVaclav Hapla ierr = PetscStrallocpy("jc", &hdf5->mataij_iname);CHKERRQ(ierr); 452058bd781SVaclav Hapla ierr = PetscStrallocpy("ir", &hdf5->mataij_jname);CHKERRQ(ierr); 453058bd781SVaclav Hapla ierr = PetscStrallocpy("data",&hdf5->mataij_aname);CHKERRQ(ierr); 454058bd781SVaclav Hapla ierr = PetscStrallocpy("MATLAB_sparse", &hdf5->mataij_cname);CHKERRQ(ierr); 455058bd781SVaclav Hapla 4560b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 457d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4580b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 45982ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4609a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 461058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 462058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4635c6c1daeSBarry Smith PetscFunctionReturn(0); 4645c6c1daeSBarry Smith } 4655c6c1daeSBarry Smith 4665c6c1daeSBarry Smith /*@C 4675c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4685c6c1daeSBarry Smith 4695c6c1daeSBarry Smith Collective on MPI_Comm 4705c6c1daeSBarry Smith 4715c6c1daeSBarry Smith Input Parameters: 4725c6c1daeSBarry Smith + comm - MPI communicator 4735c6c1daeSBarry Smith . name - name of file 4745c6c1daeSBarry Smith - type - type of file 4755c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4765c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4775c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4785c6c1daeSBarry Smith 4795c6c1daeSBarry Smith Output Parameter: 4805c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4815c6c1daeSBarry Smith 48282ea9b62SBarry Smith Options Database: 48382ea9b62SBarry 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 4849a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 48582ea9b62SBarry Smith 4865c6c1daeSBarry Smith Level: beginner 4875c6c1daeSBarry Smith 4885c6c1daeSBarry Smith Note: 4895c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4905c6c1daeSBarry Smith 4915c6c1daeSBarry Smith Concepts: HDF5 files 4925c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4935c6c1daeSBarry Smith 4946a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4959a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 496a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 4975c6c1daeSBarry Smith @*/ 4985c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 4995c6c1daeSBarry Smith { 5005c6c1daeSBarry Smith PetscErrorCode ierr; 5015c6c1daeSBarry Smith 5025c6c1daeSBarry Smith PetscFunctionBegin; 5035c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5045c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5055c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5065c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 5075c6c1daeSBarry Smith PetscFunctionReturn(0); 5085c6c1daeSBarry Smith } 5095c6c1daeSBarry Smith 5105c6c1daeSBarry Smith /*@C 5115c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5125c6c1daeSBarry Smith 5135c6c1daeSBarry Smith Not collective 5145c6c1daeSBarry Smith 5155c6c1daeSBarry Smith Input Parameter: 5165c6c1daeSBarry Smith . viewer - the PetscViewer 5175c6c1daeSBarry Smith 5185c6c1daeSBarry Smith Output Parameter: 5195c6c1daeSBarry Smith . file_id - The file id 5205c6c1daeSBarry Smith 5215c6c1daeSBarry Smith Level: intermediate 5225c6c1daeSBarry Smith 5235c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5245c6c1daeSBarry Smith @*/ 5255c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5265c6c1daeSBarry Smith { 5275c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5285c6c1daeSBarry Smith 5295c6c1daeSBarry Smith PetscFunctionBegin; 5305c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5315c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5325c6c1daeSBarry Smith PetscFunctionReturn(0); 5335c6c1daeSBarry Smith } 5345c6c1daeSBarry Smith 5355c6c1daeSBarry Smith /*@C 5365c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5375c6c1daeSBarry Smith 5385c6c1daeSBarry Smith Not collective 5395c6c1daeSBarry Smith 5405c6c1daeSBarry Smith Input Parameters: 5415c6c1daeSBarry Smith + viewer - the PetscViewer 5425c6c1daeSBarry Smith - name - The group name 5435c6c1daeSBarry Smith 5445c6c1daeSBarry Smith Level: intermediate 5455c6c1daeSBarry Smith 546874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5475c6c1daeSBarry Smith @*/ 5485c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5495c6c1daeSBarry Smith { 5505c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5515c6c1daeSBarry Smith GroupList *groupNode; 5525c6c1daeSBarry Smith PetscErrorCode ierr; 5535c6c1daeSBarry Smith 5545c6c1daeSBarry Smith PetscFunctionBegin; 5555c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5565c6c1daeSBarry Smith PetscValidCharPointer(name,2); 55795dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5585c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 559a297a907SKarl Rupp 5605c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5615c6c1daeSBarry Smith hdf5->groups = groupNode; 5625c6c1daeSBarry Smith PetscFunctionReturn(0); 5635c6c1daeSBarry Smith } 5645c6c1daeSBarry Smith 5653ef9c667SSatish Balay /*@ 5665c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5675c6c1daeSBarry Smith 5685c6c1daeSBarry Smith Not collective 5695c6c1daeSBarry Smith 5705c6c1daeSBarry Smith Input Parameter: 5715c6c1daeSBarry Smith . viewer - the PetscViewer 5725c6c1daeSBarry Smith 5735c6c1daeSBarry Smith Level: intermediate 5745c6c1daeSBarry Smith 575874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5765c6c1daeSBarry Smith @*/ 5775c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5785c6c1daeSBarry Smith { 5795c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5805c6c1daeSBarry Smith GroupList *groupNode; 5815c6c1daeSBarry Smith PetscErrorCode ierr; 5825c6c1daeSBarry Smith 5835c6c1daeSBarry Smith PetscFunctionBegin; 5845c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 58582f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5865c6c1daeSBarry Smith groupNode = hdf5->groups; 5875c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5885c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5895c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5905c6c1daeSBarry Smith PetscFunctionReturn(0); 5915c6c1daeSBarry Smith } 5925c6c1daeSBarry Smith 5935c6c1daeSBarry Smith /*@C 594874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 595874270d9SVaclav Hapla If none has been assigned, returns NULL. 5965c6c1daeSBarry Smith 5975c6c1daeSBarry Smith Not collective 5985c6c1daeSBarry Smith 5995c6c1daeSBarry Smith Input Parameter: 6005c6c1daeSBarry Smith . viewer - the PetscViewer 6015c6c1daeSBarry Smith 6025c6c1daeSBarry Smith Output Parameter: 6035c6c1daeSBarry Smith . name - The group name 6045c6c1daeSBarry Smith 6055c6c1daeSBarry Smith Level: intermediate 6065c6c1daeSBarry Smith 607874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6085c6c1daeSBarry Smith @*/ 6095c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 6105c6c1daeSBarry Smith { 6115c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6125c6c1daeSBarry Smith 6135c6c1daeSBarry Smith PetscFunctionBegin; 6145c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 615c959eef4SJed Brown PetscValidPointer(name,2); 616a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6170298fd71SBarry Smith else *name = NULL; 6185c6c1daeSBarry Smith PetscFunctionReturn(0); 6195c6c1daeSBarry Smith } 6205c6c1daeSBarry Smith 6218c557b5aSVaclav Hapla /*@ 622874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 623874270d9SVaclav Hapla and return this group's ID and file ID. 624874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 625874270d9SVaclav Hapla 626874270d9SVaclav Hapla Not collective 627874270d9SVaclav Hapla 628874270d9SVaclav Hapla Input Parameter: 629874270d9SVaclav Hapla . viewer - the PetscViewer 630874270d9SVaclav Hapla 631874270d9SVaclav Hapla Output Parameter: 632874270d9SVaclav Hapla + fileId - The HDF5 file ID 633874270d9SVaclav Hapla - groupId - The HDF5 group ID 634874270d9SVaclav Hapla 635874270d9SVaclav Hapla Level: intermediate 636874270d9SVaclav Hapla 637874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 638874270d9SVaclav Hapla @*/ 63954dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 64054dbf706SVaclav Hapla { 641*90e03537SVaclav Hapla hid_t file_id; 642*90e03537SVaclav Hapla H5O_type_t type; 64354dbf706SVaclav Hapla const char *groupName = NULL; 64454dbf706SVaclav Hapla PetscErrorCode ierr; 64554dbf706SVaclav Hapla 64654dbf706SVaclav Hapla PetscFunctionBegin; 64754dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 64854dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 649*90e03537SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, PETSC_TRUE, NULL, &type);CHKERRQ(ierr); 650*90e03537SVaclav 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); 651*90e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 65254dbf706SVaclav Hapla *fileId = file_id; 65354dbf706SVaclav Hapla PetscFunctionReturn(0); 65454dbf706SVaclav Hapla } 65554dbf706SVaclav Hapla 6563ef9c667SSatish Balay /*@ 6575c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith Not collective 6605c6c1daeSBarry Smith 6615c6c1daeSBarry Smith Input Parameter: 6625c6c1daeSBarry Smith . viewer - the PetscViewer 6635c6c1daeSBarry Smith 6645c6c1daeSBarry Smith Level: intermediate 6655c6c1daeSBarry Smith 6665c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6675c6c1daeSBarry Smith @*/ 6685c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6695c6c1daeSBarry Smith { 6705c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6715c6c1daeSBarry Smith 6725c6c1daeSBarry Smith PetscFunctionBegin; 6735c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6745c6c1daeSBarry Smith ++hdf5->timestep; 6755c6c1daeSBarry Smith PetscFunctionReturn(0); 6765c6c1daeSBarry Smith } 6775c6c1daeSBarry Smith 6783ef9c667SSatish Balay /*@ 6795c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6805c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith Not collective 6835c6c1daeSBarry Smith 6845c6c1daeSBarry Smith Input Parameters: 6855c6c1daeSBarry Smith + viewer - the PetscViewer 6865c6c1daeSBarry Smith - timestep - The timestep number 6875c6c1daeSBarry Smith 6885c6c1daeSBarry Smith Level: intermediate 6895c6c1daeSBarry Smith 6905c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 6915c6c1daeSBarry Smith @*/ 6925c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 6935c6c1daeSBarry Smith { 6945c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6955c6c1daeSBarry Smith 6965c6c1daeSBarry Smith PetscFunctionBegin; 6975c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6985c6c1daeSBarry Smith hdf5->timestep = timestep; 6995c6c1daeSBarry Smith PetscFunctionReturn(0); 7005c6c1daeSBarry Smith } 7015c6c1daeSBarry Smith 7023ef9c667SSatish Balay /*@ 7035c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7045c6c1daeSBarry Smith 7055c6c1daeSBarry Smith Not collective 7065c6c1daeSBarry Smith 7075c6c1daeSBarry Smith Input Parameter: 7085c6c1daeSBarry Smith . viewer - the PetscViewer 7095c6c1daeSBarry Smith 7105c6c1daeSBarry Smith Output Parameter: 7115c6c1daeSBarry Smith . timestep - The timestep number 7125c6c1daeSBarry Smith 7135c6c1daeSBarry Smith Level: intermediate 7145c6c1daeSBarry Smith 7155c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7165c6c1daeSBarry Smith @*/ 7175c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7185c6c1daeSBarry Smith { 7195c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7205c6c1daeSBarry Smith 7215c6c1daeSBarry Smith PetscFunctionBegin; 7225c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7235c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7245c6c1daeSBarry Smith *timestep = hdf5->timestep; 7255c6c1daeSBarry Smith PetscFunctionReturn(0); 7265c6c1daeSBarry Smith } 7275c6c1daeSBarry Smith 72836bce6e8SMatthew G. Knepley /*@C 72936bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 73036bce6e8SMatthew G. Knepley 73136bce6e8SMatthew G. Knepley Not collective 73236bce6e8SMatthew G. Knepley 73336bce6e8SMatthew G. Knepley Input Parameter: 73436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 73536bce6e8SMatthew G. Knepley 73636bce6e8SMatthew G. Knepley Output Parameter: 73736bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 73836bce6e8SMatthew G. Knepley 73936bce6e8SMatthew G. Knepley Level: advanced 74036bce6e8SMatthew G. Knepley 74136bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 74236bce6e8SMatthew G. Knepley @*/ 74336bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 74436bce6e8SMatthew G. Knepley { 74536bce6e8SMatthew G. Knepley PetscFunctionBegin; 74636bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 74736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 74836bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 74936bce6e8SMatthew G. Knepley #else 75036bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 75136bce6e8SMatthew G. Knepley #endif 75236bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 75336bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 75436bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 75536bce6e8SMatthew G. Knepley else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_DOUBLE; 75636bce6e8SMatthew G. Knepley else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_DOUBLE; 75736bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 75836bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 75936bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7607e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 76136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 76236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 76336bce6e8SMatthew G. Knepley } 76436bce6e8SMatthew G. Knepley 76536bce6e8SMatthew G. Knepley /*@C 76636bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 76736bce6e8SMatthew G. Knepley 76836bce6e8SMatthew G. Knepley Not collective 76936bce6e8SMatthew G. Knepley 77036bce6e8SMatthew G. Knepley Input Parameter: 77136bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 77236bce6e8SMatthew G. Knepley 77336bce6e8SMatthew G. Knepley Output Parameter: 77436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 77536bce6e8SMatthew G. Knepley 77636bce6e8SMatthew G. Knepley Level: advanced 77736bce6e8SMatthew G. Knepley 77836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 77936bce6e8SMatthew G. Knepley @*/ 78036bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 78136bce6e8SMatthew G. Knepley { 78236bce6e8SMatthew G. Knepley PetscFunctionBegin; 78336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 78436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 78536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 78636bce6e8SMatthew G. Knepley #else 78736bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 78836bce6e8SMatthew G. Knepley #endif 78936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 79036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 79136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 79236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 79336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 79436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 7957e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 79636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 79736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 79836bce6e8SMatthew G. Knepley } 79936bce6e8SMatthew G. Knepley 800df863907SAlex Fikl /*@C 801b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 80236bce6e8SMatthew G. Knepley 80336bce6e8SMatthew G. Knepley Input Parameters: 80436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 80557d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 80636bce6e8SMatthew G. Knepley . name - The attribute name 80736bce6e8SMatthew G. Knepley . datatype - The attribute type 80836bce6e8SMatthew G. Knepley - value - The attribute value 80936bce6e8SMatthew G. Knepley 81036bce6e8SMatthew G. Knepley Level: advanced 81136bce6e8SMatthew G. Knepley 812577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 81336bce6e8SMatthew G. Knepley @*/ 81457d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 81536bce6e8SMatthew G. Knepley { 81657d22f7dSVaclav Hapla char *parent; 81760568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 818080f144cSVaclav Hapla PetscBool has; 81936bce6e8SMatthew G. Knepley PetscErrorCode ierr; 82036bce6e8SMatthew G. Knepley 82136bce6e8SMatthew G. Knepley PetscFunctionBegin; 8225cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 82357d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 824c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 825b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 82657d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 827bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 828b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 82936bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8307e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8317e97c476SMatthew G. Knepley size_t len; 8327e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 833729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8347e97c476SMatthew G. Knepley } 83536bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 836729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 83760568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 838080f144cSVaclav Hapla if (has) { 839080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 840080f144cSVaclav Hapla } else { 84160568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 842080f144cSVaclav Hapla } 843729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 844729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 845729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 84660568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 847729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 84857d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 84936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 85036bce6e8SMatthew G. Knepley } 85136bce6e8SMatthew G. Knepley 852df863907SAlex Fikl /*@C 853577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 854577e0e04SVaclav Hapla 855577e0e04SVaclav Hapla Input Parameters: 856577e0e04SVaclav Hapla + viewer - The HDF5 viewer 857577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 858577e0e04SVaclav Hapla . name - The attribute name 859577e0e04SVaclav Hapla . datatype - The attribute type 860577e0e04SVaclav Hapla - value - The attribute value 861577e0e04SVaclav Hapla 862577e0e04SVaclav Hapla Notes: 863577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 864577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 865577e0e04SVaclav Hapla 866577e0e04SVaclav Hapla Level: advanced 867577e0e04SVaclav Hapla 868577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 869577e0e04SVaclav Hapla @*/ 870577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 871577e0e04SVaclav Hapla { 872577e0e04SVaclav Hapla PetscErrorCode ierr; 873577e0e04SVaclav Hapla 874577e0e04SVaclav Hapla PetscFunctionBegin; 875577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 876577e0e04SVaclav Hapla PetscValidHeader(obj,2); 877577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 878b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 879577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 880577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 881577e0e04SVaclav Hapla PetscFunctionReturn(0); 882577e0e04SVaclav Hapla } 883577e0e04SVaclav Hapla 884577e0e04SVaclav Hapla /*@C 885b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 886ad0c61feSMatthew G. Knepley 887ad0c61feSMatthew G. Knepley Input Parameters: 888ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 88957d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 890ad0c61feSMatthew G. Knepley . name - The attribute name 891ad0c61feSMatthew G. Knepley - datatype - The attribute type 892ad0c61feSMatthew G. Knepley 893ad0c61feSMatthew G. Knepley Output Parameter: 894ad0c61feSMatthew G. Knepley . value - The attribute value 895ad0c61feSMatthew G. Knepley 896ad0c61feSMatthew G. Knepley Level: advanced 897ad0c61feSMatthew G. Knepley 898577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 899ad0c61feSMatthew G. Knepley @*/ 90057d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 901ad0c61feSMatthew G. Knepley { 90257d22f7dSVaclav Hapla char *parent; 903f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 904969748fdSVaclav Hapla PetscBool has; 905ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 906ad0c61feSMatthew G. Knepley 907ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9085cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 90957d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 910c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 911b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 91257d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 913969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 914969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 915969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 916ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 917ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 91860568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 91960568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 920f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 921f0b58479SMatthew G. Knepley size_t len; 922e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 923f0b58479SMatthew G. Knepley PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 924f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 925f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 926f0b58479SMatthew G. Knepley } 92770efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 928729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 929e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 930e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 93157d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 932ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 933ad0c61feSMatthew G. Knepley } 934ad0c61feSMatthew G. Knepley 935577e0e04SVaclav Hapla /*@C 936577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 937577e0e04SVaclav Hapla 938577e0e04SVaclav Hapla Input Parameters: 939577e0e04SVaclav Hapla + viewer - The HDF5 viewer 940577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 941577e0e04SVaclav Hapla . name - The attribute name 942577e0e04SVaclav Hapla - datatype - The attribute type 943577e0e04SVaclav Hapla 944577e0e04SVaclav Hapla Output Parameter: 945577e0e04SVaclav Hapla . value - The attribute value 946577e0e04SVaclav Hapla 947577e0e04SVaclav Hapla Notes: 948577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 949577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 950577e0e04SVaclav Hapla 951577e0e04SVaclav Hapla Level: advanced 952577e0e04SVaclav Hapla 953577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 954577e0e04SVaclav Hapla @*/ 955577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 956577e0e04SVaclav Hapla { 957577e0e04SVaclav Hapla PetscErrorCode ierr; 958577e0e04SVaclav Hapla 959577e0e04SVaclav Hapla PetscFunctionBegin; 960577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 961577e0e04SVaclav Hapla PetscValidHeader(obj,2); 962577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 963b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 964577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 965577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 966577e0e04SVaclav Hapla PetscFunctionReturn(0); 967577e0e04SVaclav Hapla } 968577e0e04SVaclav Hapla 969a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 970a75c3fd4SVaclav Hapla { 971a75c3fd4SVaclav Hapla htri_t exists; 972a75c3fd4SVaclav Hapla hid_t group; 973a75c3fd4SVaclav Hapla 974a75c3fd4SVaclav Hapla PetscFunctionBegin; 975a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 976a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 977a75c3fd4SVaclav Hapla if (!exists && createGroup) { 978a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 979a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 980a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 981a75c3fd4SVaclav Hapla } 982a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 983a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 984a75c3fd4SVaclav Hapla } 985a75c3fd4SVaclav Hapla 986bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 9875cdeb108SMatthew G. Knepley { 988*90e03537SVaclav Hapla const char rootGroupName[] = "/"; 9895cdeb108SMatthew G. Knepley hid_t h5; 990e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 99185688ae2SVaclav Hapla PetscInt i,n; 99285688ae2SVaclav Hapla char **hierarchy; 99385688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 9945cdeb108SMatthew G. Knepley PetscErrorCode ierr; 9955cdeb108SMatthew G. Knepley 9965cdeb108SMatthew G. Knepley PetscFunctionBegin; 9975cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 998*90e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 999*90e03537SVaclav Hapla else name = rootGroupName; 1000ccd66a83SVaclav Hapla if (has) { 100156cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10025cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1003ccd66a83SVaclav Hapla } 1004ccd66a83SVaclav Hapla if (otype) { 1005ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 100656cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1007ccd66a83SVaclav Hapla } 1008c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 100985688ae2SVaclav Hapla 101085688ae2SVaclav Hapla /* 101185688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 101285688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 101385688ae2SVaclav Hapla 1) whether it's a valid link 101485688ae2SVaclav Hapla 2) whether this link resolves to an object 101585688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 101685688ae2SVaclav Hapla */ 101785688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 101885688ae2SVaclav Hapla if (!n) { 101985688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1020ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1021ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 102285688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 102385688ae2SVaclav Hapla PetscFunctionReturn(0); 102485688ae2SVaclav Hapla } 102585688ae2SVaclav Hapla for (i=0; i<n; i++) { 102685688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 102785688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1028a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1029a75c3fd4SVaclav Hapla if (!exists) break; 103085688ae2SVaclav Hapla } 103185688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 103285688ae2SVaclav Hapla 103385688ae2SVaclav Hapla /* If the object exists, get its type */ 1034ccd66a83SVaclav Hapla if (exists && otype) { 10355cdeb108SMatthew G. Knepley H5O_info_t info; 10365cdeb108SMatthew G. Knepley 103704633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 103856cc0592SVaclav Hapla *otype = info.type; 10395cdeb108SMatthew G. Knepley } 1040ccd66a83SVaclav Hapla if (has) *has = exists; 10415cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 10425cdeb108SMatthew G. Knepley } 10435cdeb108SMatthew G. Knepley 1044ecce1506SVaclav Hapla /*@ 10450a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 10460a9f38efSVaclav Hapla 10470a9f38efSVaclav Hapla Input Parameters: 10480a9f38efSVaclav Hapla . viewer - The HDF5 viewer 10490a9f38efSVaclav Hapla 10500a9f38efSVaclav Hapla Output Parameter: 10510a9f38efSVaclav Hapla . has - Flag for group existence 10520a9f38efSVaclav Hapla 10530a9f38efSVaclav Hapla Notes: 10540a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 10550a9f38efSVaclav Hapla 10560a9f38efSVaclav Hapla Level: advanced 10570a9f38efSVaclav Hapla 10580a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 10590a9f38efSVaclav Hapla @*/ 10600a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 10610a9f38efSVaclav Hapla { 10620a9f38efSVaclav Hapla H5O_type_t type; 10630a9f38efSVaclav Hapla const char *name; 10640a9f38efSVaclav Hapla PetscErrorCode ierr; 10650a9f38efSVaclav Hapla 10660a9f38efSVaclav Hapla PetscFunctionBegin; 10670a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10680a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 10690a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 10700a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 10710a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 10720a9f38efSVaclav Hapla PetscFunctionReturn(0); 10730a9f38efSVaclav Hapla } 10740a9f38efSVaclav Hapla 10750a9f38efSVaclav Hapla /*@ 1076e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1077ecce1506SVaclav Hapla 1078ecce1506SVaclav Hapla Input Parameters: 1079ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1080ecce1506SVaclav Hapla - obj - The named object 1081ecce1506SVaclav Hapla 1082ecce1506SVaclav Hapla Output Parameter: 1083ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1084ecce1506SVaclav Hapla 1085e3f143f7SVaclav Hapla Notes: 1086e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1087e3f143f7SVaclav Hapla 1088ecce1506SVaclav Hapla Level: advanced 1089ecce1506SVaclav Hapla 1090e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1091ecce1506SVaclav Hapla @*/ 1092ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1093ecce1506SVaclav Hapla { 109456cc0592SVaclav Hapla H5O_type_t type; 10956c132bc1SVaclav Hapla char *path; 1096ecce1506SVaclav Hapla PetscErrorCode ierr; 1097ecce1506SVaclav Hapla 1098ecce1506SVaclav Hapla PetscFunctionBegin; 1099c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1100c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1101c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1102ecce1506SVaclav Hapla *has = PETSC_FALSE; 1103e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11046c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11056c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 110656cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11076c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1108ecce1506SVaclav Hapla PetscFunctionReturn(0); 1109ecce1506SVaclav Hapla } 1110ecce1506SVaclav Hapla 1111df863907SAlex Fikl /*@C 1112ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1113e7bdbf8eSMatthew G. Knepley 1114e7bdbf8eSMatthew G. Knepley Input Parameters: 1115e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 111610e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1117e7bdbf8eSMatthew G. Knepley - name - The attribute name 1118e7bdbf8eSMatthew G. Knepley 1119e7bdbf8eSMatthew G. Knepley Output Parameter: 1120e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1121e7bdbf8eSMatthew G. Knepley 1122e7bdbf8eSMatthew G. Knepley Level: advanced 1123e7bdbf8eSMatthew G. Knepley 1124577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1125e7bdbf8eSMatthew G. Knepley @*/ 112610e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1127e7bdbf8eSMatthew G. Knepley { 112810e69e8fSVaclav Hapla char *parent; 1129e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1130e7bdbf8eSMatthew G. Knepley 1131e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 11325cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 113310e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1134c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1135c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 113610e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1137bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 113810e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 113910e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 114006db490cSVaclav Hapla PetscFunctionReturn(0); 114106db490cSVaclav Hapla } 114206db490cSVaclav Hapla 1143577e0e04SVaclav Hapla /*@C 1144577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1145577e0e04SVaclav Hapla 1146577e0e04SVaclav Hapla Input Parameters: 1147577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1148577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1149577e0e04SVaclav Hapla - name - The attribute name 1150577e0e04SVaclav Hapla 1151577e0e04SVaclav Hapla Output Parameter: 1152577e0e04SVaclav Hapla . has - Flag for attribute existence 1153577e0e04SVaclav Hapla 1154577e0e04SVaclav Hapla Notes: 1155577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1156577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1157577e0e04SVaclav Hapla 1158577e0e04SVaclav Hapla Level: advanced 1159577e0e04SVaclav Hapla 1160577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1161577e0e04SVaclav Hapla @*/ 1162577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1163577e0e04SVaclav Hapla { 1164577e0e04SVaclav Hapla PetscErrorCode ierr; 1165577e0e04SVaclav Hapla 1166577e0e04SVaclav Hapla PetscFunctionBegin; 1167577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1168577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1169577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1170577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1171577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1172577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1173577e0e04SVaclav Hapla PetscFunctionReturn(0); 1174577e0e04SVaclav Hapla } 1175577e0e04SVaclav Hapla 117606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 117706db490cSVaclav Hapla { 117806db490cSVaclav Hapla hid_t h5; 117906db490cSVaclav Hapla htri_t hhas; 118006db490cSVaclav Hapla PetscErrorCode ierr; 118106db490cSVaclav Hapla 118206db490cSVaclav Hapla PetscFunctionBegin; 118306db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11842f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 11855cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1186e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1187e7bdbf8eSMatthew G. Knepley } 1188e7bdbf8eSMatthew G. Knepley 1189eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadInitialize_Private(PetscViewer viewer, const char name[], HDF5ReadCtx *ctx) 1190a5e1feadSVaclav Hapla { 1191fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1192a5e1feadSVaclav Hapla PetscErrorCode ierr; 1193a5e1feadSVaclav Hapla 1194a5e1feadSVaclav Hapla PetscFunctionBegin; 119569a06e7bSVaclav Hapla ierr = PetscNew(&h);CHKERRQ(ierr); 119669a06e7bSVaclav Hapla ierr = PetscViewerHDF5OpenGroup(viewer, &h->file, &h->group);CHKERRQ(ierr); 119769a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataset,H5Dopen2,(h->group, name, H5P_DEFAULT)); 119869a06e7bSVaclav Hapla PetscStackCallHDF5Return(h->dataspace,H5Dget_space,(h->dataset)); 11999568d583SVaclav Hapla ierr = PetscViewerHDF5GetTimestep(viewer, &h->timestep);CHKERRQ(ierr); 1200adb363a2SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,name,"complex",&h->complexVal);CHKERRQ(ierr); 120109dabeb0SVaclav Hapla /* MATLAB stores column vectors horizontally */ 1202adb363a2SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer,name,"MATLAB_class",&h->horizontal);CHKERRQ(ierr); 1203b8ef406cSVaclav Hapla /* Create property list for collective dataset read */ 1204b8ef406cSVaclav Hapla PetscStackCallHDF5Return(h->plist,H5Pcreate,(H5P_DATASET_XFER)); 1205b8ef406cSVaclav Hapla #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 1206b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(h->plist, H5FD_MPIO_COLLECTIVE)); 1207b8ef406cSVaclav Hapla #endif 1208b8ef406cSVaclav Hapla /* To write dataset independently use H5Pset_dxpl_mpio(plist_id, H5FD_MPIO_INDEPENDENT) */ 120969a06e7bSVaclav Hapla *ctx = h; 121069a06e7bSVaclav Hapla PetscFunctionReturn(0); 121169a06e7bSVaclav Hapla } 121269a06e7bSVaclav Hapla 1213eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadFinalize_Private(PetscViewer viewer, HDF5ReadCtx *ctx) 121469a06e7bSVaclav Hapla { 121569a06e7bSVaclav Hapla HDF5ReadCtx h; 121669a06e7bSVaclav Hapla PetscErrorCode ierr; 121769a06e7bSVaclav Hapla 121869a06e7bSVaclav Hapla PetscFunctionBegin; 121969a06e7bSVaclav Hapla h = *ctx; 1220b8ef406cSVaclav Hapla PetscStackCallHDF5(H5Pclose,(h->plist)); 1221*90e03537SVaclav Hapla PetscStackCallHDF5(H5Gclose,(h->group)); 122269a06e7bSVaclav Hapla PetscStackCallHDF5(H5Sclose,(h->dataspace)); 122369a06e7bSVaclav Hapla PetscStackCallHDF5(H5Dclose,(h->dataset)); 122469a06e7bSVaclav Hapla ierr = PetscFree(*ctx);CHKERRQ(ierr); 122569a06e7bSVaclav Hapla PetscFunctionReturn(0); 122669a06e7bSVaclav Hapla } 122769a06e7bSVaclav Hapla 1228eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSizes_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout *map_) 122969a06e7bSVaclav Hapla { 123069a06e7bSVaclav Hapla int rdim, dim; 123169a06e7bSVaclav Hapla hsize_t dims[4]; 123209dabeb0SVaclav Hapla PetscInt bsInd, lenInd, bs, len, N; 12338374c777SVaclav Hapla PetscLayout map; 12348374c777SVaclav Hapla PetscErrorCode ierr; 123569a06e7bSVaclav Hapla 123669a06e7bSVaclav Hapla PetscFunctionBegin; 12378374c777SVaclav Hapla if (!(*map_)) { 12388374c777SVaclav Hapla ierr = PetscLayoutCreate(PetscObjectComm((PetscObject)viewer),map_);CHKERRQ(ierr); 12398374c777SVaclav Hapla } 12408374c777SVaclav Hapla map = *map_; 12419787f08bSVaclav Hapla /* calculate expected number of dimensions */ 1242a5e1feadSVaclav Hapla dim = 0; 12439568d583SVaclav Hapla if (ctx->timestep >= 0) ++dim; 1244a5e1feadSVaclav Hapla ++dim; /* length in blocks */ 12459568d583SVaclav Hapla if (ctx->complexVal) ++dim; 12469787f08bSVaclav Hapla /* get actual number of dimensions in dataset */ 124769a06e7bSVaclav Hapla PetscStackCallHDF5Return(rdim,H5Sget_simple_extent_dims,(ctx->dataspace, dims, NULL)); 12489787f08bSVaclav Hapla /* calculate expected dimension indices */ 12499787f08bSVaclav Hapla lenInd = 0; 12509568d583SVaclav Hapla if (ctx->timestep >= 0) ++lenInd; 12519787f08bSVaclav Hapla bsInd = lenInd + 1; 12529568d583SVaclav Hapla ctx->dim2 = PETSC_FALSE; 12539787f08bSVaclav Hapla if (rdim == dim) { 125445e21b7fSVaclav Hapla bs = 1; /* support vectors stored as 1D array */ 12559787f08bSVaclav Hapla } else if (rdim == dim+1) { 125645e21b7fSVaclav Hapla bs = (PetscInt) dims[bsInd]; 12579568d583SVaclav Hapla if (bs == 1) ctx->dim2 = PETSC_TRUE; /* vector with blocksize of 1, still stored as 2D array */ 1258a5e1feadSVaclav Hapla } else { 1259b54f1188SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Number of dimensions %d not %d as expected", rdim, dim); 1260a5e1feadSVaclav Hapla } 126109dabeb0SVaclav Hapla len = dims[lenInd]; 126209dabeb0SVaclav Hapla if (ctx->horizontal) { 1263b54f1188SVaclav 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."); 126409dabeb0SVaclav Hapla len = bs; 126509dabeb0SVaclav Hapla bs = 1; 126609dabeb0SVaclav Hapla } 126709dabeb0SVaclav Hapla N = (PetscInt) len*bs; 12688374c777SVaclav Hapla 12698374c777SVaclav Hapla /* Set Vec sizes,blocksize,and type if not already set */ 12708374c777SVaclav Hapla if (map->bs < 0) { 127145e21b7fSVaclav Hapla ierr = PetscLayoutSetBlockSize(map, bs);CHKERRQ(ierr); 127245e21b7fSVaclav 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); 12738374c777SVaclav Hapla if (map->N < 0) { 127445e21b7fSVaclav Hapla ierr = PetscLayoutSetSize(map, N);CHKERRQ(ierr); 1275b54f1188SVaclav 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); 127669a06e7bSVaclav Hapla PetscFunctionReturn(0); 127769a06e7bSVaclav Hapla } 127869a06e7bSVaclav Hapla 1279eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadSelectHyperslab_Private(PetscViewer viewer, HDF5ReadCtx ctx, PetscLayout map, hid_t *memspace) 128016f6402dSVaclav Hapla { 128116f6402dSVaclav Hapla hsize_t count[4], offset[4]; 128216f6402dSVaclav Hapla int dim; 128316f6402dSVaclav Hapla PetscInt bs, n, low; 128416f6402dSVaclav Hapla PetscErrorCode ierr; 128516f6402dSVaclav Hapla 128616f6402dSVaclav Hapla PetscFunctionBegin; 128716f6402dSVaclav Hapla /* Compute local size and ownership range */ 128816f6402dSVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 128916f6402dSVaclav Hapla ierr = PetscLayoutGetBlockSize(map, &bs);CHKERRQ(ierr); 129016f6402dSVaclav Hapla ierr = PetscLayoutGetLocalSize(map, &n);CHKERRQ(ierr); 129116f6402dSVaclav Hapla ierr = PetscLayoutGetRange(map, &low, NULL);CHKERRQ(ierr); 129216f6402dSVaclav Hapla 129316f6402dSVaclav Hapla /* Each process defines a dataset and reads it from the hyperslab in the file */ 129416f6402dSVaclav Hapla dim = 0; 129516f6402dSVaclav Hapla if (ctx->timestep >= 0) { 129616f6402dSVaclav Hapla count[dim] = 1; 129716f6402dSVaclav Hapla offset[dim] = ctx->timestep; 129816f6402dSVaclav Hapla ++dim; 129916f6402dSVaclav Hapla } 130009dabeb0SVaclav Hapla if (ctx->horizontal) { 130109dabeb0SVaclav Hapla count[dim] = 1; 130209dabeb0SVaclav Hapla offset[dim] = 0; 130309dabeb0SVaclav Hapla ++dim; 130409dabeb0SVaclav Hapla } 130516f6402dSVaclav Hapla { 130616f6402dSVaclav Hapla ierr = PetscHDF5IntCast(n/bs, &count[dim]);CHKERRQ(ierr); 130716f6402dSVaclav Hapla ierr = PetscHDF5IntCast(low/bs, &offset[dim]);CHKERRQ(ierr); 130816f6402dSVaclav Hapla ++dim; 130916f6402dSVaclav Hapla } 131016f6402dSVaclav Hapla if (bs > 1 || ctx->dim2) { 131109dabeb0SVaclav Hapla if (PetscUnlikely(ctx->horizontal)) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "cannot have horizontal array with blocksize > 1"); 131216f6402dSVaclav Hapla count[dim] = bs; 131316f6402dSVaclav Hapla offset[dim] = 0; 131416f6402dSVaclav Hapla ++dim; 131516f6402dSVaclav Hapla } 131616f6402dSVaclav Hapla if (ctx->complexVal) { 131716f6402dSVaclav Hapla count[dim] = 2; 131816f6402dSVaclav Hapla offset[dim] = 0; 131916f6402dSVaclav Hapla ++dim; 132016f6402dSVaclav Hapla } 132116f6402dSVaclav Hapla PetscStackCallHDF5Return(*memspace,H5Screate_simple,(dim, count, NULL)); 132216f6402dSVaclav Hapla PetscStackCallHDF5(H5Sselect_hyperslab,(ctx->dataspace, H5S_SELECT_SET, offset, NULL, count, NULL)); 132316f6402dSVaclav Hapla PetscFunctionReturn(0); 132416f6402dSVaclav Hapla } 132516f6402dSVaclav Hapla 1326eb5a92b4SVaclav Hapla static PetscErrorCode PetscViewerHDF5ReadArray_Private(PetscViewer viewer, HDF5ReadCtx h, hid_t datatype, hid_t memspace, void *arr) 1327ef2d82ceSVaclav Hapla { 1328ef2d82ceSVaclav Hapla PetscFunctionBegin; 1329ef2d82ceSVaclav Hapla PetscStackCallHDF5(H5Dread,(h->dataset, datatype, memspace, h->dataspace, h->plist, arr)); 1330ef2d82ceSVaclav Hapla PetscFunctionReturn(0); 1331ef2d82ceSVaclav Hapla } 1332ef2d82ceSVaclav Hapla 1333dad982a8SVaclav Hapla PetscErrorCode PetscViewerHDF5Load(PetscViewer viewer, const char *name, PetscLayout map, hid_t datatype, void **newarr) 1334a25c73c6SVaclav Hapla { 1335fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 1336fbd37863SVaclav Hapla hid_t memspace=0; 1337a25c73c6SVaclav Hapla size_t unitsize; 1338a25c73c6SVaclav Hapla void *arr; 1339a25c73c6SVaclav Hapla PetscErrorCode ierr; 1340a25c73c6SVaclav Hapla 1341a25c73c6SVaclav Hapla PetscFunctionBegin; 1342eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 13435a89fdf4SVaclav Hapla #if defined(PETSC_USE_COMPLEX) 13445a89fdf4SVaclav Hapla if (!h->complexVal) { 1345c76551afSVaclav Hapla H5T_class_t clazz = H5Tget_class(datatype); 1346c76551afSVaclav 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."); 13475a89fdf4SVaclav Hapla } 13485a89fdf4SVaclav Hapla #else 13495a89fdf4SVaclav 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."); 13505a89fdf4SVaclav Hapla #endif 1351e9e90110SVaclav Hapla 1352e9e90110SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1353e9e90110SVaclav Hapla ierr = PetscLayoutSetUp(map);CHKERRQ(ierr); 1354e9e90110SVaclav Hapla ierr = PetscViewerHDF5ReadSelectHyperslab_Private(viewer, h, map, &memspace);CHKERRQ(ierr); 1355e9e90110SVaclav Hapla 13564fc17bcdSVaclav Hapla unitsize = H5Tget_size(datatype); 13574fc17bcdSVaclav Hapla if (h->complexVal) unitsize *= 2; 1358dff35581SVaclav 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); 13594fc17bcdSVaclav Hapla ierr = PetscMalloc(map->n*unitsize, &arr);CHKERRQ(ierr); 13604fc17bcdSVaclav Hapla 1361eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadArray_Private(viewer, h, datatype, memspace, arr);CHKERRQ(ierr); 1362a25c73c6SVaclav Hapla PetscStackCallHDF5(H5Sclose,(memspace)); 1363eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 1364a25c73c6SVaclav Hapla *newarr = arr; 1365a25c73c6SVaclav Hapla PetscFunctionReturn(0); 1366a25c73c6SVaclav Hapla } 1367a25c73c6SVaclav Hapla 1368c1aaad9cSVaclav Hapla /*@C 1369c1aaad9cSVaclav Hapla PetscViewerHDF5ReadSizes - Read block size and global size of a vector (Vec or IS) stored in an HDF5 file. 1370c1aaad9cSVaclav Hapla 1371c1aaad9cSVaclav Hapla Input Parameters: 1372c1aaad9cSVaclav Hapla + viewer - The HDF5 viewer 1373c1aaad9cSVaclav Hapla - name - The vector name 1374c1aaad9cSVaclav Hapla 1375c1aaad9cSVaclav Hapla Output Parameter: 1376c1aaad9cSVaclav Hapla + bs - block size 1377c1aaad9cSVaclav Hapla - N - global size 1378c1aaad9cSVaclav Hapla 1379c1aaad9cSVaclav Hapla Note: 1380c1aaad9cSVaclav Hapla A vector is stored as an HDF5 dataspace with 1-4 dimensions in this order: 1381c1aaad9cSVaclav Hapla 1) # timesteps (optional), 2) # blocks, 3) # elements per block (optional), 4) real and imaginary part (only for complex). 1382c1aaad9cSVaclav Hapla 1383c1aaad9cSVaclav Hapla A vectors can be stored as a 2D dataspace even if its blocksize is 1; see PetscViewerHDF5SetBaseDimension2(). 1384c1aaad9cSVaclav Hapla 1385c1aaad9cSVaclav Hapla Level: advanced 1386c1aaad9cSVaclav Hapla 1387c1aaad9cSVaclav Hapla .seealso: PetscViewerHDF5Open(), VecLoad(), ISLoad(), VecGetSize(), ISGetSize(), PetscViewerHDF5SetBaseDimension2() 1388c1aaad9cSVaclav Hapla @*/ 138969a06e7bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadSizes(PetscViewer viewer, const char name[], PetscInt *bs, PetscInt *N) 139069a06e7bSVaclav Hapla { 1391fbd37863SVaclav Hapla HDF5ReadCtx h=NULL; 13928374c777SVaclav Hapla PetscLayout map=NULL; 139369a06e7bSVaclav Hapla PetscErrorCode ierr; 139469a06e7bSVaclav Hapla 139569a06e7bSVaclav Hapla PetscFunctionBegin; 1396c1aaad9cSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1397eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadInitialize_Private(viewer, name, &h);CHKERRQ(ierr); 1398eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadSizes_Private(viewer, h, &map);CHKERRQ(ierr); 1399eb5a92b4SVaclav Hapla ierr = PetscViewerHDF5ReadFinalize_Private(viewer, &h);CHKERRQ(ierr); 14008374c777SVaclav Hapla if (bs) *bs = map->bs; 14018374c777SVaclav Hapla if (N) *N = map->N; 14028374c777SVaclav Hapla ierr = PetscLayoutDestroy(&map);CHKERRQ(ierr); 1403a5e1feadSVaclav Hapla PetscFunctionReturn(0); 1404a5e1feadSVaclav Hapla } 1405a5e1feadSVaclav Hapla 1406a75e6a4aSMatthew G. Knepley /* 1407a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1408a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1409a75e6a4aSMatthew G. Knepley */ 1410d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1411a75e6a4aSMatthew G. Knepley 1412a75e6a4aSMatthew G. Knepley /*@C 1413a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1414a75e6a4aSMatthew G. Knepley 1415a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1416a75e6a4aSMatthew G. Knepley 1417a75e6a4aSMatthew G. Knepley Input Parameter: 1418a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1419a75e6a4aSMatthew G. Knepley 1420a75e6a4aSMatthew G. Knepley Level: intermediate 1421a75e6a4aSMatthew G. Knepley 1422a75e6a4aSMatthew G. Knepley Options Database Keys: 1423a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1424a75e6a4aSMatthew G. Knepley 1425a75e6a4aSMatthew G. Knepley Environmental variables: 1426a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1427a75e6a4aSMatthew G. Knepley 1428a75e6a4aSMatthew G. Knepley Notes: 1429a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1430a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1431a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1432a75e6a4aSMatthew G. Knepley 1433a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1434a75e6a4aSMatthew G. Knepley @*/ 1435a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1436a75e6a4aSMatthew G. Knepley { 1437a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1438a75e6a4aSMatthew G. Knepley PetscBool flg; 1439a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1440a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1441a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1442a75e6a4aSMatthew G. Knepley 1443a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1444a75e6a4aSMatthew 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);} 1445a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 144612801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1447a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1448a75e6a4aSMatthew G. Knepley } 144947435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1450a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1451a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1452a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1453a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1454a75e6a4aSMatthew G. Knepley if (!flg) { 1455a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1456a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1457a75e6a4aSMatthew G. Knepley } 1458a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1459a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1460a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1461a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 146247435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1463a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1464a75e6a4aSMatthew G. Knepley } 1465a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1466a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1467a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1468a75e6a4aSMatthew G. Knepley } 1469