1bc307625SVaclav Hapla #include <petsc/private/viewerimpl.h> 27d964842SVaclav Hapla #include <petsc/private/viewerhdf5impl.h> 3d70abbfaSBarry Smith #include <petscviewerhdf5.h> /*I "petscviewerhdf5.h" I*/ 45c6c1daeSBarry Smith 5bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 606db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 706db490cSVaclav Hapla 86c132bc1SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char objname[], char **fullpath) 96c132bc1SVaclav Hapla { 106c132bc1SVaclav Hapla const char *group; 116c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 126c132bc1SVaclav Hapla PetscErrorCode ierr; 136c132bc1SVaclav Hapla 146c132bc1SVaclav Hapla PetscFunctionBegin; 156c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 166c132bc1SVaclav Hapla ierr = PetscStrcat(buf, group);CHKERRQ(ierr); 176c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 186c132bc1SVaclav Hapla ierr = PetscStrcat(buf, objname);CHKERRQ(ierr); 196c132bc1SVaclav Hapla ierr = PetscStrallocpy(buf, fullpath);CHKERRQ(ierr); 206c132bc1SVaclav Hapla PetscFunctionReturn(0); 216c132bc1SVaclav Hapla } 226c132bc1SVaclav Hapla 23577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 24577e0e04SVaclav Hapla { 25577e0e04SVaclav Hapla PetscBool has; 26577e0e04SVaclav Hapla const char *group; 27577e0e04SVaclav Hapla PetscErrorCode ierr; 28577e0e04SVaclav Hapla 29577e0e04SVaclav Hapla PetscFunctionBegin; 30577e0e04SVaclav Hapla if (!obj->name) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 31577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 32577e0e04SVaclav Hapla if (!has) { 33577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 34577e0e04SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) %s not stored in group %s", obj->name, group); 35577e0e04SVaclav Hapla } 36577e0e04SVaclav Hapla PetscFunctionReturn(0); 37577e0e04SVaclav Hapla } 38577e0e04SVaclav Hapla 394416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4082ea9b62SBarry Smith { 4182ea9b62SBarry Smith PetscErrorCode ierr; 42*ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4482ea9b62SBarry Smith 4582ea9b62SBarry Smith PetscFunctionBegin; 4682ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 4782ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 489a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 49*ee8b9732SVaclav Hapla ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 50*ee8b9732SVaclav Hapla if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 5182ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5282ea9b62SBarry Smith PetscFunctionReturn(0); 5382ea9b62SBarry Smith } 5482ea9b62SBarry Smith 551b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 561b793a25SVaclav Hapla { 571b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 581b793a25SVaclav Hapla PetscErrorCode ierr; 591b793a25SVaclav Hapla 601b793a25SVaclav Hapla PetscFunctionBegin; 611b793a25SVaclav Hapla if (hdf5->filename) { 621b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 631b793a25SVaclav Hapla } 641b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 651b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 661b793a25SVaclav Hapla PetscFunctionReturn(0); 671b793a25SVaclav Hapla } 681b793a25SVaclav Hapla 695c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 705c6c1daeSBarry Smith { 715c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 725c6c1daeSBarry Smith PetscErrorCode ierr; 735c6c1daeSBarry Smith 745c6c1daeSBarry Smith PetscFunctionBegin; 755c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 76729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 775c6c1daeSBarry Smith PetscFunctionReturn(0); 785c6c1daeSBarry Smith } 795c6c1daeSBarry Smith 809b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 815c6c1daeSBarry Smith { 825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 835c6c1daeSBarry Smith PetscErrorCode ierr; 845c6c1daeSBarry Smith 855c6c1daeSBarry Smith PetscFunctionBegin; 869c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 875c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 885c6c1daeSBarry Smith while (hdf5->groups) { 897d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 905c6c1daeSBarry Smith 915c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 925c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 935c6c1daeSBarry Smith hdf5->groups = tmp; 945c6c1daeSBarry Smith } 955c6c1daeSBarry Smith ierr = PetscFree(hdf5);CHKERRQ(ierr); 960b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);CHKERRQ(ierr); 97d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);CHKERRQ(ierr); 980b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);CHKERRQ(ierr); 99058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL);CHKERRQ(ierr); 100058bd781SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL);CHKERRQ(ierr); 1015c6c1daeSBarry Smith PetscFunctionReturn(0); 1025c6c1daeSBarry Smith } 1035c6c1daeSBarry Smith 1049b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1055c6c1daeSBarry Smith { 1065c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1075c6c1daeSBarry Smith 1085c6c1daeSBarry Smith PetscFunctionBegin; 1095c6c1daeSBarry Smith hdf5->btype = type; 1105c6c1daeSBarry Smith PetscFunctionReturn(0); 1115c6c1daeSBarry Smith } 1125c6c1daeSBarry Smith 1130b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1140b62783dSVaclav Hapla { 1150b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1160b62783dSVaclav Hapla 1170b62783dSVaclav Hapla PetscFunctionBegin; 1180b62783dSVaclav Hapla *type = hdf5->btype; 1190b62783dSVaclav Hapla PetscFunctionReturn(0); 1200b62783dSVaclav Hapla } 1210b62783dSVaclav Hapla 1229b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 12382ea9b62SBarry Smith { 12482ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 12582ea9b62SBarry Smith 12682ea9b62SBarry Smith PetscFunctionBegin; 12782ea9b62SBarry Smith hdf5->basedimension2 = flg; 12882ea9b62SBarry Smith PetscFunctionReturn(0); 12982ea9b62SBarry Smith } 13082ea9b62SBarry Smith 131df863907SAlex Fikl /*@ 13282ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 13382ea9b62SBarry Smith dimension of 2. 13482ea9b62SBarry Smith 13582ea9b62SBarry Smith Logically Collective on PetscViewer 13682ea9b62SBarry Smith 13782ea9b62SBarry Smith Input Parameters: 13882ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 13982ea9b62SBarry 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 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith Options Database: 14282ea9b62SBarry 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 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith 14595452b02SPatrick Sanan Notes: 14695452b02SPatrick 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 14782ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 14882ea9b62SBarry Smith 14982ea9b62SBarry Smith Level: intermediate 15082ea9b62SBarry Smith 15182ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 15282ea9b62SBarry Smith 15382ea9b62SBarry Smith @*/ 15482ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 15582ea9b62SBarry Smith { 15682ea9b62SBarry Smith PetscErrorCode ierr; 15782ea9b62SBarry Smith 15882ea9b62SBarry Smith PetscFunctionBegin; 15982ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16082ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16182ea9b62SBarry Smith PetscFunctionReturn(0); 16282ea9b62SBarry Smith } 16382ea9b62SBarry Smith 164df863907SAlex Fikl /*@ 16582ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 16682ea9b62SBarry Smith dimension of 2. 16782ea9b62SBarry Smith 16882ea9b62SBarry Smith Logically Collective on PetscViewer 16982ea9b62SBarry Smith 17082ea9b62SBarry Smith Input Parameter: 17182ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 17282ea9b62SBarry Smith 17382ea9b62SBarry Smith Output Parameter: 17482ea9b62SBarry 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 17582ea9b62SBarry Smith 17695452b02SPatrick Sanan Notes: 17795452b02SPatrick 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 17882ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 17982ea9b62SBarry Smith 18082ea9b62SBarry Smith Level: intermediate 18182ea9b62SBarry Smith 18282ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 18382ea9b62SBarry Smith 18482ea9b62SBarry Smith @*/ 18582ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 18682ea9b62SBarry Smith { 18782ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 18882ea9b62SBarry Smith 18982ea9b62SBarry Smith PetscFunctionBegin; 19082ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19182ea9b62SBarry Smith *flg = hdf5->basedimension2; 19282ea9b62SBarry Smith PetscFunctionReturn(0); 19382ea9b62SBarry Smith } 19482ea9b62SBarry Smith 1959b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 1969a0502c6SHåkon Strandenes { 1979a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1989a0502c6SHåkon Strandenes 1999a0502c6SHåkon Strandenes PetscFunctionBegin; 2009a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2019a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2029a0502c6SHåkon Strandenes } 2039a0502c6SHåkon Strandenes 204df863907SAlex Fikl /*@ 2059a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2069a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2079a0502c6SHåkon Strandenes 2089a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2099a0502c6SHåkon Strandenes 2109a0502c6SHåkon Strandenes Input Parameters: 2119a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2129a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2139a0502c6SHåkon Strandenes 2149a0502c6SHåkon Strandenes Options Database: 2159a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2169a0502c6SHåkon Strandenes 2179a0502c6SHåkon Strandenes 21895452b02SPatrick Sanan Notes: 21995452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2209a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2219a0502c6SHåkon Strandenes 2229a0502c6SHåkon Strandenes Level: intermediate 2239a0502c6SHåkon Strandenes 2249a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2259a0502c6SHåkon Strandenes PetscReal 2269a0502c6SHåkon Strandenes 2279a0502c6SHåkon Strandenes @*/ 2289a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2299a0502c6SHåkon Strandenes { 2309a0502c6SHåkon Strandenes PetscErrorCode ierr; 2319a0502c6SHåkon Strandenes 2329a0502c6SHåkon Strandenes PetscFunctionBegin; 2339a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2349a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2359a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2369a0502c6SHåkon Strandenes } 2379a0502c6SHåkon Strandenes 238df863907SAlex Fikl /*@ 2399a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2409a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2419a0502c6SHåkon Strandenes 2429a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2439a0502c6SHåkon Strandenes 2449a0502c6SHåkon Strandenes Input Parameter: 2459a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2469a0502c6SHåkon Strandenes 2479a0502c6SHåkon Strandenes Output Parameter: 2489a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2499a0502c6SHåkon Strandenes 25095452b02SPatrick Sanan Notes: 25195452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2529a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2539a0502c6SHåkon Strandenes 2549a0502c6SHåkon Strandenes Level: intermediate 2559a0502c6SHåkon Strandenes 2569a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2579a0502c6SHåkon Strandenes PetscReal 2589a0502c6SHåkon Strandenes 2599a0502c6SHåkon Strandenes @*/ 2609a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2619a0502c6SHåkon Strandenes { 2629a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2639a0502c6SHåkon Strandenes 2649a0502c6SHåkon Strandenes PetscFunctionBegin; 2659a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2669a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2679a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2689a0502c6SHåkon Strandenes } 2699a0502c6SHåkon Strandenes 270*ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 271*ee8b9732SVaclav Hapla { 272*ee8b9732SVaclav Hapla PetscFunctionBegin; 273*ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 274*ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 275*ee8b9732SVaclav Hapla #if H5_VERSION_GE(1,10,3) 276*ee8b9732SVaclav Hapla { 277*ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 278*ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 279*ee8b9732SVaclav Hapla } 280*ee8b9732SVaclav Hapla #else 281*ee8b9732SVaclav Hapla if (flg) { 282*ee8b9732SVaclav Hapla PetscErrorCode ierr; 283*ee8b9732SVaclav Hapla ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3\n");CHKERRQ(ierr); 284*ee8b9732SVaclav Hapla } 285*ee8b9732SVaclav Hapla #endif 286*ee8b9732SVaclav Hapla PetscFunctionReturn(0); 287*ee8b9732SVaclav Hapla } 288*ee8b9732SVaclav Hapla 289*ee8b9732SVaclav Hapla /*@ 290*ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 291*ee8b9732SVaclav Hapla 292*ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 293*ee8b9732SVaclav Hapla 294*ee8b9732SVaclav Hapla Input Parameters: 295*ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 296*ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 297*ee8b9732SVaclav Hapla 298*ee8b9732SVaclav Hapla Options Database: 299*ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 300*ee8b9732SVaclav Hapla 301*ee8b9732SVaclav Hapla Notes: 302*ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 303*ee8b9732SVaclav Hapla However, this works correctly only since HDF5 1.10.3; hence, we ignore this setting for older versions. 304*ee8b9732SVaclav Hapla 305*ee8b9732SVaclav Hapla Developer notes: 306*ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 307*ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 308*ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 309*ee8b9732SVaclav Hapla 310*ee8b9732SVaclav Hapla Level: intermediate 311*ee8b9732SVaclav Hapla 312*ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 313*ee8b9732SVaclav Hapla 314*ee8b9732SVaclav Hapla @*/ 315*ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 316*ee8b9732SVaclav Hapla { 317*ee8b9732SVaclav Hapla PetscErrorCode ierr; 318*ee8b9732SVaclav Hapla 319*ee8b9732SVaclav Hapla PetscFunctionBegin; 320*ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 321*ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 322*ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 323*ee8b9732SVaclav Hapla PetscFunctionReturn(0); 324*ee8b9732SVaclav Hapla } 325*ee8b9732SVaclav Hapla 326*ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 327*ee8b9732SVaclav Hapla { 328*ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 329*ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 330*ee8b9732SVaclav Hapla 331*ee8b9732SVaclav Hapla PetscFunctionBegin; 332*ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 333*ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 334*ee8b9732SVaclav Hapla PetscFunctionReturn(0); 335*ee8b9732SVaclav Hapla } 336*ee8b9732SVaclav Hapla 337*ee8b9732SVaclav Hapla /*@ 338*ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 339*ee8b9732SVaclav Hapla 340*ee8b9732SVaclav Hapla Not Collective 341*ee8b9732SVaclav Hapla 342*ee8b9732SVaclav Hapla Input Parameters: 343*ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 344*ee8b9732SVaclav Hapla 345*ee8b9732SVaclav Hapla Output Parameters: 346*ee8b9732SVaclav Hapla . flg - the flag 347*ee8b9732SVaclav Hapla 348*ee8b9732SVaclav Hapla Level: intermediate 349*ee8b9732SVaclav Hapla 350*ee8b9732SVaclav Hapla Notes: 351*ee8b9732SVaclav Hapla This setting works correctly only since HDF5 1.10.3. For older versions, PETSC_FALSE will be always returned. 352*ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 353*ee8b9732SVaclav Hapla 354*ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 355*ee8b9732SVaclav Hapla 356*ee8b9732SVaclav Hapla @*/ 357*ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 358*ee8b9732SVaclav Hapla { 359*ee8b9732SVaclav Hapla PetscErrorCode ierr; 360*ee8b9732SVaclav Hapla 361*ee8b9732SVaclav Hapla PetscFunctionBegin; 362*ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 363*ee8b9732SVaclav Hapla PetscValidPointer(flg,2); 364*ee8b9732SVaclav Hapla 365*ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 366*ee8b9732SVaclav Hapla PetscFunctionReturn(0); 367*ee8b9732SVaclav Hapla } 368*ee8b9732SVaclav Hapla 3699b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3705c6c1daeSBarry Smith { 3715c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3725c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 3735c6c1daeSBarry Smith MPI_Info info = MPI_INFO_NULL; 3745c6c1daeSBarry Smith #endif 3755c6c1daeSBarry Smith hid_t plist_id; 3765c6c1daeSBarry Smith PetscErrorCode ierr; 3775c6c1daeSBarry Smith 3785c6c1daeSBarry Smith PetscFunctionBegin; 379f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 380f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3815c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3825c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 383729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 3845c6c1daeSBarry Smith #if defined(PETSC_HAVE_H5PSET_FAPL_MPIO) 385729ab6d9SBarry Smith PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), info)); 3865c6c1daeSBarry Smith #endif 3875c6c1daeSBarry Smith /* Create or open the file collectively */ 3885c6c1daeSBarry Smith switch (hdf5->btype) { 3895c6c1daeSBarry Smith case FILE_MODE_READ: 390729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 3915c6c1daeSBarry Smith break; 3925c6c1daeSBarry Smith case FILE_MODE_APPEND: 393729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 3945c6c1daeSBarry Smith break; 3955c6c1daeSBarry Smith case FILE_MODE_WRITE: 396729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 3975c6c1daeSBarry Smith break; 3985c6c1daeSBarry Smith default: 3995c6c1daeSBarry Smith SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4005c6c1daeSBarry Smith } 4015c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 402729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4035c6c1daeSBarry Smith PetscFunctionReturn(0); 4045c6c1daeSBarry Smith } 4055c6c1daeSBarry Smith 406d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 407d1232d7fSToby Isaac { 408d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 409d1232d7fSToby Isaac 410d1232d7fSToby Isaac PetscFunctionBegin; 411d1232d7fSToby Isaac *name = vhdf5->filename; 412d1232d7fSToby Isaac PetscFunctionReturn(0); 413d1232d7fSToby Isaac } 414d1232d7fSToby Isaac 415b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 416b723ab35SVaclav Hapla { 4175dc64a97SVaclav Hapla /* 418b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 419b723ab35SVaclav Hapla PetscErrorCode ierr; 4205dc64a97SVaclav Hapla */ 421b723ab35SVaclav Hapla 422b723ab35SVaclav Hapla PetscFunctionBegin; 423b723ab35SVaclav Hapla PetscFunctionReturn(0); 424b723ab35SVaclav Hapla } 425b723ab35SVaclav Hapla 4268556b5ebSBarry Smith /*MC 4278556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4288556b5ebSBarry Smith 4298556b5ebSBarry Smith 4308556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4318556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4328556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4338556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4348556b5ebSBarry Smith 4351b266c99SBarry Smith Level: beginner 4368556b5ebSBarry Smith M*/ 437d1232d7fSToby Isaac 4388cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4395c6c1daeSBarry Smith { 4405c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4415c6c1daeSBarry Smith PetscErrorCode ierr; 4425c6c1daeSBarry Smith 4435c6c1daeSBarry Smith PetscFunctionBegin; 444b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4455c6c1daeSBarry Smith 4465c6c1daeSBarry Smith v->data = (void*) hdf5; 4475c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 44882ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 449b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4501b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 4515c6c1daeSBarry Smith v->ops->flush = 0; 4525c6c1daeSBarry Smith hdf5->btype = (PetscFileMode) -1; 4535c6c1daeSBarry Smith hdf5->filename = 0; 4545c6c1daeSBarry Smith hdf5->timestep = -1; 4550298fd71SBarry Smith hdf5->groups = NULL; 4565c6c1daeSBarry Smith 4579c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4589c5072fbSVaclav Hapla 4590b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 460d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4610b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4620b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 46382ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4649a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 465*ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 466*ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4675c6c1daeSBarry Smith PetscFunctionReturn(0); 4685c6c1daeSBarry Smith } 4695c6c1daeSBarry Smith 4705c6c1daeSBarry Smith /*@C 4715c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4725c6c1daeSBarry Smith 4735c6c1daeSBarry Smith Collective on MPI_Comm 4745c6c1daeSBarry Smith 4755c6c1daeSBarry Smith Input Parameters: 4765c6c1daeSBarry Smith + comm - MPI communicator 4775c6c1daeSBarry Smith . name - name of file 4785c6c1daeSBarry Smith - type - type of file 4795c6c1daeSBarry Smith $ FILE_MODE_WRITE - create new file for binary output 4805c6c1daeSBarry Smith $ FILE_MODE_READ - open existing file for binary input 4815c6c1daeSBarry Smith $ FILE_MODE_APPEND - open existing file for binary output 4825c6c1daeSBarry Smith 4835c6c1daeSBarry Smith Output Parameter: 4845c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 4855c6c1daeSBarry Smith 48682ea9b62SBarry Smith Options Database: 48782ea9b62SBarry 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 4889a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 48982ea9b62SBarry Smith 4905c6c1daeSBarry Smith Level: beginner 4915c6c1daeSBarry Smith 4925c6c1daeSBarry Smith Note: 4935c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 4945c6c1daeSBarry Smith 4955c6c1daeSBarry Smith Concepts: HDF5 files 4965c6c1daeSBarry Smith Concepts: PetscViewerHDF5^creating 4975c6c1daeSBarry Smith 4986a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 4999a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 500a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5015c6c1daeSBarry Smith @*/ 5025c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5035c6c1daeSBarry Smith { 5045c6c1daeSBarry Smith PetscErrorCode ierr; 5055c6c1daeSBarry Smith 5065c6c1daeSBarry Smith PetscFunctionBegin; 5075c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5085c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5095c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5105c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 511b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5125c6c1daeSBarry Smith PetscFunctionReturn(0); 5135c6c1daeSBarry Smith } 5145c6c1daeSBarry Smith 5155c6c1daeSBarry Smith /*@C 5165c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5175c6c1daeSBarry Smith 5185c6c1daeSBarry Smith Not collective 5195c6c1daeSBarry Smith 5205c6c1daeSBarry Smith Input Parameter: 5215c6c1daeSBarry Smith . viewer - the PetscViewer 5225c6c1daeSBarry Smith 5235c6c1daeSBarry Smith Output Parameter: 5245c6c1daeSBarry Smith . file_id - The file id 5255c6c1daeSBarry Smith 5265c6c1daeSBarry Smith Level: intermediate 5275c6c1daeSBarry Smith 5285c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5295c6c1daeSBarry Smith @*/ 5305c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5315c6c1daeSBarry Smith { 5325c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5335c6c1daeSBarry Smith 5345c6c1daeSBarry Smith PetscFunctionBegin; 5355c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5365c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5375c6c1daeSBarry Smith PetscFunctionReturn(0); 5385c6c1daeSBarry Smith } 5395c6c1daeSBarry Smith 5405c6c1daeSBarry Smith /*@C 5415c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5425c6c1daeSBarry Smith 5435c6c1daeSBarry Smith Not collective 5445c6c1daeSBarry Smith 5455c6c1daeSBarry Smith Input Parameters: 5465c6c1daeSBarry Smith + viewer - the PetscViewer 5475c6c1daeSBarry Smith - name - The group name 5485c6c1daeSBarry Smith 5495c6c1daeSBarry Smith Level: intermediate 5505c6c1daeSBarry Smith 551874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5525c6c1daeSBarry Smith @*/ 5535c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char *name) 5545c6c1daeSBarry Smith { 5555c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5567d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5575c6c1daeSBarry Smith PetscErrorCode ierr; 5585c6c1daeSBarry Smith 5595c6c1daeSBarry Smith PetscFunctionBegin; 5605c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5615c6c1daeSBarry Smith PetscValidCharPointer(name,2); 56295dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 5635c6c1daeSBarry Smith ierr = PetscStrallocpy(name, (char**) &groupNode->name);CHKERRQ(ierr); 564a297a907SKarl Rupp 5655c6c1daeSBarry Smith groupNode->next = hdf5->groups; 5665c6c1daeSBarry Smith hdf5->groups = groupNode; 5675c6c1daeSBarry Smith PetscFunctionReturn(0); 5685c6c1daeSBarry Smith } 5695c6c1daeSBarry Smith 5703ef9c667SSatish Balay /*@ 5715c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 5725c6c1daeSBarry Smith 5735c6c1daeSBarry Smith Not collective 5745c6c1daeSBarry Smith 5755c6c1daeSBarry Smith Input Parameter: 5765c6c1daeSBarry Smith . viewer - the PetscViewer 5775c6c1daeSBarry Smith 5785c6c1daeSBarry Smith Level: intermediate 5795c6c1daeSBarry Smith 580874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5815c6c1daeSBarry Smith @*/ 5825c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 5835c6c1daeSBarry Smith { 5845c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5857d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5865c6c1daeSBarry Smith PetscErrorCode ierr; 5875c6c1daeSBarry Smith 5885c6c1daeSBarry Smith PetscFunctionBegin; 5895c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 59082f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 5915c6c1daeSBarry Smith groupNode = hdf5->groups; 5925c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 5935c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 5945c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 5955c6c1daeSBarry Smith PetscFunctionReturn(0); 5965c6c1daeSBarry Smith } 5975c6c1daeSBarry Smith 5985c6c1daeSBarry Smith /*@C 599874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 600874270d9SVaclav Hapla If none has been assigned, returns NULL. 6015c6c1daeSBarry Smith 6025c6c1daeSBarry Smith Not collective 6035c6c1daeSBarry Smith 6045c6c1daeSBarry Smith Input Parameter: 6055c6c1daeSBarry Smith . viewer - the PetscViewer 6065c6c1daeSBarry Smith 6075c6c1daeSBarry Smith Output Parameter: 6085c6c1daeSBarry Smith . name - The group name 6095c6c1daeSBarry Smith 6105c6c1daeSBarry Smith Level: intermediate 6115c6c1daeSBarry Smith 612874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6135c6c1daeSBarry Smith @*/ 6145c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char **name) 6155c6c1daeSBarry Smith { 6165c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6175c6c1daeSBarry Smith 6185c6c1daeSBarry Smith PetscFunctionBegin; 6195c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 620c959eef4SJed Brown PetscValidPointer(name,2); 621a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6220298fd71SBarry Smith else *name = NULL; 6235c6c1daeSBarry Smith PetscFunctionReturn(0); 6245c6c1daeSBarry Smith } 6255c6c1daeSBarry Smith 6268c557b5aSVaclav Hapla /*@ 627874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 628874270d9SVaclav Hapla and return this group's ID and file ID. 629874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 630874270d9SVaclav Hapla 631874270d9SVaclav Hapla Not collective 632874270d9SVaclav Hapla 633874270d9SVaclav Hapla Input Parameter: 634874270d9SVaclav Hapla . viewer - the PetscViewer 635874270d9SVaclav Hapla 636874270d9SVaclav Hapla Output Parameter: 637874270d9SVaclav Hapla + fileId - The HDF5 file ID 638874270d9SVaclav Hapla - groupId - The HDF5 group ID 639874270d9SVaclav Hapla 640e74616beSVaclav Hapla Notes: 641e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 642e74616beSVaclav Hapla 643874270d9SVaclav Hapla Level: intermediate 644874270d9SVaclav Hapla 645874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 646874270d9SVaclav Hapla @*/ 64754dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 64854dbf706SVaclav Hapla { 64990e03537SVaclav Hapla hid_t file_id; 65090e03537SVaclav Hapla H5O_type_t type; 65154dbf706SVaclav Hapla const char *groupName = NULL; 652e74616beSVaclav Hapla PetscBool create; 65354dbf706SVaclav Hapla PetscErrorCode ierr; 65454dbf706SVaclav Hapla 65554dbf706SVaclav Hapla PetscFunctionBegin; 656e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 65754dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 65854dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 659e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 66090e03537SVaclav 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); 66190e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 66254dbf706SVaclav Hapla *fileId = file_id; 66354dbf706SVaclav Hapla PetscFunctionReturn(0); 66454dbf706SVaclav Hapla } 66554dbf706SVaclav Hapla 6663ef9c667SSatish Balay /*@ 6675c6c1daeSBarry Smith PetscViewerHDF5IncrementTimestep - Increments the current timestep for the HDF5 output. Fields are stacked in time. 6685c6c1daeSBarry Smith 6695c6c1daeSBarry Smith Not collective 6705c6c1daeSBarry Smith 6715c6c1daeSBarry Smith Input Parameter: 6725c6c1daeSBarry Smith . viewer - the PetscViewer 6735c6c1daeSBarry Smith 6745c6c1daeSBarry Smith Level: intermediate 6755c6c1daeSBarry Smith 6765c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 6775c6c1daeSBarry Smith @*/ 6785c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 6795c6c1daeSBarry Smith { 6805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6815c6c1daeSBarry Smith 6825c6c1daeSBarry Smith PetscFunctionBegin; 6835c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6845c6c1daeSBarry Smith ++hdf5->timestep; 6855c6c1daeSBarry Smith PetscFunctionReturn(0); 6865c6c1daeSBarry Smith } 6875c6c1daeSBarry Smith 6883ef9c667SSatish Balay /*@ 6895c6c1daeSBarry Smith PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. A timestep 6905c6c1daeSBarry Smith of -1 disables blocking with timesteps. 6915c6c1daeSBarry Smith 6925c6c1daeSBarry Smith Not collective 6935c6c1daeSBarry Smith 6945c6c1daeSBarry Smith Input Parameters: 6955c6c1daeSBarry Smith + viewer - the PetscViewer 6965c6c1daeSBarry Smith - timestep - The timestep number 6975c6c1daeSBarry Smith 6985c6c1daeSBarry Smith Level: intermediate 6995c6c1daeSBarry Smith 7005c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 7015c6c1daeSBarry Smith @*/ 7025c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 7035c6c1daeSBarry Smith { 7045c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7055c6c1daeSBarry Smith 7065c6c1daeSBarry Smith PetscFunctionBegin; 7075c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7085c6c1daeSBarry Smith hdf5->timestep = timestep; 7095c6c1daeSBarry Smith PetscFunctionReturn(0); 7105c6c1daeSBarry Smith } 7115c6c1daeSBarry Smith 7123ef9c667SSatish Balay /*@ 7135c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 7145c6c1daeSBarry Smith 7155c6c1daeSBarry Smith Not collective 7165c6c1daeSBarry Smith 7175c6c1daeSBarry Smith Input Parameter: 7185c6c1daeSBarry Smith . viewer - the PetscViewer 7195c6c1daeSBarry Smith 7205c6c1daeSBarry Smith Output Parameter: 7215c6c1daeSBarry Smith . timestep - The timestep number 7225c6c1daeSBarry Smith 7235c6c1daeSBarry Smith Level: intermediate 7245c6c1daeSBarry Smith 7255c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 7265c6c1daeSBarry Smith @*/ 7275c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 7285c6c1daeSBarry Smith { 7295c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7305c6c1daeSBarry Smith 7315c6c1daeSBarry Smith PetscFunctionBegin; 7325c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7335c6c1daeSBarry Smith PetscValidPointer(timestep,2); 7345c6c1daeSBarry Smith *timestep = hdf5->timestep; 7355c6c1daeSBarry Smith PetscFunctionReturn(0); 7365c6c1daeSBarry Smith } 7375c6c1daeSBarry Smith 73836bce6e8SMatthew G. Knepley /*@C 73936bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 74036bce6e8SMatthew G. Knepley 74136bce6e8SMatthew G. Knepley Not collective 74236bce6e8SMatthew G. Knepley 74336bce6e8SMatthew G. Knepley Input Parameter: 74436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 74536bce6e8SMatthew G. Knepley 74636bce6e8SMatthew G. Knepley Output Parameter: 74736bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 74836bce6e8SMatthew G. Knepley 74936bce6e8SMatthew G. Knepley Level: advanced 75036bce6e8SMatthew G. Knepley 75136bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 75236bce6e8SMatthew G. Knepley @*/ 75336bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 75436bce6e8SMatthew G. Knepley { 75536bce6e8SMatthew G. Knepley PetscFunctionBegin; 75636bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 75736bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 75836bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 75936bce6e8SMatthew G. Knepley #else 76036bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 76136bce6e8SMatthew G. Knepley #endif 76236bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 76336bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 76436bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 765c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 766de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 76736bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 76836bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 76936bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 7707e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 77136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 77236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 77336bce6e8SMatthew G. Knepley } 77436bce6e8SMatthew G. Knepley 77536bce6e8SMatthew G. Knepley /*@C 77636bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 77736bce6e8SMatthew G. Knepley 77836bce6e8SMatthew G. Knepley Not collective 77936bce6e8SMatthew G. Knepley 78036bce6e8SMatthew G. Knepley Input Parameter: 78136bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 78236bce6e8SMatthew G. Knepley 78336bce6e8SMatthew G. Knepley Output Parameter: 78436bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 78536bce6e8SMatthew G. Knepley 78636bce6e8SMatthew G. Knepley Level: advanced 78736bce6e8SMatthew G. Knepley 78836bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 78936bce6e8SMatthew G. Knepley @*/ 79036bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 79136bce6e8SMatthew G. Knepley { 79236bce6e8SMatthew G. Knepley PetscFunctionBegin; 79336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 79436bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 79536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 79636bce6e8SMatthew G. Knepley #else 79736bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 79836bce6e8SMatthew G. Knepley #endif 79936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 80036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 80136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 80236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 80336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 80436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 8057e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 80636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 80736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 80836bce6e8SMatthew G. Knepley } 80936bce6e8SMatthew G. Knepley 810df863907SAlex Fikl /*@C 811b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 81236bce6e8SMatthew G. Knepley 81336bce6e8SMatthew G. Knepley Input Parameters: 81436bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 81557d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 81636bce6e8SMatthew G. Knepley . name - The attribute name 81736bce6e8SMatthew G. Knepley . datatype - The attribute type 81836bce6e8SMatthew G. Knepley - value - The attribute value 81936bce6e8SMatthew G. Knepley 82036bce6e8SMatthew G. Knepley Level: advanced 82136bce6e8SMatthew G. Knepley 822577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 82336bce6e8SMatthew G. Knepley @*/ 82457d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, const void *value) 82536bce6e8SMatthew G. Knepley { 82657d22f7dSVaclav Hapla char *parent; 82760568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 828080f144cSVaclav Hapla PetscBool has; 82936bce6e8SMatthew G. Knepley PetscErrorCode ierr; 83036bce6e8SMatthew G. Knepley 83136bce6e8SMatthew G. Knepley PetscFunctionBegin; 8325cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 83357d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 834c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 835b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 83657d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 837bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 838b67695ecSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr); 83936bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 8407e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 8417e97c476SMatthew G. Knepley size_t len; 8427e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 843729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 8447e97c476SMatthew G. Knepley } 84536bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 846729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 84760568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 848080f144cSVaclav Hapla if (has) { 849080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 850080f144cSVaclav Hapla } else { 85160568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 852080f144cSVaclav Hapla } 853729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 854729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 855729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 85660568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 857729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 85857d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 85936bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 86036bce6e8SMatthew G. Knepley } 86136bce6e8SMatthew G. Knepley 862df863907SAlex Fikl /*@C 863577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 864577e0e04SVaclav Hapla 865577e0e04SVaclav Hapla Input Parameters: 866577e0e04SVaclav Hapla + viewer - The HDF5 viewer 867577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 868577e0e04SVaclav Hapla . name - The attribute name 869577e0e04SVaclav Hapla . datatype - The attribute type 870577e0e04SVaclav Hapla - value - The attribute value 871577e0e04SVaclav Hapla 872577e0e04SVaclav Hapla Notes: 873577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 874577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 875577e0e04SVaclav Hapla 876577e0e04SVaclav Hapla Level: advanced 877577e0e04SVaclav Hapla 878577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 879577e0e04SVaclav Hapla @*/ 880577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 881577e0e04SVaclav Hapla { 882577e0e04SVaclav Hapla PetscErrorCode ierr; 883577e0e04SVaclav Hapla 884577e0e04SVaclav Hapla PetscFunctionBegin; 885577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 886577e0e04SVaclav Hapla PetscValidHeader(obj,2); 887577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 888b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 889577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 890577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 891577e0e04SVaclav Hapla PetscFunctionReturn(0); 892577e0e04SVaclav Hapla } 893577e0e04SVaclav Hapla 894577e0e04SVaclav Hapla /*@C 895b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 896ad0c61feSMatthew G. Knepley 897ad0c61feSMatthew G. Knepley Input Parameters: 898ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 89957d22f7dSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 900ad0c61feSMatthew G. Knepley . name - The attribute name 901ad0c61feSMatthew G. Knepley - datatype - The attribute type 902ad0c61feSMatthew G. Knepley 903ad0c61feSMatthew G. Knepley Output Parameter: 904ad0c61feSMatthew G. Knepley . value - The attribute value 905ad0c61feSMatthew G. Knepley 906ad0c61feSMatthew G. Knepley Level: advanced 907ad0c61feSMatthew G. Knepley 908577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 909ad0c61feSMatthew G. Knepley @*/ 91057d22f7dSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscDataType datatype, void *value) 911ad0c61feSMatthew G. Knepley { 91257d22f7dSVaclav Hapla char *parent; 913f0b58479SMatthew G. Knepley hid_t h5, obj, attribute, atype, dtype; 914969748fdSVaclav Hapla PetscBool has; 915ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 916ad0c61feSMatthew G. Knepley 917ad0c61feSMatthew G. Knepley PetscFunctionBegin; 9185cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 91957d22f7dSVaclav Hapla if (dataset) PetscValidCharPointer(dataset, 2); 920c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 921b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 92257d22f7dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 923969748fdSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 924969748fdSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, &has);CHKERRQ(ierr);} 925969748fdSVaclav Hapla if (!has) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist", parent, name); 926ad0c61feSMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 927ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 92860568592SMatthew G. Knepley PetscStackCallHDF5Return(obj,H5Oopen,(h5, parent, H5P_DEFAULT)); 92960568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 930f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 931f0b58479SMatthew G. Knepley size_t len; 932e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 93315b861d2SVaclav Hapla PetscStackCall("H5Tget_size",len = H5Tget_size(atype)); 934f0b58479SMatthew G. Knepley PetscStackCallHDF5(H5Tclose,(atype)); 935f0b58479SMatthew G. Knepley ierr = PetscMalloc((len+1) * sizeof(char *), &value);CHKERRQ(ierr); 936f0b58479SMatthew G. Knepley } 93770efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 938729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 939e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 940e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 94157d22f7dSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 942ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 943ad0c61feSMatthew G. Knepley } 944ad0c61feSMatthew G. Knepley 945577e0e04SVaclav Hapla /*@C 946577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 947577e0e04SVaclav Hapla 948577e0e04SVaclav Hapla Input Parameters: 949577e0e04SVaclav Hapla + viewer - The HDF5 viewer 950577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 951577e0e04SVaclav Hapla . name - The attribute name 952577e0e04SVaclav Hapla - datatype - The attribute type 953577e0e04SVaclav Hapla 954577e0e04SVaclav Hapla Output Parameter: 955577e0e04SVaclav Hapla . value - The attribute value 956577e0e04SVaclav Hapla 957577e0e04SVaclav Hapla Notes: 958577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 959577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 960577e0e04SVaclav Hapla 961577e0e04SVaclav Hapla Level: advanced 962577e0e04SVaclav Hapla 963577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 964577e0e04SVaclav Hapla @*/ 965577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *value) 966577e0e04SVaclav Hapla { 967577e0e04SVaclav Hapla PetscErrorCode ierr; 968577e0e04SVaclav Hapla 969577e0e04SVaclav Hapla PetscFunctionBegin; 970577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 971577e0e04SVaclav Hapla PetscValidHeader(obj,2); 972577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 973b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 974577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 975577e0e04SVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 976577e0e04SVaclav Hapla PetscFunctionReturn(0); 977577e0e04SVaclav Hapla } 978577e0e04SVaclav Hapla 979a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 980a75c3fd4SVaclav Hapla { 981a75c3fd4SVaclav Hapla htri_t exists; 982a75c3fd4SVaclav Hapla hid_t group; 983a75c3fd4SVaclav Hapla 984a75c3fd4SVaclav Hapla PetscFunctionBegin; 985a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 986a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 987a75c3fd4SVaclav Hapla if (!exists && createGroup) { 988a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 989a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 990a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 991a75c3fd4SVaclav Hapla } 992a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 993a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 994a75c3fd4SVaclav Hapla } 995a75c3fd4SVaclav Hapla 996bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 9975cdeb108SMatthew G. Knepley { 99890e03537SVaclav Hapla const char rootGroupName[] = "/"; 9995cdeb108SMatthew G. Knepley hid_t h5; 1000e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 100115b861d2SVaclav Hapla PetscInt i; 100215b861d2SVaclav Hapla int n; 100385688ae2SVaclav Hapla char **hierarchy; 100485688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 10055cdeb108SMatthew G. Knepley PetscErrorCode ierr; 10065cdeb108SMatthew G. Knepley 10075cdeb108SMatthew G. Knepley PetscFunctionBegin; 10085cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 100990e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 101090e03537SVaclav Hapla else name = rootGroupName; 1011ccd66a83SVaclav Hapla if (has) { 101256cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 10135cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1014ccd66a83SVaclav Hapla } 1015ccd66a83SVaclav Hapla if (otype) { 1016ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 101756cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1018ccd66a83SVaclav Hapla } 1019c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 102085688ae2SVaclav Hapla 102185688ae2SVaclav Hapla /* 102285688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 102385688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 102485688ae2SVaclav Hapla 1) whether it's a valid link 102585688ae2SVaclav Hapla 2) whether this link resolves to an object 102685688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 102785688ae2SVaclav Hapla */ 102885688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 102985688ae2SVaclav Hapla if (!n) { 103085688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1031ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1032ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 103385688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 103485688ae2SVaclav Hapla PetscFunctionReturn(0); 103585688ae2SVaclav Hapla } 103685688ae2SVaclav Hapla for (i=0; i<n; i++) { 103785688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 103885688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1039a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1040a75c3fd4SVaclav Hapla if (!exists) break; 104185688ae2SVaclav Hapla } 104285688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 104385688ae2SVaclav Hapla 104485688ae2SVaclav Hapla /* If the object exists, get its type */ 1045ccd66a83SVaclav Hapla if (exists && otype) { 10465cdeb108SMatthew G. Knepley H5O_info_t info; 10475cdeb108SMatthew G. Knepley 104876276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 104904633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 105056cc0592SVaclav Hapla *otype = info.type; 10515cdeb108SMatthew G. Knepley } 1052ccd66a83SVaclav Hapla if (has) *has = exists; 10535cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 10545cdeb108SMatthew G. Knepley } 10555cdeb108SMatthew G. Knepley 1056ecce1506SVaclav Hapla /*@ 10570a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 10580a9f38efSVaclav Hapla 10590a9f38efSVaclav Hapla Input Parameters: 10600a9f38efSVaclav Hapla . viewer - The HDF5 viewer 10610a9f38efSVaclav Hapla 10620a9f38efSVaclav Hapla Output Parameter: 10630a9f38efSVaclav Hapla . has - Flag for group existence 10640a9f38efSVaclav Hapla 10650a9f38efSVaclav Hapla Notes: 10660a9f38efSVaclav Hapla If the path exists but is not a group, this returns PETSC_FALSE as well. 10670a9f38efSVaclav Hapla 10680a9f38efSVaclav Hapla Level: advanced 10690a9f38efSVaclav Hapla 10700a9f38efSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5OpenGroup() 10710a9f38efSVaclav Hapla @*/ 10720a9f38efSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, PetscBool *has) 10730a9f38efSVaclav Hapla { 10740a9f38efSVaclav Hapla H5O_type_t type; 10750a9f38efSVaclav Hapla const char *name; 10760a9f38efSVaclav Hapla PetscErrorCode ierr; 10770a9f38efSVaclav Hapla 10780a9f38efSVaclav Hapla PetscFunctionBegin; 10790a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10800a9f38efSVaclav Hapla PetscValidIntPointer(has,2); 10810a9f38efSVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &name);CHKERRQ(ierr); 10820a9f38efSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, name, PETSC_FALSE, has, &type);CHKERRQ(ierr); 10830a9f38efSVaclav Hapla *has = (type == H5O_TYPE_GROUP) ? PETSC_TRUE : PETSC_FALSE; 10840a9f38efSVaclav Hapla PetscFunctionReturn(0); 10850a9f38efSVaclav Hapla } 10860a9f38efSVaclav Hapla 10870a9f38efSVaclav Hapla /*@ 1088e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1089ecce1506SVaclav Hapla 1090ecce1506SVaclav Hapla Input Parameters: 1091ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1092ecce1506SVaclav Hapla - obj - The named object 1093ecce1506SVaclav Hapla 1094ecce1506SVaclav Hapla Output Parameter: 1095ecce1506SVaclav Hapla . has - Flag for dataset existence; PETSC_FALSE for unnamed object 1096ecce1506SVaclav Hapla 1097e3f143f7SVaclav Hapla Notes: 1098e3f143f7SVaclav Hapla If the path exists but is not a dataset, this returns PETSC_FALSE as well. 1099e3f143f7SVaclav Hapla 1100ecce1506SVaclav Hapla Level: advanced 1101ecce1506SVaclav Hapla 1102e3f143f7SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1103ecce1506SVaclav Hapla @*/ 1104ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1105ecce1506SVaclav Hapla { 110656cc0592SVaclav Hapla H5O_type_t type; 11076c132bc1SVaclav Hapla char *path; 1108ecce1506SVaclav Hapla PetscErrorCode ierr; 1109ecce1506SVaclav Hapla 1110ecce1506SVaclav Hapla PetscFunctionBegin; 1111c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1112c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 1113c57a1a9aSVaclav Hapla PetscValidIntPointer(has,3); 1114ecce1506SVaclav Hapla *has = PETSC_FALSE; 1115e3f143f7SVaclav Hapla if (!obj->name) PetscFunctionReturn(0); 11166c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, obj->name, &path);CHKERRQ(ierr); 11176c132bc1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, path, PETSC_FALSE, has, &type);CHKERRQ(ierr); 111856cc0592SVaclav Hapla *has = (type == H5O_TYPE_DATASET) ? PETSC_TRUE : PETSC_FALSE; 11196c132bc1SVaclav Hapla ierr = PetscFree(path);CHKERRQ(ierr); 1120ecce1506SVaclav Hapla PetscFunctionReturn(0); 1121ecce1506SVaclav Hapla } 1122ecce1506SVaclav Hapla 1123df863907SAlex Fikl /*@C 1124ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1125e7bdbf8eSMatthew G. Knepley 1126e7bdbf8eSMatthew G. Knepley Input Parameters: 1127e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 112810e69e8fSVaclav Hapla . dataset - The parent dataset name, relative to the current group. NULL means a group-wise attribute. 1129e7bdbf8eSMatthew G. Knepley - name - The attribute name 1130e7bdbf8eSMatthew G. Knepley 1131e7bdbf8eSMatthew G. Knepley Output Parameter: 1132e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1133e7bdbf8eSMatthew G. Knepley 1134e7bdbf8eSMatthew G. Knepley Level: advanced 1135e7bdbf8eSMatthew G. Knepley 1136577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1137e7bdbf8eSMatthew G. Knepley @*/ 113810e69e8fSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char dataset[], const char name[], PetscBool *has) 1139e7bdbf8eSMatthew G. Knepley { 114010e69e8fSVaclav Hapla char *parent; 1141e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1142e7bdbf8eSMatthew G. Knepley 1143e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 11445cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 114510e69e8fSVaclav Hapla if (dataset) PetscValidCharPointer(dataset,2); 1146c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 1147c57a1a9aSVaclav Hapla PetscValidIntPointer(has,4); 114810e69e8fSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, dataset, &parent);CHKERRQ(ierr); 1149bb286ee1SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parent, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 115010e69e8fSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parent, name, has);CHKERRQ(ierr);} 115110e69e8fSVaclav Hapla ierr = PetscFree(parent);CHKERRQ(ierr); 115206db490cSVaclav Hapla PetscFunctionReturn(0); 115306db490cSVaclav Hapla } 115406db490cSVaclav Hapla 1155577e0e04SVaclav Hapla /*@C 1156577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1157577e0e04SVaclav Hapla 1158577e0e04SVaclav Hapla Input Parameters: 1159577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1160577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1161577e0e04SVaclav Hapla - name - The attribute name 1162577e0e04SVaclav Hapla 1163577e0e04SVaclav Hapla Output Parameter: 1164577e0e04SVaclav Hapla . has - Flag for attribute existence 1165577e0e04SVaclav Hapla 1166577e0e04SVaclav Hapla Notes: 1167577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1168577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1169577e0e04SVaclav Hapla 1170577e0e04SVaclav Hapla Level: advanced 1171577e0e04SVaclav Hapla 1172577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1173577e0e04SVaclav Hapla @*/ 1174577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1175577e0e04SVaclav Hapla { 1176577e0e04SVaclav Hapla PetscErrorCode ierr; 1177577e0e04SVaclav Hapla 1178577e0e04SVaclav Hapla PetscFunctionBegin; 1179577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1180577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1181577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1182577e0e04SVaclav Hapla PetscValidIntPointer(has,4); 1183577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1184577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1185577e0e04SVaclav Hapla PetscFunctionReturn(0); 1186577e0e04SVaclav Hapla } 1187577e0e04SVaclav Hapla 118806db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 118906db490cSVaclav Hapla { 119006db490cSVaclav Hapla hid_t h5; 119106db490cSVaclav Hapla htri_t hhas; 119206db490cSVaclav Hapla PetscErrorCode ierr; 119306db490cSVaclav Hapla 119406db490cSVaclav Hapla PetscFunctionBegin; 119506db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11962f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 11975cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1198e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1199e7bdbf8eSMatthew G. Knepley } 1200e7bdbf8eSMatthew G. Knepley 1201a75e6a4aSMatthew G. Knepley /* 1202a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1203a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1204a75e6a4aSMatthew G. Knepley */ 1205d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1206a75e6a4aSMatthew G. Knepley 1207a75e6a4aSMatthew G. Knepley /*@C 1208a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1209a75e6a4aSMatthew G. Knepley 1210a75e6a4aSMatthew G. Knepley Collective on MPI_Comm 1211a75e6a4aSMatthew G. Knepley 1212a75e6a4aSMatthew G. Knepley Input Parameter: 1213a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1214a75e6a4aSMatthew G. Knepley 1215a75e6a4aSMatthew G. Knepley Level: intermediate 1216a75e6a4aSMatthew G. Knepley 1217a75e6a4aSMatthew G. Knepley Options Database Keys: 1218a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1219a75e6a4aSMatthew G. Knepley 1220a75e6a4aSMatthew G. Knepley Environmental variables: 1221a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1222a75e6a4aSMatthew G. Knepley 1223a75e6a4aSMatthew G. Knepley Notes: 1224a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1225a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1226a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1227a75e6a4aSMatthew G. Knepley 1228a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1229a75e6a4aSMatthew G. Knepley @*/ 1230a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1231a75e6a4aSMatthew G. Knepley { 1232a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1233a75e6a4aSMatthew G. Knepley PetscBool flg; 1234a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1235a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1236a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1237a75e6a4aSMatthew G. Knepley 1238a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1239a75e6a4aSMatthew 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);} 1240a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 124112801b39SBarry Smith ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,0); 1242a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1243a75e6a4aSMatthew G. Knepley } 124447435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1245a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1246a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1247a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 1248a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1249a75e6a4aSMatthew G. Knepley if (!flg) { 1250a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 1251a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1252a75e6a4aSMatthew G. Knepley } 1253a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 1254a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1255a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 1256a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 125747435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1258a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1259a75e6a4aSMatthew G. Knepley } 1260a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 1261a75e6a4aSMatthew G. Knepley if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(0);} 1262a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1263a75e6a4aSMatthew G. Knepley } 1264