1*bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h> 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; 27b723ab35SVaclav Hapla PetscBool mataij_names_set; 285c6c1daeSBarry Smith } PetscViewer_HDF5; 295c6c1daeSBarry Smith 306c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 316c132bc1SVaclav Hapla { 326c132bc1SVaclav Hapla const char *group; 336c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 346c132bc1SVaclav Hapla PetscErrorCode ierr; 356c132bc1SVaclav Hapla 366c132bc1SVaclav Hapla PetscFunctionBegin; 376c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 386c132bc1SVaclav Hapla ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 396c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 406c132bc1SVaclav Hapla ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 416c132bc1SVaclav Hapla ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 426c132bc1SVaclav Hapla PetscFunctionReturn(0); 436c132bc1SVaclav Hapla } 446c132bc1SVaclav Hapla 45577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 46577e0e04SVaclav Hapla { 47577e0e04SVaclav Hapla PetscBool has; 48577e0e04SVaclav Hapla const char *group; 49577e0e04SVaclav Hapla PetscErrorCode ierr; 50577e0e04SVaclav Hapla 51577e0e04SVaclav Hapla PetscFunctionBegin; 52577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 53577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 54577e0e04SVaclav Hapla if (!has) { 55577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 56577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 57577e0e04SVaclav Hapla } 58577e0e04SVaclav Hapla PetscFunctionReturn(0); 59577e0e04SVaclav Hapla } 60577e0e04SVaclav Hapla 614416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 6282ea9b62SBarry Smith { 6382ea9b62SBarry Smith PetscErrorCode ierr; 6482ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6582ea9b62SBarry Smith 6682ea9b62SBarry Smith PetscFunctionBegin; 6782ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 6882ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 699a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 7082ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7182ea9b62SBarry Smith PetscFunctionReturn(0); 7282ea9b62SBarry Smith } 7382ea9b62SBarry Smith 745c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 755c6c1daeSBarry Smith { 765c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 775c6c1daeSBarry Smith PetscErrorCode ierr; 785c6c1daeSBarry Smith 795c6c1daeSBarry Smith PetscFunctionBegin; 805c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 81729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 825c6c1daeSBarry Smith PetscFunctionReturn(0); 835c6c1daeSBarry Smith } 845c6c1daeSBarry Smith 859b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 865c6c1daeSBarry Smith { 875c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 885c6c1daeSBarry Smith PetscErrorCode ierr; 895c6c1daeSBarry Smith 905c6c1daeSBarry Smith PetscFunctionBegin; 915c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 925c6c1daeSBarry Smith while (hdf5->groups) { 935c6c1daeSBarry Smith GroupList *tmp = hdf5->groups->next; 945c6c1daeSBarry Smith 955c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 965c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 975c6c1daeSBarry Smith hdf5->groups = tmp; 985c6c1daeSBarry Smith } 9961912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 10061912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 10161912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 10261912a0cSVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 1035c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 1040b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 105d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 1060b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 107058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 108058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 109058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetAIJNames_C",NULL);CHKERRQ(ierr); 110058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetAIJNames_C",NULL);CHKERRQ(ierr); 1115c6c1daeSBarry Smith PetscFunctionReturn(0); 1125c6c1daeSBarry Smith } 1135c6c1daeSBarry Smith 1149b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1155c6c1daeSBarry Smith { 1165c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1175c6c1daeSBarry Smith 1185c6c1daeSBarry Smith PetscFunctionBegin; 1195c6c1daeSBarry Smith hdf5->btype = type; 1205c6c1daeSBarry Smith PetscFunctionReturn(0); 1215c6c1daeSBarry Smith } 1225c6c1daeSBarry Smith 1230b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1240b62783dSVaclav Hapla { 1250b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1260b62783dSVaclav Hapla 1270b62783dSVaclav Hapla PetscFunctionBegin; 1280b62783dSVaclav Hapla *type = hdf5->btype; 1290b62783dSVaclav Hapla PetscFunctionReturn(0); 1300b62783dSVaclav Hapla } 1310b62783dSVaclav Hapla 1329b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13382ea9b62SBarry Smith { 13482ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13582ea9b62SBarry Smith 13682ea9b62SBarry Smith PetscFunctionBegin; 13782ea9b62SBarry Smith hdf5->basedimension2 = flg; 13882ea9b62SBarry Smith PetscFunctionReturn(0); 13982ea9b62SBarry Smith } 14082ea9b62SBarry Smith 141df863907SAlex Fikl /*@ 14282ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14382ea9b62SBarry Smith dimension of 2. 14482ea9b62SBarry Smith 14582ea9b62SBarry Smith Logically Collective on PetscViewer 14682ea9b62SBarry Smith 14782ea9b62SBarry Smith Input Parameters: 14882ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14982ea9b62SBarry 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 15082ea9b62SBarry Smith 15182ea9b62SBarry Smith Options Database: 15282ea9b62SBarry 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 15382ea9b62SBarry Smith 15482ea9b62SBarry Smith 15595452b02SPatrick Sanan Notes: 15695452b02SPatrick 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 15782ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith Level: intermediate 16082ea9b62SBarry Smith 16182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16282ea9b62SBarry Smith 16382ea9b62SBarry Smith @*/ 16482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16582ea9b62SBarry Smith { 16682ea9b62SBarry Smith PetscErrorCode ierr; 16782ea9b62SBarry Smith 16882ea9b62SBarry Smith PetscFunctionBegin; 16982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 17082ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 17182ea9b62SBarry Smith PetscFunctionReturn(0); 17282ea9b62SBarry Smith } 17382ea9b62SBarry Smith 174df863907SAlex Fikl /*@ 17582ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17682ea9b62SBarry Smith dimension of 2. 17782ea9b62SBarry Smith 17882ea9b62SBarry Smith Logically Collective on PetscViewer 17982ea9b62SBarry Smith 18082ea9b62SBarry Smith Input Parameter: 18182ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18282ea9b62SBarry Smith 18382ea9b62SBarry Smith Output Parameter: 18482ea9b62SBarry 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 18582ea9b62SBarry Smith 18695452b02SPatrick Sanan Notes: 18795452b02SPatrick 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 18882ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18982ea9b62SBarry Smith 19082ea9b62SBarry Smith Level: intermediate 19182ea9b62SBarry Smith 19282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 19382ea9b62SBarry Smith 19482ea9b62SBarry Smith @*/ 19582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19682ea9b62SBarry Smith { 19782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19882ea9b62SBarry Smith 19982ea9b62SBarry Smith PetscFunctionBegin; 20082ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 20182ea9b62SBarry Smith *flg = hdf5->basedimension2; 20282ea9b62SBarry Smith PetscFunctionReturn(0); 20382ea9b62SBarry Smith } 20482ea9b62SBarry Smith 2059b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2069a0502c6SHåkon Strandenes { 2079a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2089a0502c6SHåkon Strandenes 2099a0502c6SHåkon Strandenes PetscFunctionBegin; 2109a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2119a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2129a0502c6SHåkon Strandenes } 2139a0502c6SHåkon Strandenes 214df863907SAlex Fikl /*@ 2159a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2169a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2179a0502c6SHåkon Strandenes 2189a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2199a0502c6SHåkon Strandenes 2209a0502c6SHåkon Strandenes Input Parameters: 2219a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2229a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes Options Database: 2259a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2269a0502c6SHåkon Strandenes 2279a0502c6SHåkon Strandenes 22895452b02SPatrick Sanan Notes: 22995452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2309a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2319a0502c6SHåkon Strandenes 2329a0502c6SHåkon Strandenes Level: intermediate 2339a0502c6SHåkon Strandenes 2349a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2359a0502c6SHåkon Strandenes PetscReal 2369a0502c6SHåkon Strandenes 2379a0502c6SHåkon Strandenes @*/ 2389a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2399a0502c6SHåkon Strandenes { 2409a0502c6SHåkon Strandenes PetscErrorCode ierr; 2419a0502c6SHåkon Strandenes 2429a0502c6SHåkon Strandenes PetscFunctionBegin; 2439a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2449a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2459a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2469a0502c6SHåkon Strandenes } 2479a0502c6SHåkon Strandenes 248df863907SAlex Fikl /*@ 2499a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2509a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2519a0502c6SHåkon Strandenes 2529a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2539a0502c6SHåkon Strandenes 2549a0502c6SHåkon Strandenes Input Parameter: 2559a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2569a0502c6SHåkon Strandenes 2579a0502c6SHåkon Strandenes Output Parameter: 2589a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2599a0502c6SHåkon Strandenes 26095452b02SPatrick Sanan Notes: 26195452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2629a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2639a0502c6SHåkon Strandenes 2649a0502c6SHåkon Strandenes Level: intermediate 2659a0502c6SHåkon Strandenes 2669a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2679a0502c6SHåkon Strandenes PetscReal 2689a0502c6SHåkon Strandenes 2699a0502c6SHåkon Strandenes @*/ 2709a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2719a0502c6SHåkon Strandenes { 2729a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2739a0502c6SHåkon Strandenes 2749a0502c6SHåkon Strandenes PetscFunctionBegin; 2759a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2769a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2779a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2789a0502c6SHåkon Strandenes } 2799a0502c6SHåkon Strandenes 2809b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 2815c6c1daeSBarry Smith { 2825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2835c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 2845c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 2855c6c1daeSBarry Smith #endif 2865c6c1daeSBarry Smith hid_t plist_id; 2875c6c1daeSBarry Smith PetscErrorCode ierr; 2885c6c1daeSBarry Smith 2895c6c1daeSBarry Smith PetscFunctionBegin; 290f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 291f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 2925c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 2935c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 294729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 2955c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 296729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 2975c6c1daeSBarry Smith #endif 2985c6c1daeSBarry Smith /* Create or open the file collectively */ 2995c6c1daeSBarry Smith switch (hdf5->btype) { 3005c6c1daeSBarry Smith case FILE_MODE_READ: 301729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3025c6c1daeSBarry Smith break; 3035c6c1daeSBarry Smith case FILE_MODE_APPEND: 304729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 3055c6c1daeSBarry Smith break; 3065c6c1daeSBarry Smith case FILE_MODE_WRITE: 307729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 3085c6c1daeSBarry Smith break; 3095c6c1daeSBarry Smith default: 3105c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 3115c6c1daeSBarry Smith } 3125c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 313729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 3145c6c1daeSBarry Smith PetscFunctionReturn(0); 3155c6c1daeSBarry Smith } 3165c6c1daeSBarry Smith 317d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 318d1232d7fSToby Isaac { 319d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 320d1232d7fSToby Isaac 321d1232d7fSToby Isaac PetscFunctionBegin; 322d1232d7fSToby Isaac *name = vhdf5->filename; 323d1232d7fSToby Isaac PetscFunctionReturn(0); 324d1232d7fSToby Isaac } 325d1232d7fSToby Isaac 3269b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetAIJNames_HDF5(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 327058bd781SVaclav Hapla { 328058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 329058bd781SVaclav Hapla PetscErrorCode ierr; 330058bd781SVaclav Hapla 331058bd781SVaclav Hapla PetscFunctionBegin; 332058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_iname);CHKERRQ(ierr); 333058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_jname);CHKERRQ(ierr); 334058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_aname);CHKERRQ(ierr); 335058bd781SVaclav Hapla ierr = PetscFree(hdf5->mataij_cname);CHKERRQ(ierr); 336058bd781SVaclav Hapla ierr = PetscStrallocpy(iname,&hdf5->mataij_iname);CHKERRQ(ierr); 337058bd781SVaclav Hapla ierr = PetscStrallocpy(jname,&hdf5->mataij_jname);CHKERRQ(ierr); 338058bd781SVaclav Hapla ierr = PetscStrallocpy(aname,&hdf5->mataij_aname);CHKERRQ(ierr); 339058bd781SVaclav Hapla ierr = PetscStrallocpy(cname,&hdf5->mataij_cname);CHKERRQ(ierr); 340b723ab35SVaclav Hapla hdf5->mataij_names_set = PETSC_TRUE; 341058bd781SVaclav Hapla PetscFunctionReturn(0); 342058bd781SVaclav Hapla } 343058bd781SVaclav Hapla 344058bd781SVaclav Hapla /*@C 345058bd781SVaclav 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. 346058bd781SVaclav Hapla 347058bd781SVaclav Hapla Collective on PetscViewer 348058bd781SVaclav Hapla 349058bd781SVaclav Hapla Input Parameters: 350058bd781SVaclav Hapla + viewer - the PetscViewer; either ASCII or binary 351058bd781SVaclav 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 352058bd781SVaclav Hapla . jname - name of dataset j representing column indices 353058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 354058bd781SVaclav Hapla - cname - name of attribute stoting column count 355058bd781SVaclav Hapla 356058bd781SVaclav Hapla Level: advanced 357058bd781SVaclav Hapla 358058bd781SVaclav Hapla Notes: 359bd7fcbe1SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 360058bd781SVaclav Hapla 361058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5GetAIJNames() 362058bd781SVaclav Hapla @*/ 363bd7fcbe1SVaclav Hapla /* TODO Once corresponding MatView is implemented, add this: 364bd7fcbe1SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols"). 365bd7fcbe1SVaclav Hapla For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded. 366bd7fcbe1SVaclav Hapla */ 367058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5SetAIJNames(PetscViewer viewer, const char iname[], const char jname[], const char aname[], const char cname[]) 368058bd781SVaclav Hapla { 369058bd781SVaclav Hapla PetscErrorCode ierr; 370058bd781SVaclav Hapla 371058bd781SVaclav Hapla PetscFunctionBegin; 372058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 373058bd781SVaclav Hapla PetscValidCharPointer(iname,2); 374058bd781SVaclav Hapla PetscValidCharPointer(jname,3); 375058bd781SVaclav Hapla PetscValidCharPointer(aname,4); 376058bd781SVaclav Hapla PetscValidCharPointer(cname,5); 377058bd781SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetAIJNames_C",(PetscViewer,const char[],const char[],const char[],const char[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 378058bd781SVaclav Hapla PetscFunctionReturn(0); 379058bd781SVaclav Hapla } 380058bd781SVaclav Hapla 3819b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAIJNames_HDF5(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 382058bd781SVaclav Hapla { 383058bd781SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 384058bd781SVaclav Hapla 385058bd781SVaclav Hapla PetscFunctionBegin; 386058bd781SVaclav Hapla *iname = hdf5->mataij_iname; 387058bd781SVaclav Hapla *jname = hdf5->mataij_jname; 388058bd781SVaclav Hapla *aname = hdf5->mataij_aname; 389058bd781SVaclav Hapla *cname = hdf5->mataij_cname; 390058bd781SVaclav Hapla PetscFunctionReturn(0); 391058bd781SVaclav Hapla } 392058bd781SVaclav Hapla 393058bd781SVaclav Hapla /*@C 394058bd781SVaclav 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. 395058bd781SVaclav Hapla 396058bd781SVaclav Hapla Collective on PetscViewer 397058bd781SVaclav Hapla 398058bd781SVaclav Hapla Input Parameters: 399058bd781SVaclav Hapla . viewer - the PetscViewer; either ASCII or binary 400058bd781SVaclav Hapla 401058bd781SVaclav Hapla Output Parameters: 402058bd781SVaclav 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 403058bd781SVaclav Hapla . jname - name of dataset j representing column indices 404058bd781SVaclav Hapla . aname - name of dataset a representing matrix values 405058bd781SVaclav Hapla - cname - name of attribute stoting column count 406058bd781SVaclav Hapla 407058bd781SVaclav Hapla Level: advanced 408058bd781SVaclav Hapla 409058bd781SVaclav Hapla Notes: 410bd7fcbe1SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be readily loaded. 411058bd781SVaclav Hapla 412058bd781SVaclav Hapla .seealso: MatLoad(), PetscViewerCreate(), PetscViewerSetType(), PETSCVIEWERHDF5, PetscViewerHDF5SetAIJNames() 413058bd781SVaclav Hapla @*/ 414bd7fcbe1SVaclav Hapla /* TODO Once corresponding MatView is implemented, add this: 415bd7fcbe1SVaclav Hapla Current defaults are (iname, jname, aname, cname) = ("i", "j", "a", "ncols"). 416bd7fcbe1SVaclav Hapla For PetscViewerFormat PETSC_VIEWER_HDF5_MAT they are ("jc", "ir", "data", "MATLAB_sparse") so that MAT files can be loaded. 417bd7fcbe1SVaclav Hapla */ 418058bd781SVaclav Hapla PetscErrorCode PetscViewerHDF5GetAIJNames(PetscViewer viewer, const char *iname[], const char *jname[], const char *aname[], const char *cname[]) 419058bd781SVaclav Hapla { 420058bd781SVaclav Hapla PetscErrorCode ierr; 421058bd781SVaclav Hapla 422058bd781SVaclav Hapla PetscFunctionBegin; 423058bd781SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 424058bd781SVaclav Hapla PetscValidPointer(iname,2); 425058bd781SVaclav Hapla PetscValidPointer(jname,3); 426058bd781SVaclav Hapla PetscValidPointer(aname,4); 427058bd781SVaclav Hapla PetscValidPointer(cname,5); 428058bd781SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetAIJNames_C",(PetscViewer,const char*[],const char*[],const char*[],const char*[]),(viewer,iname,jname,aname,cname));CHKERRQ(ierr); 429058bd781SVaclav Hapla PetscFunctionReturn(0); 430058bd781SVaclav Hapla } 431058bd781SVaclav Hapla 432b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 433b723ab35SVaclav Hapla { 434b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 435b723ab35SVaclav Hapla PetscErrorCode ierr; 436b723ab35SVaclav Hapla 437b723ab35SVaclav Hapla PetscFunctionBegin; 438b723ab35SVaclav Hapla if (!hdf5->mataij_names_set) { 439bd7fcbe1SVaclav Hapla /* TODO Once corresponding MatView is implemented, uncomment */ 440bd7fcbe1SVaclav Hapla #if 0 441cbb4c999SVaclav Hapla if (viewer->format == PETSC_VIEWER_HDF5_MAT) { 442bd7fcbe1SVaclav Hapla #endif 443b723ab35SVaclav Hapla ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"jc","ir","data","MATLAB_sparse");CHKERRQ(ierr); 444bd7fcbe1SVaclav Hapla #if 0 445cbb4c999SVaclav Hapla } else { 446cbb4c999SVaclav Hapla ierr = PetscViewerHDF5SetAIJNames_HDF5(viewer,"i","j","a","ncols");CHKERRQ(ierr); 447cbb4c999SVaclav Hapla } 448bd7fcbe1SVaclav Hapla #endif 449b723ab35SVaclav Hapla } 450b723ab35SVaclav Hapla PetscFunctionReturn(0); 451b723ab35SVaclav Hapla } 452b723ab35SVaclav Hapla 4538556b5ebSBarry Smith /*MC 4548556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4558556b5ebSBarry Smith 4568556b5ebSBarry Smith 4578556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4588556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4598556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4608556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4618556b5ebSBarry Smith 4621b266c99SBarry Smith Level: beginner 4638556b5ebSBarry Smith M*/ 464d1232d7fSToby Isaac 4658cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4665c6c1daeSBarry Smith { 4675c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4685c6c1daeSBarry Smith PetscErrorCode ierr; 4695c6c1daeSBarry Smith 4705c6c1daeSBarry Smith PetscFunctionBegin; 471b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4725c6c1daeSBarry Smith 4735c6c1daeSBarry Smith v->data = (void*) hdf5; 4745c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 47582ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 476b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4775c6c1daeSBarry Smith v->ops->flush = 0; 4785c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4795c6c1daeSBarry Smith hdf5->filename = 0; 4805c6c1daeSBarry Smith hdf5->timestep = -1; 4810298fd71SBarry Smith hdf5->groups = NULL; 4825c6c1daeSBarry Smith 483b723ab35SVaclav Hapla hdf5->mataij_iname = NULL; 484b723ab35SVaclav Hapla hdf5->mataij_jname = NULL; 485b723ab35SVaclav Hapla hdf5->mataij_aname = NULL; 486b723ab35SVaclav Hapla hdf5->mataij_cname = NULL; 487b723ab35SVaclav Hapla hdf5->mataij_names_set = PETSC_FALSE; 488058bd781SVaclav Hapla 4890b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 490d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4910b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4920b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 49382ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4949a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 495058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetAIJNames_C",PetscViewerHDF5SetAIJNames_HDF5);CHKERRQ(ierr); 496058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetAIJNames_C",PetscViewerHDF5GetAIJNames_HDF5);CHKERRQ(ierr); 4975c6c1daeSBarry Smith PetscFunctionReturn(0); 4985c6c1daeSBarry Smith } 4995c6c1daeSBarry Smith 5005c6c1daeSBarry Smith /*@C 5015c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 5025c6c1daeSBarry Smith 5035c6c1daeSBarry Smith Collective on MPI_Comm 5045c6c1daeSBarry Smith 5055c6c1daeSBarry Smith Input Parameters: 5065c6c1daeSBarry Smith + comm - MPI communicator 5075c6c1daeSBarry Smith . name - name of file 5085c6c1daeSBarry Smith - type - type of file 5095c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 5105c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 5115c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 5125c6c1daeSBarry Smith 5135c6c1daeSBarry Smith Output Parameter: 5145c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5155c6c1daeSBarry Smith 51682ea9b62SBarry Smith Options Database: 51782ea9b62SBarry 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 5189a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 51982ea9b62SBarry Smith 5205c6c1daeSBarry Smith Level: beginner 5215c6c1daeSBarry Smith 5225c6c1daeSBarry Smith Note: 5235c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5245c6c1daeSBarry Smith 5255c6c1daeSBarry Smith Concepts: HDF5 files 5265c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 5275c6c1daeSBarry Smith 5286a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5299a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 530a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5315c6c1daeSBarry Smith @*/ 5325c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5335c6c1daeSBarry Smith { 5345c6c1daeSBarry Smith PetscErrorCode ierr; 5355c6c1daeSBarry Smith 5365c6c1daeSBarry Smith PetscFunctionBegin; 5375c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5385c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5395c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5405c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 5415c6c1daeSBarry Smith PetscFunctionReturn(0); 5425c6c1daeSBarry Smith } 5435c6c1daeSBarry Smith 5445c6c1daeSBarry Smith /*@C 5455c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5465c6c1daeSBarry Smith 5475c6c1daeSBarry Smith Not collective 5485c6c1daeSBarry Smith 5495c6c1daeSBarry Smith Input Parameter: 5505c6c1daeSBarry Smith . viewer - the PetscViewer 5515c6c1daeSBarry Smith 5525c6c1daeSBarry Smith Output Parameter: 5535c6c1daeSBarry Smith . file_id - The file id 5545c6c1daeSBarry Smith 5555c6c1daeSBarry Smith Level: intermediate 5565c6c1daeSBarry Smith 5575c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5585c6c1daeSBarry Smith @*/ 5595c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5605c6c1daeSBarry Smith { 5615c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5625c6c1daeSBarry Smith 5635c6c1daeSBarry Smith PetscFunctionBegin; 5645c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5655c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5665c6c1daeSBarry Smith PetscFunctionReturn(0); 5675c6c1daeSBarry Smith } 5685c6c1daeSBarry Smith 5695c6c1daeSBarry Smith /*@C 5705c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5715c6c1daeSBarry Smith 5725c6c1daeSBarry Smith Not collective 5735c6c1daeSBarry Smith 5745c6c1daeSBarry Smith Input Parameters: 5755c6c1daeSBarry Smith + viewer - the PetscViewer 5765c6c1daeSBarry Smith - name - The group name 5775c6c1daeSBarry Smith 5785c6c1daeSBarry Smith Level: intermediate 5795c6c1daeSBarry Smith 580874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5815c6c1daeSBarry Smith @*/ 5825c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5835c6c1daeSBarry Smith { 5845c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5855c6c1daeSBarry Smith GroupList *groupNode; 5865c6c1daeSBarry Smith PetscErrorCode ierr; 5875c6c1daeSBarry Smith 5885c6c1daeSBarry Smith PetscFunctionBegin; 5895c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5905c6c1daeSBarry Smith PetscValidCharPointer(name,2); 59195dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5925c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 593a297a907SKarl Rupp 5945c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5955c6c1daeSBarry Smith hdf5->groups = groupNode; 5965c6c1daeSBarry Smith PetscFunctionReturn(0); 5975c6c1daeSBarry Smith } 5985c6c1daeSBarry Smith 5993ef9c667SSatish Balay /*@ 6005c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6015c6c1daeSBarry Smith 6025c6c1daeSBarry Smith Not collective 6035c6c1daeSBarry Smith 6045c6c1daeSBarry Smith Input Parameter: 6055c6c1daeSBarry Smith . viewer - the PetscViewer 6065c6c1daeSBarry Smith 6075c6c1daeSBarry Smith Level: intermediate 6085c6c1daeSBarry Smith 609874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6105c6c1daeSBarry Smith @*/ 6115c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6125c6c1daeSBarry Smith { 6135c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6145c6c1daeSBarry Smith GroupList *groupNode; 6155c6c1daeSBarry Smith PetscErrorCode ierr; 6165c6c1daeSBarry Smith 6175c6c1daeSBarry Smith PetscFunctionBegin; 6185c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 61982f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6205c6c1daeSBarry Smith groupNode = hdf5->groups; 6215c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6225c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6235c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6245c6c1daeSBarry Smith PetscFunctionReturn(0); 6255c6c1daeSBarry Smith } 6265c6c1daeSBarry Smith 6275c6c1daeSBarry Smith /*@C 628874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 629874270d9SVaclav Hapla If none has been assigned, returns NULL. 6305c6c1daeSBarry Smith 6315c6c1daeSBarry Smith Not collective 6325c6c1daeSBarry Smith 6335c6c1daeSBarry Smith Input Parameter: 6345c6c1daeSBarry Smith . viewer - the PetscViewer 6355c6c1daeSBarry Smith 6365c6c1daeSBarry Smith Output Parameter: 6375c6c1daeSBarry Smith . name - The group name 6385c6c1daeSBarry Smith 6395c6c1daeSBarry Smith Level: intermediate 6405c6c1daeSBarry Smith 641874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6425c6c1daeSBarry Smith @*/ 6435c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 6445c6c1daeSBarry Smith { 6455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6465c6c1daeSBarry Smith 6475c6c1daeSBarry Smith PetscFunctionBegin; 6485c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 649c959eef4SJed Brown PetscValidPointer(name,2); 650a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6510298fd71SBarry Smith else *name = NULL; 6525c6c1daeSBarry Smith PetscFunctionReturn(0); 6535c6c1daeSBarry Smith } 6545c6c1daeSBarry Smith 6558c557b5aSVaclav Hapla /*@ 656874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 657874270d9SVaclav Hapla and return this group's ID and file ID. 658874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 659874270d9SVaclav Hapla 660874270d9SVaclav Hapla Not collective 661874270d9SVaclav Hapla 662874270d9SVaclav Hapla Input Parameter: 663874270d9SVaclav Hapla . viewer - the PetscViewer 664874270d9SVaclav Hapla 665874270d9SVaclav Hapla Output Parameter: 666874270d9SVaclav Hapla + fileId - The HDF5 file ID 667874270d9SVaclav Hapla - groupId - The HDF5 group ID 668874270d9SVaclav Hapla 669e74616beSVaclav Hapla Notes: 670e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 671e74616beSVaclav Hapla 672874270d9SVaclav Hapla Level: intermediate 673874270d9SVaclav Hapla 674874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 675874270d9SVaclav Hapla @*/ 67654dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 67754dbf706SVaclav Hapla { 67890e03537SVaclav Hapla hid_t file_id; 67990e03537SVaclav Hapla H5O_type_t type; 68054dbf706SVaclav Hapla const char *groupName = NULL; 681e74616beSVaclav Hapla PetscBool create; 68254dbf706SVaclav Hapla PetscErrorCode ierr; 68354dbf706SVaclav Hapla 68454dbf706SVaclav Hapla PetscFunctionBegin; 685e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 68654dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 68754dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 688e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 68990e03537SVaclav 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); 69090e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 69154dbf706SVaclav Hapla *fileId = file_id; 69254dbf706SVaclav Hapla PetscFunctionReturn(0); 69354dbf706SVaclav Hapla } 69454dbf706SVaclav Hapla 6953ef9c667SSatish Balay /*@ 6965c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6975c6c1daeSBarry Smith 6985c6c1daeSBarry Smith Not collective 6995c6c1daeSBarry Smith 7005c6c1daeSBarry Smith Input Parameter: 7015c6c1daeSBarry Smith . viewer - the PetscViewer 7025c6c1daeSBarry Smith 7035c6c1daeSBarry Smith Level: intermediate 7045c6c1daeSBarry Smith 7055c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 7065c6c1daeSBarry Smith @*/ 7075c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 7085c6c1daeSBarry Smith { 7095c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7105c6c1daeSBarry Smith 7115c6c1daeSBarry Smith PetscFunctionBegin; 7125c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7135c6c1daeSBarry Smith ++hdf5->timestep; 7145c6c1daeSBarry Smith PetscFunctionReturn(0); 7155c6c1daeSBarry Smith } 7165c6c1daeSBarry Smith 7173ef9c667SSatish Balay /*@ 7185c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 7195c6c1daeSBarry Smith of -1 disables blocking with timesteps. 7205c6c1daeSBarry Smith 7215c6c1daeSBarry Smith Not collective 7225c6c1daeSBarry Smith 7235c6c1daeSBarry Smith Input Parameters: 7245c6c1daeSBarry Smith + viewer - the PetscViewer 7255c6c1daeSBarry Smith - timestep - The timestep number 7265c6c1daeSBarry Smith 7275c6c1daeSBarry Smith Level: intermediate 7285c6c1daeSBarry Smith 7295c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 7305c6c1daeSBarry Smith @*/ 7315c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 7325c6c1daeSBarry Smith { 7335c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7345c6c1daeSBarry Smith 7355c6c1daeSBarry Smith PetscFunctionBegin; 7365c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7375c6c1daeSBarry Smith hdf5->timestep = timestep; 7385c6c1daeSBarry Smith PetscFunctionReturn(0); 7395c6c1daeSBarry Smith } 7405c6c1daeSBarry Smith 7413ef9c667SSatish Balay /*@ 7425c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7435c6c1daeSBarry Smith 7445c6c1daeSBarry Smith Not collective 7455c6c1daeSBarry Smith 7465c6c1daeSBarry Smith Input Parameter: 7475c6c1daeSBarry Smith . viewer - the PetscViewer 7485c6c1daeSBarry Smith 7495c6c1daeSBarry Smith Output Parameter: 7505c6c1daeSBarry Smith . timestep - The timestep number 7515c6c1daeSBarry Smith 7525c6c1daeSBarry Smith Level: intermediate 7535c6c1daeSBarry Smith 7545c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7555c6c1daeSBarry Smith @*/ 7565c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7575c6c1daeSBarry Smith { 7585c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7595c6c1daeSBarry Smith 7605c6c1daeSBarry Smith PetscFunctionBegin; 7615c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7625c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7635c6c1daeSBarry Smith *timestep = hdf5->timestep; 7645c6c1daeSBarry Smith PetscFunctionReturn(0); 7655c6c1daeSBarry Smith } 7665c6c1daeSBarry Smith 76736bce6e8SMatthew G. Knepley /*@C 76836bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 76936bce6e8SMatthew G. Knepley 77036bce6e8SMatthew G. Knepley Not collective 77136bce6e8SMatthew G. Knepley 77236bce6e8SMatthew G. Knepley Input Parameter: 77336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 77436bce6e8SMatthew G. Knepley 77536bce6e8SMatthew G. Knepley Output Parameter: 77636bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 77736bce6e8SMatthew G. Knepley 77836bce6e8SMatthew G. Knepley Level: advanced 77936bce6e8SMatthew G. Knepley 78036bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 78136bce6e8SMatthew G. Knepley @*/ 78236bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 78336bce6e8SMatthew G. Knepley { 78436bce6e8SMatthew G. Knepley PetscFunctionBegin; 78536bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 78636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 78736bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 78836bce6e8SMatthew G. Knepley #else 78936bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 79036bce6e8SMatthew G. Knepley #endif 79136bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 79236bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 79336bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 794c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 795de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 79636bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 79736bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 79836bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7997e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 80036bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 80136bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 80236bce6e8SMatthew G. Knepley } 80336bce6e8SMatthew G. Knepley 80436bce6e8SMatthew G. Knepley /*@C 80536bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 80636bce6e8SMatthew G. Knepley 80736bce6e8SMatthew G. Knepley Not collective 80836bce6e8SMatthew G. Knepley 80936bce6e8SMatthew G. Knepley Input Parameter: 81036bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 81136bce6e8SMatthew G. Knepley 81236bce6e8SMatthew G. Knepley Output Parameter: 81336bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 81436bce6e8SMatthew G. Knepley 81536bce6e8SMatthew G. Knepley Level: advanced 81636bce6e8SMatthew G. Knepley 81736bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 81836bce6e8SMatthew G. Knepley @*/ 81936bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 82036bce6e8SMatthew G. Knepley { 82136bce6e8SMatthew G. Knepley PetscFunctionBegin; 82236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 82336bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 82436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 82536bce6e8SMatthew G. Knepley #else 82636bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 82736bce6e8SMatthew G. Knepley #endif 82836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 82936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 83036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 83136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 83236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 83336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 8347e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 83536bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 83636bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 83736bce6e8SMatthew G. Knepley } 83836bce6e8SMatthew G. Knepley 839df863907SAlex Fikl /*@C 840b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 84136bce6e8SMatthew G. Knepley 84236bce6e8SMatthew G. Knepley Input Parameters: 84336bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 84457d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 84536bce6e8SMatthew G. Knepley . name - The attribute name 84636bce6e8SMatthew G. Knepley . datatype - The attribute type 84736bce6e8SMatthew G. Knepley - value - The attribute value 84836bce6e8SMatthew G. Knepley 84936bce6e8SMatthew G. Knepley Level: advanced 85036bce6e8SMatthew G. Knepley 851577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 85236bce6e8SMatthew G. Knepley @*/ 85357d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 85436bce6e8SMatthew G. Knepley { 85557d22f7dSVaclav Hapla char *parent; 85660568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 857080f144cSVaclav Hapla PetscBool has; 85836bce6e8SMatthew G. Knepley PetscErrorCode ierr; 85936bce6e8SMatthew G. Knepley 86036bce6e8SMatthew G. Knepley PetscFunctionBegin; 8615cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 86257d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 863c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 864b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 86557d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 866bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 867b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 86836bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8697e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8707e97c476SMatthew G. Knepley size_t len; 8717e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 872729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8737e97c476SMatthew G. Knepley } 87436bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 875729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 87660568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 877080f144cSVaclav Hapla if (has) { 878080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 879080f144cSVaclav Hapla } else { 88060568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 881080f144cSVaclav Hapla } 882729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 883729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 884729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 88560568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 886729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 88757d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 88836bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 88936bce6e8SMatthew G. Knepley } 89036bce6e8SMatthew G. Knepley 891df863907SAlex Fikl /*@C 892577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 893577e0e04SVaclav Hapla 894577e0e04SVaclav Hapla Input Parameters: 895577e0e04SVaclav Hapla + viewer - The HDF5 viewer 896577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 897577e0e04SVaclav Hapla . name - The attribute name 898577e0e04SVaclav Hapla . datatype - The attribute type 899577e0e04SVaclav Hapla - value - The attribute value 900577e0e04SVaclav Hapla 901577e0e04SVaclav Hapla Notes: 902577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 903577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 904577e0e04SVaclav Hapla 905577e0e04SVaclav Hapla Level: advanced 906577e0e04SVaclav Hapla 907577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 908577e0e04SVaclav Hapla @*/ 909577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 910577e0e04SVaclav Hapla { 911577e0e04SVaclav Hapla PetscErrorCode ierr; 912577e0e04SVaclav Hapla 913577e0e04SVaclav Hapla PetscFunctionBegin; 914577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 915577e0e04SVaclav Hapla PetscValidHeader(obj,2); 916577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 917b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 918577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 919577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 920577e0e04SVaclav Hapla PetscFunctionReturn(0); 921577e0e04SVaclav Hapla } 922577e0e04SVaclav Hapla 923577e0e04SVaclav Hapla /*@C 924b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 925ad0c61feSMatthew G. Knepley 926ad0c61feSMatthew G. Knepley Input Parameters: 927ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 92857d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 929ad0c61feSMatthew G. Knepley . name - The attribute name 930ad0c61feSMatthew G. Knepley - datatype - The attribute type 931ad0c61feSMatthew G. Knepley 932ad0c61feSMatthew G. Knepley Output Parameter: 933ad0c61feSMatthew G. Knepley . value - The attribute value 934ad0c61feSMatthew G. Knepley 935ad0c61feSMatthew G. Knepley Level: advanced 936ad0c61feSMatthew G. Knepley 937577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 938ad0c61feSMatthew G. Knepley @*/ 93957d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 940ad0c61feSMatthew G. Knepley { 94157d22f7dSVaclav Hapla char *parent; 942f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 943969748fdSVaclav Hapla PetscBool has; 944ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 945ad0c61feSMatthew G. Knepley 946ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9475cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 94857d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 949c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 950b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 95157d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 952969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 953969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 954969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 955ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 956ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 95760568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 95860568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 959f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 960f0b58479SMatthew G. Knepley size_t len; 961e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 96215b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 963f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 964f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 965f0b58479SMatthew G. Knepley } 96670efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 967729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 968e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 969e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 97057d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 971ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 972ad0c61feSMatthew G. Knepley } 973ad0c61feSMatthew G. Knepley 974577e0e04SVaclav Hapla /*@C 975577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 976577e0e04SVaclav Hapla 977577e0e04SVaclav Hapla Input Parameters: 978577e0e04SVaclav Hapla + viewer - The HDF5 viewer 979577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 980577e0e04SVaclav Hapla . name - The attribute name 981577e0e04SVaclav Hapla - datatype - The attribute type 982577e0e04SVaclav Hapla 983577e0e04SVaclav Hapla Output Parameter: 984577e0e04SVaclav Hapla . value - The attribute value 985577e0e04SVaclav Hapla 986577e0e04SVaclav Hapla Notes: 987577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 988577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 989577e0e04SVaclav Hapla 990577e0e04SVaclav Hapla Level: advanced 991577e0e04SVaclav Hapla 992577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 993577e0e04SVaclav Hapla @*/ 994577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 995577e0e04SVaclav Hapla { 996577e0e04SVaclav Hapla PetscErrorCode ierr; 997577e0e04SVaclav Hapla 998577e0e04SVaclav Hapla PetscFunctionBegin; 999577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1000577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1001577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1002b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 1003577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1004577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1005577e0e04SVaclav Hapla PetscFunctionReturn(0); 1006577e0e04SVaclav Hapla } 1007577e0e04SVaclav Hapla 1008a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1009a75c3fd4SVaclav Hapla { 1010a75c3fd4SVaclav Hapla htri_t exists; 1011a75c3fd4SVaclav Hapla hid_t group; 1012a75c3fd4SVaclav Hapla 1013a75c3fd4SVaclav Hapla PetscFunctionBegin; 1014a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1015a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1016a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1017a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1018a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1019a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1020a75c3fd4SVaclav Hapla } 1021a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1022a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1023a75c3fd4SVaclav Hapla } 1024a75c3fd4SVaclav Hapla 1025bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 10265cdeb108SMatthew G. Knepley { 102790e03537SVaclav Hapla const char rootGroupName[] = "/"; 10285cdeb108SMatthew G. Knepley hid_t h5; 1029e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 103015b861d2SVaclav Hapla PetscInt i; 103115b861d2SVaclav Hapla int n; 103285688ae2SVaclav Hapla char **hierarchy; 103385688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 10345cdeb108SMatthew G. Knepley PetscErrorCode ierr; 10355cdeb108SMatthew G. Knepley 10365cdeb108SMatthew G. Knepley PetscFunctionBegin; 10375cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 103890e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 103990e03537SVaclav Hapla else name = rootGroupName; 1040ccd66a83SVaclav Hapla if (has) { 104156cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10425cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1043ccd66a83SVaclav Hapla } 1044ccd66a83SVaclav Hapla if (otype) { 1045ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 104656cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1047ccd66a83SVaclav Hapla } 1048c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 104985688ae2SVaclav Hapla 105085688ae2SVaclav Hapla /* 105185688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 105285688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 105385688ae2SVaclav Hapla 1) whether it's a valid link 105485688ae2SVaclav Hapla 2) whether this link resolves to an object 105585688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 105685688ae2SVaclav Hapla */ 105785688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 105885688ae2SVaclav Hapla if (!n) { 105985688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1060ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1061ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 106285688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 106385688ae2SVaclav Hapla PetscFunctionReturn(0); 106485688ae2SVaclav Hapla } 106585688ae2SVaclav Hapla for (i=0; i<n; i++) { 106685688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 106785688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1068a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1069a75c3fd4SVaclav Hapla if (!exists) break; 107085688ae2SVaclav Hapla } 107185688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 107285688ae2SVaclav Hapla 107385688ae2SVaclav Hapla /* If the object exists, get its type */ 1074ccd66a83SVaclav Hapla if (exists && otype) { 10755cdeb108SMatthew G. Knepley H5O_info_t info; 10765cdeb108SMatthew G. Knepley 107776276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 107804633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 107956cc0592SVaclav Hapla *otype = info.type; 10805cdeb108SMatthew G. Knepley } 1081ccd66a83SVaclav Hapla if (has) *has = exists; 10825cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 10835cdeb108SMatthew G. Knepley } 10845cdeb108SMatthew G. Knepley 1085ecce1506SVaclav Hapla /*@ 10860a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 10870a9f38efSVaclav Hapla 10880a9f38efSVaclav Hapla Input Parameters: 10890a9f38efSVaclav Hapla . viewer - The HDF5 viewer 10900a9f38efSVaclav Hapla 10910a9f38efSVaclav Hapla Output Parameter: 10920a9f38efSVaclav Hapla . has - Flag for group existence 10930a9f38efSVaclav Hapla 10940a9f38efSVaclav Hapla Notes: 10950a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 10960a9f38efSVaclav Hapla 10970a9f38efSVaclav Hapla Level: advanced 10980a9f38efSVaclav Hapla 10990a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 11000a9f38efSVaclav Hapla @*/ 11010a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 11020a9f38efSVaclav Hapla { 11030a9f38efSVaclav Hapla H5O_type_t type; 11040a9f38efSVaclav Hapla const char *name; 11050a9f38efSVaclav Hapla PetscErrorCode ierr; 11060a9f38efSVaclav Hapla 11070a9f38efSVaclav Hapla PetscFunctionBegin; 11080a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11090a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 11100a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 11110a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 11120a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 11130a9f38efSVaclav Hapla PetscFunctionReturn(0); 11140a9f38efSVaclav Hapla } 11150a9f38efSVaclav Hapla 11160a9f38efSVaclav Hapla /*@ 1117e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1118ecce1506SVaclav Hapla 1119ecce1506SVaclav Hapla Input Parameters: 1120ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1121ecce1506SVaclav Hapla - obj - The named object 1122ecce1506SVaclav Hapla 1123ecce1506SVaclav Hapla Output Parameter: 1124ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1125ecce1506SVaclav Hapla 1126e3f143f7SVaclav Hapla Notes: 1127e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1128e3f143f7SVaclav Hapla 1129ecce1506SVaclav Hapla Level: advanced 1130ecce1506SVaclav Hapla 1131e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1132ecce1506SVaclav Hapla @*/ 1133ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1134ecce1506SVaclav Hapla { 113556cc0592SVaclav Hapla H5O_type_t type; 11366c132bc1SVaclav Hapla char *path; 1137ecce1506SVaclav Hapla PetscErrorCode ierr; 1138ecce1506SVaclav Hapla 1139ecce1506SVaclav Hapla PetscFunctionBegin; 1140c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1141c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1142c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1143ecce1506SVaclav Hapla *has = PETSC_FALSE; 1144e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11456c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11466c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 114756cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11486c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1149ecce1506SVaclav Hapla PetscFunctionReturn(0); 1150ecce1506SVaclav Hapla } 1151ecce1506SVaclav Hapla 1152df863907SAlex Fikl /*@C 1153ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1154e7bdbf8eSMatthew G. Knepley 1155e7bdbf8eSMatthew G. Knepley Input Parameters: 1156e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 115710e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1158e7bdbf8eSMatthew G. Knepley - name - The attribute name 1159e7bdbf8eSMatthew G. Knepley 1160e7bdbf8eSMatthew G. Knepley Output Parameter: 1161e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1162e7bdbf8eSMatthew G. Knepley 1163e7bdbf8eSMatthew G. Knepley Level: advanced 1164e7bdbf8eSMatthew G. Knepley 1165577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1166e7bdbf8eSMatthew G. Knepley @*/ 116710e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1168e7bdbf8eSMatthew G. Knepley { 116910e69e8fSVaclav Hapla char *parent; 1170e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1171e7bdbf8eSMatthew G. Knepley 1172e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 11735cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 117410e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1175c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1176c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 117710e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1178bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 117910e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 118010e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 118106db490cSVaclav Hapla PetscFunctionReturn(0); 118206db490cSVaclav Hapla } 118306db490cSVaclav Hapla 1184577e0e04SVaclav Hapla /*@C 1185577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1186577e0e04SVaclav Hapla 1187577e0e04SVaclav Hapla Input Parameters: 1188577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1189577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1190577e0e04SVaclav Hapla - name - The attribute name 1191577e0e04SVaclav Hapla 1192577e0e04SVaclav Hapla Output Parameter: 1193577e0e04SVaclav Hapla . has - Flag for attribute existence 1194577e0e04SVaclav Hapla 1195577e0e04SVaclav Hapla Notes: 1196577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1197577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1198577e0e04SVaclav Hapla 1199577e0e04SVaclav Hapla Level: advanced 1200577e0e04SVaclav Hapla 1201577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1202577e0e04SVaclav Hapla @*/ 1203577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1204577e0e04SVaclav Hapla { 1205577e0e04SVaclav Hapla PetscErrorCode ierr; 1206577e0e04SVaclav Hapla 1207577e0e04SVaclav Hapla PetscFunctionBegin; 1208577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1209577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1210577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1211577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1212577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1213577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1214577e0e04SVaclav Hapla PetscFunctionReturn(0); 1215577e0e04SVaclav Hapla } 1216577e0e04SVaclav Hapla 121706db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 121806db490cSVaclav Hapla { 121906db490cSVaclav Hapla hid_t h5; 122006db490cSVaclav Hapla htri_t hhas; 122106db490cSVaclav Hapla PetscErrorCode ierr; 122206db490cSVaclav Hapla 122306db490cSVaclav Hapla PetscFunctionBegin; 122406db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 12252f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 12265cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1227e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1228e7bdbf8eSMatthew G. Knepley } 1229e7bdbf8eSMatthew G. Knepley 1230a75e6a4aSMatthew G. Knepley /* 1231a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1232a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1233a75e6a4aSMatthew G. Knepley */ 1234d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1235a75e6a4aSMatthew G. Knepley 1236a75e6a4aSMatthew G. Knepley /*@C 1237a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1238a75e6a4aSMatthew G. Knepley 1239a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1240a75e6a4aSMatthew G. Knepley 1241a75e6a4aSMatthew G. Knepley Input Parameter: 1242a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1243a75e6a4aSMatthew G. Knepley 1244a75e6a4aSMatthew G. Knepley Level: intermediate 1245a75e6a4aSMatthew G. Knepley 1246a75e6a4aSMatthew G. Knepley Options Database Keys: 1247a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1248a75e6a4aSMatthew G. Knepley 1249a75e6a4aSMatthew G. Knepley Environmental variables: 1250a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1251a75e6a4aSMatthew G. Knepley 1252a75e6a4aSMatthew G. Knepley Notes: 1253a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1254a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1255a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1256a75e6a4aSMatthew G. Knepley 1257a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1258a75e6a4aSMatthew G. Knepley @*/ 1259a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1260a75e6a4aSMatthew G. Knepley { 1261a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1262a75e6a4aSMatthew G. Knepley PetscBool flg; 1263a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1264a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1265a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1266a75e6a4aSMatthew G. Knepley 1267a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1268a75e6a4aSMatthew 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);} 1269a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 127012801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1271a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1272a75e6a4aSMatthew G. Knepley } 127347435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1274a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1275a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1276a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1277a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1278a75e6a4aSMatthew G. Knepley if (!flg) { 1279a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1280a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1281a75e6a4aSMatthew G. Knepley } 1282a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1283a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1284a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1285a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 128647435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1287a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1288a75e6a4aSMatthew G. Knepley } 1289a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1290a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1291a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1292a75e6a4aSMatthew G. Knepley } 1293