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 84302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) 96c132bc1SVaclav Hapla { 104302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 116c132bc1SVaclav Hapla const char *group; 126c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 136c132bc1SVaclav Hapla PetscErrorCode ierr; 146c132bc1SVaclav Hapla 156c132bc1SVaclav Hapla PetscFunctionBegin; 166c132bc1SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 174302210dSVaclav Hapla ierr = PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative);CHKERRQ(ierr); 184302210dSVaclav Hapla if (relative) { 194302210dSVaclav Hapla ierr = PetscStrcpy(buf, group);CHKERRQ(ierr); 206c132bc1SVaclav Hapla ierr = PetscStrcat(buf, "/");CHKERRQ(ierr); 214302210dSVaclav Hapla ierr = PetscStrcat(buf, path);CHKERRQ(ierr); 224302210dSVaclav Hapla ierr = PetscStrallocpy(buf, abspath);CHKERRQ(ierr); 234302210dSVaclav Hapla } else { 244302210dSVaclav Hapla ierr = PetscStrallocpy(path, abspath);CHKERRQ(ierr); 254302210dSVaclav Hapla } 266c132bc1SVaclav Hapla PetscFunctionReturn(0); 276c132bc1SVaclav Hapla } 286c132bc1SVaclav Hapla 29577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 30577e0e04SVaclav Hapla { 31577e0e04SVaclav Hapla PetscBool has; 32577e0e04SVaclav Hapla PetscErrorCode ierr; 33577e0e04SVaclav Hapla 34577e0e04SVaclav Hapla PetscFunctionBegin; 35577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasObject(viewer, obj, &has);CHKERRQ(ierr); 36577e0e04SVaclav Hapla if (!has) { 3789e0ef10SVaclav Hapla const char *group; 38577e0e04SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &group);CHKERRQ(ierr); 3989e0ef10SVaclav Hapla SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/"); 40577e0e04SVaclav Hapla } 41577e0e04SVaclav Hapla PetscFunctionReturn(0); 42577e0e04SVaclav Hapla } 43577e0e04SVaclav Hapla 444416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4582ea9b62SBarry Smith { 4682ea9b62SBarry Smith PetscErrorCode ierr; 47ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4882ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4982ea9b62SBarry Smith 5082ea9b62SBarry Smith PetscFunctionBegin; 5182ea9b62SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"HDF5 PetscViewer Options");CHKERRQ(ierr); 5282ea9b62SBarry Smith ierr = PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL);CHKERRQ(ierr); 539a0502c6SHåkon Strandenes ierr = PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL);CHKERRQ(ierr); 54ee8b9732SVaclav Hapla ierr = PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set);CHKERRQ(ierr); 55ee8b9732SVaclav Hapla if (set) {ierr = PetscViewerHDF5SetCollective(v,flg);CHKERRQ(ierr);} 5682ea9b62SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 5782ea9b62SBarry Smith PetscFunctionReturn(0); 5882ea9b62SBarry Smith } 5982ea9b62SBarry Smith 601b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 611b793a25SVaclav Hapla { 621b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6303fa8834SVaclav Hapla PetscBool flg; 641b793a25SVaclav Hapla PetscErrorCode ierr; 651b793a25SVaclav Hapla 661b793a25SVaclav Hapla PetscFunctionBegin; 671b793a25SVaclav Hapla if (hdf5->filename) { 681b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename);CHKERRQ(ierr); 691b793a25SVaclav Hapla } 701b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2]);CHKERRQ(ierr); 711b793a25SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput]);CHKERRQ(ierr); 7203fa8834SVaclav Hapla ierr = PetscViewerHDF5GetCollective(v,&flg);CHKERRQ(ierr); 7303fa8834SVaclav Hapla ierr = PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent");CHKERRQ(ierr); 741b793a25SVaclav Hapla PetscFunctionReturn(0); 751b793a25SVaclav Hapla } 761b793a25SVaclav Hapla 775c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 785c6c1daeSBarry Smith { 795c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 805c6c1daeSBarry Smith PetscErrorCode ierr; 815c6c1daeSBarry Smith 825c6c1daeSBarry Smith PetscFunctionBegin; 835c6c1daeSBarry Smith ierr = PetscFree(hdf5->filename);CHKERRQ(ierr); 84729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 855c6c1daeSBarry Smith PetscFunctionReturn(0); 865c6c1daeSBarry Smith } 875c6c1daeSBarry Smith 889b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 895c6c1daeSBarry Smith { 905c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 915c6c1daeSBarry Smith PetscErrorCode ierr; 925c6c1daeSBarry Smith 935c6c1daeSBarry Smith PetscFunctionBegin; 949c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 955c6c1daeSBarry Smith ierr = PetscViewerFileClose_HDF5(viewer);CHKERRQ(ierr); 965c6c1daeSBarry Smith while (hdf5->groups) { 977d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 985c6c1daeSBarry Smith 995c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups->name);CHKERRQ(ierr); 1005c6c1daeSBarry Smith ierr = PetscFree(hdf5->groups);CHKERRQ(ierr); 1015c6c1daeSBarry Smith hdf5->groups = tmp; 1025c6c1daeSBarry Smith } 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); 1095c6c1daeSBarry Smith PetscFunctionReturn(0); 1105c6c1daeSBarry Smith } 1115c6c1daeSBarry Smith 1129b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1135c6c1daeSBarry Smith { 1145c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1155c6c1daeSBarry Smith 1165c6c1daeSBarry Smith PetscFunctionBegin; 1175c6c1daeSBarry Smith hdf5->btype = type; 1185c6c1daeSBarry Smith PetscFunctionReturn(0); 1195c6c1daeSBarry Smith } 1205c6c1daeSBarry Smith 1210b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1220b62783dSVaclav Hapla { 1230b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1240b62783dSVaclav Hapla 1250b62783dSVaclav Hapla PetscFunctionBegin; 1260b62783dSVaclav Hapla *type = hdf5->btype; 1270b62783dSVaclav Hapla PetscFunctionReturn(0); 1280b62783dSVaclav Hapla } 1290b62783dSVaclav Hapla 1309b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13182ea9b62SBarry Smith { 13282ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 13382ea9b62SBarry Smith 13482ea9b62SBarry Smith PetscFunctionBegin; 13582ea9b62SBarry Smith hdf5->basedimension2 = flg; 13682ea9b62SBarry Smith PetscFunctionReturn(0); 13782ea9b62SBarry Smith } 13882ea9b62SBarry Smith 139df863907SAlex Fikl /*@ 14082ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14182ea9b62SBarry Smith dimension of 2. 14282ea9b62SBarry Smith 14382ea9b62SBarry Smith Logically Collective on PetscViewer 14482ea9b62SBarry Smith 14582ea9b62SBarry Smith Input Parameters: 14682ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 14782ea9b62SBarry Smith - flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 14882ea9b62SBarry Smith 14982ea9b62SBarry Smith Options Database: 15082ea9b62SBarry Smith . -viewer_hdf5_base_dimension2 - turns on (true) or off (false) using a dimension of 2 in the HDF5 file even if the bs/dof of the vector is 1 15182ea9b62SBarry Smith 15282ea9b62SBarry Smith 15395452b02SPatrick Sanan Notes: 15495452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 15582ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 15682ea9b62SBarry Smith 15782ea9b62SBarry Smith Level: intermediate 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 16082ea9b62SBarry Smith 16182ea9b62SBarry Smith @*/ 16282ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16382ea9b62SBarry Smith { 16482ea9b62SBarry Smith PetscErrorCode ierr; 16582ea9b62SBarry Smith 16682ea9b62SBarry Smith PetscFunctionBegin; 16782ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 16882ea9b62SBarry Smith ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 16982ea9b62SBarry Smith PetscFunctionReturn(0); 17082ea9b62SBarry Smith } 17182ea9b62SBarry Smith 172df863907SAlex Fikl /*@ 17382ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17482ea9b62SBarry Smith dimension of 2. 17582ea9b62SBarry Smith 17682ea9b62SBarry Smith Logically Collective on PetscViewer 17782ea9b62SBarry Smith 17882ea9b62SBarry Smith Input Parameter: 17982ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18082ea9b62SBarry Smith 18182ea9b62SBarry Smith Output Parameter: 18282ea9b62SBarry Smith . flg - if PETSC_TRUE the vector will always have at least a dimension of 2 even if that first dimension is of size 1 18382ea9b62SBarry Smith 18495452b02SPatrick Sanan Notes: 18595452b02SPatrick Sanan Setting this option allegedly makes code that reads the HDF5 in easier since they do not have a "special case" of a bs/dof 18682ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 18782ea9b62SBarry Smith 18882ea9b62SBarry Smith Level: intermediate 18982ea9b62SBarry Smith 19082ea9b62SBarry Smith .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen() 19182ea9b62SBarry Smith 19282ea9b62SBarry Smith @*/ 19382ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19482ea9b62SBarry Smith { 19582ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 19682ea9b62SBarry Smith 19782ea9b62SBarry Smith PetscFunctionBegin; 19882ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 19982ea9b62SBarry Smith *flg = hdf5->basedimension2; 20082ea9b62SBarry Smith PetscFunctionReturn(0); 20182ea9b62SBarry Smith } 20282ea9b62SBarry Smith 2039b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2049a0502c6SHåkon Strandenes { 2059a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2069a0502c6SHåkon Strandenes 2079a0502c6SHåkon Strandenes PetscFunctionBegin; 2089a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2099a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2109a0502c6SHåkon Strandenes } 2119a0502c6SHåkon Strandenes 212df863907SAlex Fikl /*@ 2139a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2149a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2159a0502c6SHåkon Strandenes 2169a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2179a0502c6SHåkon Strandenes 2189a0502c6SHåkon Strandenes Input Parameters: 2199a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2209a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2219a0502c6SHåkon Strandenes 2229a0502c6SHåkon Strandenes Options Database: 2239a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2249a0502c6SHåkon Strandenes 2259a0502c6SHåkon Strandenes 22695452b02SPatrick Sanan Notes: 22795452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2289a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2299a0502c6SHåkon Strandenes 2309a0502c6SHåkon Strandenes Level: intermediate 2319a0502c6SHåkon Strandenes 2329a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2339a0502c6SHåkon Strandenes PetscReal 2349a0502c6SHåkon Strandenes 2359a0502c6SHåkon Strandenes @*/ 2369a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2379a0502c6SHåkon Strandenes { 2389a0502c6SHåkon Strandenes PetscErrorCode ierr; 2399a0502c6SHåkon Strandenes 2409a0502c6SHåkon Strandenes PetscFunctionBegin; 2419a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2429a0502c6SHåkon Strandenes ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 2439a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2449a0502c6SHåkon Strandenes } 2459a0502c6SHåkon Strandenes 246df863907SAlex Fikl /*@ 2479a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2489a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2499a0502c6SHåkon Strandenes 2509a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2519a0502c6SHåkon Strandenes 2529a0502c6SHåkon Strandenes Input Parameter: 2539a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2549a0502c6SHåkon Strandenes 2559a0502c6SHåkon Strandenes Output Parameter: 2569a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2579a0502c6SHåkon Strandenes 25895452b02SPatrick Sanan Notes: 25995452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2609a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2619a0502c6SHåkon Strandenes 2629a0502c6SHåkon Strandenes Level: intermediate 2639a0502c6SHåkon Strandenes 2649a0502c6SHåkon Strandenes .seealso: PetscViewerFileSetMode(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerBinaryOpen(), 2659a0502c6SHåkon Strandenes PetscReal 2669a0502c6SHåkon Strandenes 2679a0502c6SHåkon Strandenes @*/ 2689a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2699a0502c6SHåkon Strandenes { 2709a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2719a0502c6SHåkon Strandenes 2729a0502c6SHåkon Strandenes PetscFunctionBegin; 2739a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2749a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2759a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2769a0502c6SHåkon Strandenes } 2779a0502c6SHåkon Strandenes 278ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 279ee8b9732SVaclav Hapla { 280ee8b9732SVaclav Hapla PetscFunctionBegin; 281ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 282ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 283c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 284ee8b9732SVaclav Hapla { 285ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 286ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 287ee8b9732SVaclav Hapla } 288ee8b9732SVaclav Hapla #else 289ee8b9732SVaclav Hapla if (flg) { 290ee8b9732SVaclav Hapla PetscErrorCode ierr; 291c796909eSBarry Smith ierr = PetscPrintf(PetscObjectComm((PetscObject)viewer), "Warning: PetscViewerHDF5SetCollective(viewer,PETSC_TRUE) is ignored for HDF5 versions prior to 1.10.3 or if built without MPI support\n");CHKERRQ(ierr); 292ee8b9732SVaclav Hapla } 293ee8b9732SVaclav Hapla #endif 294ee8b9732SVaclav Hapla PetscFunctionReturn(0); 295ee8b9732SVaclav Hapla } 296ee8b9732SVaclav Hapla 297ee8b9732SVaclav Hapla /*@ 298ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 299ee8b9732SVaclav Hapla 300ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 301ee8b9732SVaclav Hapla 302ee8b9732SVaclav Hapla Input Parameters: 303ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 304ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 305ee8b9732SVaclav Hapla 306ee8b9732SVaclav Hapla Options Database: 307ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 308ee8b9732SVaclav Hapla 309ee8b9732SVaclav Hapla Notes: 310ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 31153effed7SBarry Smith However, this works correctly only since HDF5 1.10.3 and if HDF5 is installed for MPI; hence, we ignore this setting for older versions. 312ee8b9732SVaclav Hapla 313ee8b9732SVaclav Hapla Developer notes: 314ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 315ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 316ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 317ee8b9732SVaclav Hapla 318ee8b9732SVaclav Hapla Level: intermediate 319ee8b9732SVaclav Hapla 320ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5GetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 321ee8b9732SVaclav Hapla 322ee8b9732SVaclav Hapla @*/ 323ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 324ee8b9732SVaclav Hapla { 325ee8b9732SVaclav Hapla PetscErrorCode ierr; 326ee8b9732SVaclav Hapla 327ee8b9732SVaclav Hapla PetscFunctionBegin; 328ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 329ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 330ee8b9732SVaclav Hapla ierr = PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg));CHKERRQ(ierr); 331ee8b9732SVaclav Hapla PetscFunctionReturn(0); 332ee8b9732SVaclav Hapla } 333ee8b9732SVaclav Hapla 334ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 335ee8b9732SVaclav Hapla { 336c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 337ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 338ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 339c796909eSBarry Smith #endif 340ee8b9732SVaclav Hapla 341ee8b9732SVaclav Hapla PetscFunctionBegin; 342c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 343c796909eSBarry Smith *flg = PETSC_FALSE; 344c796909eSBarry Smith #else 345ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 346ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 347c796909eSBarry Smith #endif 348ee8b9732SVaclav Hapla PetscFunctionReturn(0); 349ee8b9732SVaclav Hapla } 350ee8b9732SVaclav Hapla 351ee8b9732SVaclav Hapla /*@ 352ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 353ee8b9732SVaclav Hapla 354ee8b9732SVaclav Hapla Not Collective 355ee8b9732SVaclav Hapla 356ee8b9732SVaclav Hapla Input Parameters: 357ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 358ee8b9732SVaclav Hapla 359ee8b9732SVaclav Hapla Output Parameters: 360ee8b9732SVaclav Hapla . flg - the flag 361ee8b9732SVaclav Hapla 362ee8b9732SVaclav Hapla Level: intermediate 363ee8b9732SVaclav Hapla 364ee8b9732SVaclav Hapla Notes: 365c796909eSBarry Smith This setting works correctly only since HDF5 1.10.3 and if HDF5 was installed for MPI. For older versions, PETSC_FALSE will be always returned. 366ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 367ee8b9732SVaclav Hapla 368ee8b9732SVaclav Hapla .seealso: PetscViewerHDF5SetCollective(), PetscViewerCreate(), PetscViewerSetType(), PetscViewerHDF5Open() 369ee8b9732SVaclav Hapla 370ee8b9732SVaclav Hapla @*/ 371ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 372ee8b9732SVaclav Hapla { 373ee8b9732SVaclav Hapla PetscErrorCode ierr; 374ee8b9732SVaclav Hapla 375ee8b9732SVaclav Hapla PetscFunctionBegin; 376ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 377534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 378ee8b9732SVaclav Hapla 379ee8b9732SVaclav Hapla ierr = PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg));CHKERRQ(ierr); 380ee8b9732SVaclav Hapla PetscFunctionReturn(0); 381ee8b9732SVaclav Hapla } 382ee8b9732SVaclav Hapla 3839b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3845c6c1daeSBarry Smith { 3855c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3865c6c1daeSBarry Smith hid_t plist_id; 3875c6c1daeSBarry Smith PetscErrorCode ierr; 3885c6c1daeSBarry Smith 3895c6c1daeSBarry Smith PetscFunctionBegin; 390f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 391f683cc58SToby Isaac if (hdf5->filename) {ierr = PetscFree(hdf5->filename);CHKERRQ(ierr);} 3925c6c1daeSBarry Smith ierr = PetscStrallocpy(name, &hdf5->filename);CHKERRQ(ierr); 3935c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 394729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 395c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 396d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 397c796909eSBarry Smith #endif 3985c6c1daeSBarry Smith /* Create or open the file collectively */ 3995c6c1daeSBarry Smith switch (hdf5->btype) { 4005c6c1daeSBarry Smith case FILE_MODE_READ: 401729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4025c6c1daeSBarry Smith break; 4035c6c1daeSBarry Smith case FILE_MODE_APPEND: 4047e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4057e4fd573SVaclav Hapla { 4067e4fd573SVaclav Hapla PetscBool flg; 4077e4fd573SVaclav Hapla ierr = PetscTestFile(hdf5->filename, 'r', &flg);CHKERRQ(ierr); 4087e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4097e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4105c6c1daeSBarry Smith break; 4117e4fd573SVaclav Hapla } 4125c6c1daeSBarry Smith case FILE_MODE_WRITE: 413729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4145c6c1daeSBarry Smith break; 4157e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4167e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4175c6c1daeSBarry Smith default: 4187e4fd573SVaclav Hapla SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4195c6c1daeSBarry Smith } 4205c6c1daeSBarry Smith if (hdf5->file_id < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 421729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4225c6c1daeSBarry Smith PetscFunctionReturn(0); 4235c6c1daeSBarry Smith } 4245c6c1daeSBarry Smith 425d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 426d1232d7fSToby Isaac { 427d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 428d1232d7fSToby Isaac 429d1232d7fSToby Isaac PetscFunctionBegin; 430d1232d7fSToby Isaac *name = vhdf5->filename; 431d1232d7fSToby Isaac PetscFunctionReturn(0); 432d1232d7fSToby Isaac } 433d1232d7fSToby Isaac 434b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 435b723ab35SVaclav Hapla { 4365dc64a97SVaclav Hapla /* 437b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 438b723ab35SVaclav Hapla PetscErrorCode ierr; 4395dc64a97SVaclav Hapla */ 440b723ab35SVaclav Hapla 441b723ab35SVaclav Hapla PetscFunctionBegin; 442b723ab35SVaclav Hapla PetscFunctionReturn(0); 443b723ab35SVaclav Hapla } 444b723ab35SVaclav Hapla 4458556b5ebSBarry Smith /*MC 4468556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 4478556b5ebSBarry Smith 4488556b5ebSBarry Smith 4498556b5ebSBarry Smith .seealso: PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET, 4508556b5ebSBarry Smith PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING, 4518556b5ebSBarry Smith PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB, 4528556b5ebSBarry Smith PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType() 4538556b5ebSBarry Smith 4541b266c99SBarry Smith Level: beginner 4558556b5ebSBarry Smith M*/ 456d1232d7fSToby Isaac 4578cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 4585c6c1daeSBarry Smith { 4595c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 4605c6c1daeSBarry Smith PetscErrorCode ierr; 4615c6c1daeSBarry Smith 4625c6c1daeSBarry Smith PetscFunctionBegin; 463b00a9115SJed Brown ierr = PetscNewLog(v,&hdf5);CHKERRQ(ierr); 4645c6c1daeSBarry Smith 4655c6c1daeSBarry Smith v->data = (void*) hdf5; 4665c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 46782ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 468b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 4691b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 470908793a3SLisandro Dalcin v->ops->flush = NULL; 4717e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 472908793a3SLisandro Dalcin hdf5->filename = NULL; 4735c6c1daeSBarry Smith hdf5->timestep = -1; 4740298fd71SBarry Smith hdf5->groups = NULL; 4755c6c1daeSBarry Smith 4769c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 4779c5072fbSVaclav Hapla 4780b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5);CHKERRQ(ierr); 479d1232d7fSToby Isaac ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5);CHKERRQ(ierr); 4800b062f91SJed Brown ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5);CHKERRQ(ierr); 4810b62783dSVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5);CHKERRQ(ierr); 48282ea9b62SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5);CHKERRQ(ierr); 4839a0502c6SHåkon Strandenes ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5);CHKERRQ(ierr); 484ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5);CHKERRQ(ierr); 485ee8b9732SVaclav Hapla ierr = PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5);CHKERRQ(ierr); 4865c6c1daeSBarry Smith PetscFunctionReturn(0); 4875c6c1daeSBarry Smith } 4885c6c1daeSBarry Smith 4895c6c1daeSBarry Smith /*@C 4905c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 4915c6c1daeSBarry Smith 492d083f849SBarry Smith Collective 4935c6c1daeSBarry Smith 4945c6c1daeSBarry Smith Input Parameters: 4955c6c1daeSBarry Smith + comm - MPI communicator 4965c6c1daeSBarry Smith . name - name of file 4975c6c1daeSBarry Smith - type - type of file 4985c6c1daeSBarry Smith 4995c6c1daeSBarry Smith Output Parameter: 5005c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5015c6c1daeSBarry Smith 50282ea9b62SBarry Smith Options Database: 503a2b725a8SWilliam Gropp + -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 504a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 50582ea9b62SBarry Smith 5065c6c1daeSBarry Smith Level: beginner 5075c6c1daeSBarry Smith 5087e4fd573SVaclav Hapla Notes: 5097e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5107e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5117e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5127e4fd573SVaclav Hapla . FILE_MODE_APPEND - if file exists, keep existing contents [H5Fopen() with H5F_ACC_RDWR], else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_EXCL] 5137e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5147e4fd573SVaclav Hapla 5157e4fd573SVaclav Hapla In case of FILE_MODE_APPEND / FILE_MODE_UPDATE, any stored object (dataset, attribute) can be selectively ovewritten if the same fully qualified name (/group/path/to/object) is specified. 5167e4fd573SVaclav Hapla 5175c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 5185c6c1daeSBarry Smith 5195c6c1daeSBarry Smith 5206a9046bcSBarry Smith .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(), PetscViewerHDF5SetBaseDimension2(), 5219a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput(), PetscViewerHDF5GetBaseDimension2(), VecView(), MatView(), VecLoad(), 522a56f64adSBarry Smith MatLoad(), PetscFileMode, PetscViewer, PetscViewerSetType(), PetscViewerFileSetMode(), PetscViewerFileSetName() 5235c6c1daeSBarry Smith @*/ 5245c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 5255c6c1daeSBarry Smith { 5265c6c1daeSBarry Smith PetscErrorCode ierr; 5275c6c1daeSBarry Smith 5285c6c1daeSBarry Smith PetscFunctionBegin; 5295c6c1daeSBarry Smith ierr = PetscViewerCreate(comm, hdf5v);CHKERRQ(ierr); 5305c6c1daeSBarry Smith ierr = PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5);CHKERRQ(ierr); 5315c6c1daeSBarry Smith ierr = PetscViewerFileSetMode(*hdf5v, type);CHKERRQ(ierr); 5325c6c1daeSBarry Smith ierr = PetscViewerFileSetName(*hdf5v, name);CHKERRQ(ierr); 533b663d2d2SVaclav Hapla ierr = PetscViewerSetFromOptions(*hdf5v);CHKERRQ(ierr); 5345c6c1daeSBarry Smith PetscFunctionReturn(0); 5355c6c1daeSBarry Smith } 5365c6c1daeSBarry Smith 5375c6c1daeSBarry Smith /*@C 5385c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 5395c6c1daeSBarry Smith 5405c6c1daeSBarry Smith Not collective 5415c6c1daeSBarry Smith 5425c6c1daeSBarry Smith Input Parameter: 5435c6c1daeSBarry Smith . viewer - the PetscViewer 5445c6c1daeSBarry Smith 5455c6c1daeSBarry Smith Output Parameter: 5465c6c1daeSBarry Smith . file_id - The file id 5475c6c1daeSBarry Smith 5485c6c1daeSBarry Smith Level: intermediate 5495c6c1daeSBarry Smith 5505c6c1daeSBarry Smith .seealso: PetscViewerHDF5Open() 5515c6c1daeSBarry Smith @*/ 5525c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 5535c6c1daeSBarry Smith { 5545c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5555c6c1daeSBarry Smith 5565c6c1daeSBarry Smith PetscFunctionBegin; 5575c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 5585c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 5595c6c1daeSBarry Smith PetscFunctionReturn(0); 5605c6c1daeSBarry Smith } 5615c6c1daeSBarry Smith 5625c6c1daeSBarry Smith /*@C 5635c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 5645c6c1daeSBarry Smith 5655c6c1daeSBarry Smith Not collective 5665c6c1daeSBarry Smith 5675c6c1daeSBarry Smith Input Parameters: 5685c6c1daeSBarry Smith + viewer - the PetscViewer 5695c6c1daeSBarry Smith - name - The group name 5705c6c1daeSBarry Smith 5715c6c1daeSBarry Smith Level: intermediate 5725c6c1daeSBarry Smith 5732e361344SVaclav Hapla Notes: 5742e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 5752e361344SVaclav Hapla + If name begins with '/', it is interpreted as an absolute path fully replacing current group, otherwise it is taken as relative to the current group. 5762e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 5772e361344SVaclav Hapla - "." means the current group is pushed again. 5782e361344SVaclav Hapla 5792e361344SVaclav Hapla Example: 5802e361344SVaclav Hapla Suppose the current group is "/a". 5812e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 5822e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 5832e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 5842e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 5852e361344SVaclav Hapla 5862e361344SVaclav Hapla Developer Notes: 5872e361344SVaclav Hapla The root group "/" is internally stored as NULL. 588820fc2d1SVaclav Hapla 589874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 5905c6c1daeSBarry Smith @*/ 591be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 5925c6c1daeSBarry Smith { 5935c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 5947d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 5952e361344SVaclav Hapla size_t i,len; 5962e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 5972e361344SVaclav Hapla const char *gname; 5985c6c1daeSBarry Smith PetscErrorCode ierr; 5995c6c1daeSBarry Smith 6005c6c1daeSBarry Smith PetscFunctionBegin; 6015c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 602820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 603820fc2d1SVaclav Hapla ierr = PetscStrlen(name, &len);CHKERRQ(ierr); 6042e361344SVaclav Hapla gname = NULL; 6052e361344SVaclav Hapla if (len) { 6062e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6072e361344SVaclav Hapla /* use current name */ 6082e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6092e361344SVaclav Hapla } else if (name[0] == '/') { 6102e361344SVaclav Hapla /* absolute */ 6112e361344SVaclav Hapla for (i=1; i<len; i++) { 6122e361344SVaclav Hapla if (name[i] != '/') { 6132e361344SVaclav Hapla gname = name; 6142e361344SVaclav Hapla break; 6152e361344SVaclav Hapla } 6162e361344SVaclav Hapla } 6172e361344SVaclav Hapla } else { 6182e361344SVaclav Hapla /* relative */ 6192e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6202e361344SVaclav Hapla ierr = PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name);CHKERRQ(ierr); 6212e361344SVaclav Hapla gname = buf; 6222e361344SVaclav Hapla } 6232e361344SVaclav Hapla } 62495dccacaSBarry Smith ierr = PetscNew(&groupNode);CHKERRQ(ierr); 6252e361344SVaclav Hapla ierr = PetscStrallocpy(gname, (char**) &groupNode->name);CHKERRQ(ierr); 6265c6c1daeSBarry Smith groupNode->next = hdf5->groups; 6275c6c1daeSBarry Smith hdf5->groups = groupNode; 6285c6c1daeSBarry Smith PetscFunctionReturn(0); 6295c6c1daeSBarry Smith } 6305c6c1daeSBarry Smith 6313ef9c667SSatish Balay /*@ 6325c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 6335c6c1daeSBarry Smith 6345c6c1daeSBarry Smith Not collective 6355c6c1daeSBarry Smith 6365c6c1daeSBarry Smith Input Parameter: 6375c6c1daeSBarry Smith . viewer - the PetscViewer 6385c6c1daeSBarry Smith 6395c6c1daeSBarry Smith Level: intermediate 6405c6c1daeSBarry Smith 641874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5GetGroup(),PetscViewerHDF5OpenGroup() 6425c6c1daeSBarry Smith @*/ 6435c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 6445c6c1daeSBarry Smith { 6455c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6467d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6475c6c1daeSBarry Smith PetscErrorCode ierr; 6485c6c1daeSBarry Smith 6495c6c1daeSBarry Smith PetscFunctionBegin; 6505c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 65182f516ccSBarry Smith if (!hdf5->groups) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 6525c6c1daeSBarry Smith groupNode = hdf5->groups; 6535c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 6545c6c1daeSBarry Smith ierr = PetscFree(groupNode->name);CHKERRQ(ierr); 6555c6c1daeSBarry Smith ierr = PetscFree(groupNode);CHKERRQ(ierr); 6565c6c1daeSBarry Smith PetscFunctionReturn(0); 6575c6c1daeSBarry Smith } 6585c6c1daeSBarry Smith 6595c6c1daeSBarry Smith /*@C 660874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 661874270d9SVaclav Hapla If none has been assigned, returns NULL. 6625c6c1daeSBarry Smith 6635c6c1daeSBarry Smith Not collective 6645c6c1daeSBarry Smith 6655c6c1daeSBarry Smith Input Parameter: 6665c6c1daeSBarry Smith . viewer - the PetscViewer 6675c6c1daeSBarry Smith 6685c6c1daeSBarry Smith Output Parameter: 6695c6c1daeSBarry Smith . name - The group name 6705c6c1daeSBarry Smith 6715c6c1daeSBarry Smith Level: intermediate 6725c6c1daeSBarry Smith 673874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5OpenGroup() 6745c6c1daeSBarry Smith @*/ 675be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 6765c6c1daeSBarry Smith { 6775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 6785c6c1daeSBarry Smith 6795c6c1daeSBarry Smith PetscFunctionBegin; 6805c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 681c959eef4SJed Brown PetscValidPointer(name,2); 682a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 6830298fd71SBarry Smith else *name = NULL; 6845c6c1daeSBarry Smith PetscFunctionReturn(0); 6855c6c1daeSBarry Smith } 6865c6c1daeSBarry Smith 6878c557b5aSVaclav Hapla /*@ 688874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 689874270d9SVaclav Hapla and return this group's ID and file ID. 690874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 691874270d9SVaclav Hapla 692874270d9SVaclav Hapla Not collective 693874270d9SVaclav Hapla 694874270d9SVaclav Hapla Input Parameter: 695874270d9SVaclav Hapla . viewer - the PetscViewer 696874270d9SVaclav Hapla 697874270d9SVaclav Hapla Output Parameter: 698874270d9SVaclav Hapla + fileId - The HDF5 file ID 699874270d9SVaclav Hapla - groupId - The HDF5 group ID 700874270d9SVaclav Hapla 701e74616beSVaclav Hapla Notes: 702e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 703e74616beSVaclav Hapla 704874270d9SVaclav Hapla Level: intermediate 705874270d9SVaclav Hapla 706874270d9SVaclav Hapla .seealso: PetscViewerHDF5Open(),PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 707874270d9SVaclav Hapla @*/ 70854dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 70954dbf706SVaclav Hapla { 71090e03537SVaclav Hapla hid_t file_id; 71190e03537SVaclav Hapla H5O_type_t type; 71254dbf706SVaclav Hapla const char *groupName = NULL; 713e74616beSVaclav Hapla PetscBool create; 71454dbf706SVaclav Hapla PetscErrorCode ierr; 71554dbf706SVaclav Hapla 71654dbf706SVaclav Hapla PetscFunctionBegin; 717e74616beSVaclav Hapla ierr = PetscViewerWritable(viewer, &create);CHKERRQ(ierr); 71854dbf706SVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &file_id);CHKERRQ(ierr); 71954dbf706SVaclav Hapla ierr = PetscViewerHDF5GetGroup(viewer, &groupName);CHKERRQ(ierr); 720e74616beSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, groupName, create, NULL, &type);CHKERRQ(ierr); 72190e03537SVaclav 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); 72290e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 72354dbf706SVaclav Hapla *fileId = file_id; 72454dbf706SVaclav Hapla PetscFunctionReturn(0); 72554dbf706SVaclav Hapla } 72654dbf706SVaclav Hapla 7273ef9c667SSatish Balay /*@ 728d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 7295c6c1daeSBarry Smith 7305c6c1daeSBarry Smith Not collective 7315c6c1daeSBarry Smith 7325c6c1daeSBarry Smith Input Parameter: 7335c6c1daeSBarry Smith . viewer - the PetscViewer 7345c6c1daeSBarry Smith 7355c6c1daeSBarry Smith Level: intermediate 7365c6c1daeSBarry Smith 737d7dd068bSVaclav Hapla Notes: 738d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 739d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 740d7dd068bSVaclav Hapla Current timestep value determines which timestep is read from or written to any dataset on the next HDF5 I/O operation [e.g. VecView()]. 741d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 742d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 743d7dd068bSVaclav Hapla 744d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 745d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 746d7dd068bSVaclav Hapla 747d7dd068bSVaclav Hapla Developer notes: 748d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 749d7dd068bSVaclav Hapla 750d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 751d7dd068bSVaclav Hapla @*/ 752d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 753d7dd068bSVaclav Hapla { 754d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 755d7dd068bSVaclav Hapla 756d7dd068bSVaclav Hapla PetscFunctionBegin; 757d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 758d7dd068bSVaclav Hapla if (hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 759d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 760d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 761d7dd068bSVaclav Hapla PetscFunctionReturn(0); 762d7dd068bSVaclav Hapla } 763d7dd068bSVaclav Hapla 764d7dd068bSVaclav Hapla /*@ 765d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 766d7dd068bSVaclav Hapla 767d7dd068bSVaclav Hapla Not collective 768d7dd068bSVaclav Hapla 769d7dd068bSVaclav Hapla Input Parameter: 770d7dd068bSVaclav Hapla . viewer - the PetscViewer 771d7dd068bSVaclav Hapla 772d7dd068bSVaclav Hapla Level: intermediate 773d7dd068bSVaclav Hapla 774d7dd068bSVaclav Hapla Notes: 775d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 776d7dd068bSVaclav Hapla 777d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IsTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 778d7dd068bSVaclav Hapla @*/ 779d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 780d7dd068bSVaclav Hapla { 781d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 782d7dd068bSVaclav Hapla 783d7dd068bSVaclav Hapla PetscFunctionBegin; 784d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 785d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 786d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 787d7dd068bSVaclav Hapla PetscFunctionReturn(0); 788d7dd068bSVaclav Hapla } 789d7dd068bSVaclav Hapla 790d7dd068bSVaclav Hapla /*@ 791d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 792d7dd068bSVaclav Hapla 793d7dd068bSVaclav Hapla Not collective 794d7dd068bSVaclav Hapla 795d7dd068bSVaclav Hapla Input Parameter: 796d7dd068bSVaclav Hapla . viewer - the PetscViewer 797d7dd068bSVaclav Hapla 798d7dd068bSVaclav Hapla Output Parameter: 799d7dd068bSVaclav Hapla . flg - is timestepping active? 800d7dd068bSVaclav Hapla 801d7dd068bSVaclav Hapla Level: intermediate 802d7dd068bSVaclav Hapla 803d7dd068bSVaclav Hapla Notes: 804d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 805d7dd068bSVaclav Hapla 806d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5PopTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 807d7dd068bSVaclav Hapla @*/ 808d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 809d7dd068bSVaclav Hapla { 810d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 811d7dd068bSVaclav Hapla 812d7dd068bSVaclav Hapla PetscFunctionBegin; 813d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 814d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 815d7dd068bSVaclav Hapla PetscFunctionReturn(0); 816d7dd068bSVaclav Hapla } 817d7dd068bSVaclav Hapla 818d7dd068bSVaclav Hapla /*@ 819d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 820d7dd068bSVaclav Hapla 821d7dd068bSVaclav Hapla Not collective 822d7dd068bSVaclav Hapla 823d7dd068bSVaclav Hapla Input Parameter: 824d7dd068bSVaclav Hapla . viewer - the PetscViewer 825d7dd068bSVaclav Hapla 826d7dd068bSVaclav Hapla Level: intermediate 827d7dd068bSVaclav Hapla 828d7dd068bSVaclav Hapla Notes: 829d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 830d7dd068bSVaclav Hapla 831d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5SetTimestep(), PetscViewerHDF5GetTimestep() 8325c6c1daeSBarry Smith @*/ 8335c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 8345c6c1daeSBarry Smith { 8355c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8365c6c1daeSBarry Smith 8375c6c1daeSBarry Smith PetscFunctionBegin; 8385c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 839d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8405c6c1daeSBarry Smith ++hdf5->timestep; 8415c6c1daeSBarry Smith PetscFunctionReturn(0); 8425c6c1daeSBarry Smith } 8435c6c1daeSBarry Smith 8443ef9c667SSatish Balay /*@ 845d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 8465c6c1daeSBarry Smith 847d7dd068bSVaclav Hapla Logically collective 8485c6c1daeSBarry Smith 8495c6c1daeSBarry Smith Input Parameters: 8505c6c1daeSBarry Smith + viewer - the PetscViewer 851d7dd068bSVaclav Hapla - timestep - The timestep 8525c6c1daeSBarry Smith 8535c6c1daeSBarry Smith Level: intermediate 8545c6c1daeSBarry Smith 855d7dd068bSVaclav Hapla Notes: 856d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 857d7dd068bSVaclav Hapla 858d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5GetTimestep() 8595c6c1daeSBarry Smith @*/ 8605c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 8615c6c1daeSBarry Smith { 8625c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8635c6c1daeSBarry Smith 8645c6c1daeSBarry Smith PetscFunctionBegin; 8655c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 866d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 867d7dd068bSVaclav Hapla if (timestep < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %D is negative", timestep); 868d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8695c6c1daeSBarry Smith hdf5->timestep = timestep; 8705c6c1daeSBarry Smith PetscFunctionReturn(0); 8715c6c1daeSBarry Smith } 8725c6c1daeSBarry Smith 8733ef9c667SSatish Balay /*@ 8745c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 8755c6c1daeSBarry Smith 8765c6c1daeSBarry Smith Not collective 8775c6c1daeSBarry Smith 8785c6c1daeSBarry Smith Input Parameter: 8795c6c1daeSBarry Smith . viewer - the PetscViewer 8805c6c1daeSBarry Smith 8815c6c1daeSBarry Smith Output Parameter: 882d7dd068bSVaclav Hapla . timestep - The timestep 8835c6c1daeSBarry Smith 8845c6c1daeSBarry Smith Level: intermediate 8855c6c1daeSBarry Smith 886d7dd068bSVaclav Hapla Notes: 887d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 888d7dd068bSVaclav Hapla 889d7dd068bSVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5PushTimestepping(), PetscViewerHDF5IncrementTimestep(), PetscViewerHDF5SetTimestep() 8905c6c1daeSBarry Smith @*/ 8915c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 8925c6c1daeSBarry Smith { 8935c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 8945c6c1daeSBarry Smith 8955c6c1daeSBarry Smith PetscFunctionBegin; 8965c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 897d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 898d7dd068bSVaclav Hapla if (!hdf5->timestepping) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 8995c6c1daeSBarry Smith *timestep = hdf5->timestep; 9005c6c1daeSBarry Smith PetscFunctionReturn(0); 9015c6c1daeSBarry Smith } 9025c6c1daeSBarry Smith 90336bce6e8SMatthew G. Knepley /*@C 90436bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 90536bce6e8SMatthew G. Knepley 90636bce6e8SMatthew G. Knepley Not collective 90736bce6e8SMatthew G. Knepley 90836bce6e8SMatthew G. Knepley Input Parameter: 90936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 91036bce6e8SMatthew G. Knepley 91136bce6e8SMatthew G. Knepley Output Parameter: 91236bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 91336bce6e8SMatthew G. Knepley 91436bce6e8SMatthew G. Knepley Level: advanced 91536bce6e8SMatthew G. Knepley 91636bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 91736bce6e8SMatthew G. Knepley @*/ 91836bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 91936bce6e8SMatthew G. Knepley { 92036bce6e8SMatthew G. Knepley PetscFunctionBegin; 92136bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 92236bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 92336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 92436bce6e8SMatthew G. Knepley #else 92536bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 92636bce6e8SMatthew G. Knepley #endif 92736bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 92836bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 92936bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 930c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 931de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 93236bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 93336bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 93436bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 9357e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 93636bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 93736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 93836bce6e8SMatthew G. Knepley } 93936bce6e8SMatthew G. Knepley 94036bce6e8SMatthew G. Knepley /*@C 94136bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 94236bce6e8SMatthew G. Knepley 94336bce6e8SMatthew G. Knepley Not collective 94436bce6e8SMatthew G. Knepley 94536bce6e8SMatthew G. Knepley Input Parameter: 94636bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 94736bce6e8SMatthew G. Knepley 94836bce6e8SMatthew G. Knepley Output Parameter: 94936bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 95036bce6e8SMatthew G. Knepley 95136bce6e8SMatthew G. Knepley Level: advanced 95236bce6e8SMatthew G. Knepley 95336bce6e8SMatthew G. Knepley .seealso: PetscDataType, PetscHDF5DataTypeToPetscDataType() 95436bce6e8SMatthew G. Knepley @*/ 95536bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 95636bce6e8SMatthew G. Knepley { 95736bce6e8SMatthew G. Knepley PetscFunctionBegin; 95836bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 95936bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 96036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 96136bce6e8SMatthew G. Knepley #else 96236bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 96336bce6e8SMatthew G. Knepley #endif 96436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 96536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 96636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 96736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 96836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 96936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 9707e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 97136bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 97236bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 97336bce6e8SMatthew G. Knepley } 97436bce6e8SMatthew G. Knepley 975df863907SAlex Fikl /*@C 976b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 97736bce6e8SMatthew G. Knepley 978343a20b2SVaclav Hapla Collective 979343a20b2SVaclav Hapla 98036bce6e8SMatthew G. Knepley Input Parameters: 98136bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 9824302210dSVaclav Hapla . parent - The parent dataset/group name 98336bce6e8SMatthew G. Knepley . name - The attribute name 98436bce6e8SMatthew G. Knepley . datatype - The attribute type 98536bce6e8SMatthew G. Knepley - value - The attribute value 98636bce6e8SMatthew G. Knepley 98736bce6e8SMatthew G. Knepley Level: advanced 98836bce6e8SMatthew G. Knepley 9894302210dSVaclav Hapla Notes: 990343a20b2SVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 9914302210dSVaclav Hapla 992577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 99336bce6e8SMatthew G. Knepley @*/ 9944302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 99536bce6e8SMatthew G. Knepley { 9964302210dSVaclav Hapla char *parentAbsPath; 99760568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 998080f144cSVaclav Hapla PetscBool has; 99936bce6e8SMatthew G. Knepley PetscErrorCode ierr; 100036bce6e8SMatthew G. Knepley 100136bce6e8SMatthew G. Knepley PetscFunctionBegin; 10025cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10034302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1004c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 10054302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1006b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 10074302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 10084302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL);CHKERRQ(ierr); 10094302210dSVaclav Hapla ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr); 101036bce6e8SMatthew G. Knepley ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 10117e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 10127e97c476SMatthew G. Knepley size_t len; 10137e97c476SMatthew G. Knepley ierr = PetscStrlen((const char *) value, &len);CHKERRQ(ierr); 1014729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 10157e97c476SMatthew G. Knepley } 101636bce6e8SMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 1017729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 10184302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1019080f144cSVaclav Hapla if (has) { 1020080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1021080f144cSVaclav Hapla } else { 102260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1023080f144cSVaclav Hapla } 1024729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1025729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1026729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 102760568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1028729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 10294302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 103036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 103136bce6e8SMatthew G. Knepley } 103236bce6e8SMatthew G. Knepley 1033df863907SAlex Fikl /*@C 1034577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1035577e0e04SVaclav Hapla 1036343a20b2SVaclav Hapla Collective 1037343a20b2SVaclav Hapla 1038577e0e04SVaclav Hapla Input Parameters: 1039577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1040577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1041577e0e04SVaclav Hapla . name - The attribute name 1042577e0e04SVaclav Hapla . datatype - The attribute type 1043577e0e04SVaclav Hapla - value - The attribute value 1044577e0e04SVaclav Hapla 1045577e0e04SVaclav Hapla Notes: 1046343a20b2SVaclav Hapla This fails if the path current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1047577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1048577e0e04SVaclav Hapla 1049577e0e04SVaclav Hapla Level: advanced 1050577e0e04SVaclav Hapla 1051577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1052577e0e04SVaclav Hapla @*/ 1053577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1054577e0e04SVaclav Hapla { 1055577e0e04SVaclav Hapla PetscErrorCode ierr; 1056577e0e04SVaclav Hapla 1057577e0e04SVaclav Hapla PetscFunctionBegin; 1058577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1059577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1060577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1061b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 1062577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1063577e0e04SVaclav Hapla ierr = PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value);CHKERRQ(ierr); 1064577e0e04SVaclav Hapla PetscFunctionReturn(0); 1065577e0e04SVaclav Hapla } 1066577e0e04SVaclav Hapla 1067577e0e04SVaclav Hapla /*@C 1068b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1069ad0c61feSMatthew G. Knepley 1070343a20b2SVaclav Hapla Collective 1071343a20b2SVaclav Hapla 1072ad0c61feSMatthew G. Knepley Input Parameters: 1073ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 10744302210dSVaclav Hapla . parent - The parent dataset/group name 1075ad0c61feSMatthew G. Knepley . name - The attribute name 1076a2d6be1bSVaclav Hapla . datatype - The attribute type 1077a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1078ad0c61feSMatthew G. Knepley 1079ad0c61feSMatthew G. Knepley Output Parameter: 1080a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1081ad0c61feSMatthew G. Knepley 1082a2d6be1bSVaclav Hapla Notes: 1083a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1084a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1085a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1086a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 1087a2d6be1bSVaclav Hapla $ ierr = PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg);CHKERRQ(ierr); 1088a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1089a2d6be1bSVaclav Hapla 1090a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1091708d4cb3SBarry Smith 1092343a20b2SVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 1093ad0c61feSMatthew G. Knepley 1094343a20b2SVaclav Hapla Level: advanced 10954302210dSVaclav Hapla 1096577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1097ad0c61feSMatthew G. Knepley @*/ 1098d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1099ad0c61feSMatthew G. Knepley { 11004302210dSVaclav Hapla char *parentAbsPath; 1101a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1102969748fdSVaclav Hapla PetscBool has; 1103ad0c61feSMatthew G. Knepley PetscErrorCode ierr; 1104ad0c61feSMatthew G. Knepley 1105ad0c61feSMatthew G. Knepley PetscFunctionBegin; 11065cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11074302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1108c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1109a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1110a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 1111a2d6be1bSVaclav Hapla ierr = PetscDataTypeToHDF5DataType(datatype, &dtype);CHKERRQ(ierr); 11124302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 11134302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL);CHKERRQ(ierr); 11144302210dSVaclav Hapla if (has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has);CHKERRQ(ierr);} 1115a2d6be1bSVaclav Hapla if (!has) { 1116a2d6be1bSVaclav Hapla if (defaultValue) { 1117a2d6be1bSVaclav Hapla if (defaultValue != value) { 1118a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 1119a2d6be1bSVaclav Hapla ierr = PetscStrallocpy(*(char**)defaultValue, (char**)value);CHKERRQ(ierr); 1120a2d6be1bSVaclav Hapla } else { 1121a2d6be1bSVaclav Hapla size_t len; 1122a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(dtype)); 1123a2d6be1bSVaclav Hapla ierr = PetscMemcpy(value, defaultValue, len);CHKERRQ(ierr); 1124a2d6be1bSVaclav Hapla } 1125a2d6be1bSVaclav Hapla } 1126a2d6be1bSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1127a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 1128a2d6be1bSVaclav Hapla } else SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1129a2d6be1bSVaclav Hapla } 1130ad0c61feSMatthew G. Knepley ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 11314302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 113260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1133f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1134f0b58479SMatthew G. Knepley size_t len; 1135a2d6be1bSVaclav Hapla hid_t atype; 1136e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1137a2d6be1bSVaclav Hapla PetscStackCallHDF5Return(len,H5Tget_size,(atype)); 1138708d4cb3SBarry Smith ierr = PetscMalloc((len+1) * sizeof(char), value);CHKERRQ(ierr); 1139708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1140708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1141708d4cb3SBarry Smith } else { 114270efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1143708d4cb3SBarry Smith } 1144729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1145e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1146e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 11474302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 1148ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1149ad0c61feSMatthew G. Knepley } 1150ad0c61feSMatthew G. Knepley 1151577e0e04SVaclav Hapla /*@C 1152577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1153577e0e04SVaclav Hapla 1154343a20b2SVaclav Hapla Collective 1155343a20b2SVaclav Hapla 1156577e0e04SVaclav Hapla Input Parameters: 1157577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1158577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1159577e0e04SVaclav Hapla . name - The attribute name 1160577e0e04SVaclav Hapla - datatype - The attribute type 1161577e0e04SVaclav Hapla 1162577e0e04SVaclav Hapla Output Parameter: 1163577e0e04SVaclav Hapla . value - The attribute value 1164577e0e04SVaclav Hapla 1165577e0e04SVaclav Hapla Notes: 1166577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1167577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1168577e0e04SVaclav Hapla 1169577e0e04SVaclav Hapla Level: advanced 1170577e0e04SVaclav Hapla 1171577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5ReadAttribute() PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1172577e0e04SVaclav Hapla @*/ 1173a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1174577e0e04SVaclav Hapla { 1175577e0e04SVaclav Hapla PetscErrorCode ierr; 1176577e0e04SVaclav Hapla 1177577e0e04SVaclav Hapla PetscFunctionBegin; 1178577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1179577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1180577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1181b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 1182577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1183a2d6be1bSVaclav Hapla ierr = PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value);CHKERRQ(ierr); 1184577e0e04SVaclav Hapla PetscFunctionReturn(0); 1185577e0e04SVaclav Hapla } 1186577e0e04SVaclav Hapla 1187a75c3fd4SVaclav Hapla PETSC_STATIC_INLINE PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1188a75c3fd4SVaclav Hapla { 1189a75c3fd4SVaclav Hapla htri_t exists; 1190a75c3fd4SVaclav Hapla hid_t group; 1191a75c3fd4SVaclav Hapla 1192a75c3fd4SVaclav Hapla PetscFunctionBegin; 1193a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1194a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1195a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1196a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1197a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1198a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1199a75c3fd4SVaclav Hapla } 1200a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1201a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1202a75c3fd4SVaclav Hapla } 1203a75c3fd4SVaclav Hapla 1204bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 12055cdeb108SMatthew G. Knepley { 120690e03537SVaclav Hapla const char rootGroupName[] = "/"; 12075cdeb108SMatthew G. Knepley hid_t h5; 1208e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 120915b861d2SVaclav Hapla PetscInt i; 121015b861d2SVaclav Hapla int n; 121185688ae2SVaclav Hapla char **hierarchy; 121285688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 12135cdeb108SMatthew G. Knepley PetscErrorCode ierr; 12145cdeb108SMatthew G. Knepley 12155cdeb108SMatthew G. Knepley PetscFunctionBegin; 12165cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 121790e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 121890e03537SVaclav Hapla else name = rootGroupName; 1219ccd66a83SVaclav Hapla if (has) { 122056cc0592SVaclav Hapla PetscValidIntPointer(has, 3); 12215cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1222ccd66a83SVaclav Hapla } 1223ccd66a83SVaclav Hapla if (otype) { 1224ccd66a83SVaclav Hapla PetscValidIntPointer(otype, 4); 122556cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1226ccd66a83SVaclav Hapla } 1227c0cd0301SJed Brown ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 122885688ae2SVaclav Hapla 122985688ae2SVaclav Hapla /* 123085688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 123185688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 123285688ae2SVaclav Hapla 1) whether it's a valid link 123385688ae2SVaclav Hapla 2) whether this link resolves to an object 123485688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 123585688ae2SVaclav Hapla */ 123685688ae2SVaclav Hapla ierr = PetscStrToArray(name,'/',&n,&hierarchy);CHKERRQ(ierr); 123785688ae2SVaclav Hapla if (!n) { 123885688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1239ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1240ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 124185688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 124285688ae2SVaclav Hapla PetscFunctionReturn(0); 124385688ae2SVaclav Hapla } 124485688ae2SVaclav Hapla for (i=0; i<n; i++) { 124585688ae2SVaclav Hapla ierr = PetscStrcat(buf,"/");CHKERRQ(ierr); 124685688ae2SVaclav Hapla ierr = PetscStrcat(buf,hierarchy[i]);CHKERRQ(ierr); 1247a75c3fd4SVaclav Hapla ierr = PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists);CHKERRQ(ierr); 1248a75c3fd4SVaclav Hapla if (!exists) break; 124985688ae2SVaclav Hapla } 125085688ae2SVaclav Hapla ierr = PetscStrToArrayDestroy(n,hierarchy);CHKERRQ(ierr); 125185688ae2SVaclav Hapla 125285688ae2SVaclav Hapla /* If the object exists, get its type */ 1253ccd66a83SVaclav Hapla if (exists && otype) { 12545cdeb108SMatthew G. Knepley H5O_info_t info; 12555cdeb108SMatthew G. Knepley 125676276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 125704633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 125856cc0592SVaclav Hapla *otype = info.type; 12595cdeb108SMatthew G. Knepley } 1260ccd66a83SVaclav Hapla if (has) *has = exists; 12615cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 12625cdeb108SMatthew G. Knepley } 12635cdeb108SMatthew G. Knepley 12644302210dSVaclav Hapla /*@C 12650a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 12660a9f38efSVaclav Hapla 1267343a20b2SVaclav Hapla Collective 1268343a20b2SVaclav Hapla 12690a9f38efSVaclav Hapla Input Parameters: 12704302210dSVaclav Hapla + viewer - The HDF5 viewer 12714302210dSVaclav Hapla - path - The group path 12720a9f38efSVaclav Hapla 12730a9f38efSVaclav Hapla Output Parameter: 12740a9f38efSVaclav Hapla . has - Flag for group existence 12750a9f38efSVaclav Hapla 12760a9f38efSVaclav Hapla Level: advanced 12770a9f38efSVaclav Hapla 12784302210dSVaclav Hapla Notes: 1279*785c443eSVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 1280*785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1281343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 12824302210dSVaclav Hapla 128389e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasDataset(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup(), PetscViewerHDF5OpenGroup() 12840a9f38efSVaclav Hapla @*/ 12854302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 12860a9f38efSVaclav Hapla { 12870a9f38efSVaclav Hapla H5O_type_t type; 12884302210dSVaclav Hapla char *abspath; 12890a9f38efSVaclav Hapla PetscErrorCode ierr; 12900a9f38efSVaclav Hapla 12910a9f38efSVaclav Hapla PetscFunctionBegin; 12920a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 12934302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 12944302210dSVaclav Hapla PetscValidBoolPointer(has,3); 12954302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 12964302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 12974302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 12984302210dSVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 12990a9f38efSVaclav Hapla PetscFunctionReturn(0); 13000a9f38efSVaclav Hapla } 13010a9f38efSVaclav Hapla 130289e0ef10SVaclav Hapla /*@C 130389e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 130489e0ef10SVaclav Hapla 1305343a20b2SVaclav Hapla Collective 1306343a20b2SVaclav Hapla 130789e0ef10SVaclav Hapla Input Parameters: 130889e0ef10SVaclav Hapla + viewer - The HDF5 viewer 130989e0ef10SVaclav Hapla - path - The dataset path 131089e0ef10SVaclav Hapla 131189e0ef10SVaclav Hapla Output Parameter: 131289e0ef10SVaclav Hapla . has - Flag whether dataset exists 131389e0ef10SVaclav Hapla 131489e0ef10SVaclav Hapla Level: advanced 131589e0ef10SVaclav Hapla 131689e0ef10SVaclav Hapla Notes: 1317343a20b2SVaclav Hapla If path starts with '/', it is taken as an absolute path overriding currently pushed group, else path is relative to the current pushed group. 131889e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1319343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 132089e0ef10SVaclav Hapla 132189e0ef10SVaclav Hapla .seealso: PetscViewerHDF5HasObject(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5HasGroup(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 132289e0ef10SVaclav Hapla @*/ 132389e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 132489e0ef10SVaclav Hapla { 132589e0ef10SVaclav Hapla H5O_type_t type; 132689e0ef10SVaclav Hapla char *abspath; 132789e0ef10SVaclav Hapla PetscErrorCode ierr; 132889e0ef10SVaclav Hapla 132989e0ef10SVaclav Hapla PetscFunctionBegin; 133089e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 133189e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 133289e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 133389e0ef10SVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath);CHKERRQ(ierr); 133489e0ef10SVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type);CHKERRQ(ierr); 133589e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 133689e0ef10SVaclav Hapla ierr = PetscFree(abspath);CHKERRQ(ierr); 133789e0ef10SVaclav Hapla PetscFunctionReturn(0); 133889e0ef10SVaclav Hapla } 133989e0ef10SVaclav Hapla 13400a9f38efSVaclav Hapla /*@ 1341e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1342ecce1506SVaclav Hapla 1343343a20b2SVaclav Hapla Collective 1344343a20b2SVaclav Hapla 1345ecce1506SVaclav Hapla Input Parameters: 1346ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1347ecce1506SVaclav Hapla - obj - The named object 1348ecce1506SVaclav Hapla 1349ecce1506SVaclav Hapla Output Parameter: 135089e0ef10SVaclav Hapla . has - Flag for dataset existence 1351ecce1506SVaclav Hapla 1352e3f143f7SVaclav Hapla Notes: 135389e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1354343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1355e3f143f7SVaclav Hapla 1356ecce1506SVaclav Hapla Level: advanced 1357ecce1506SVaclav Hapla 135889e0ef10SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasDataset(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5PushGroup(), PetscViewerHDF5PopGroup(), PetscViewerHDF5GetGroup() 1359ecce1506SVaclav Hapla @*/ 1360ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1361ecce1506SVaclav Hapla { 136289e0ef10SVaclav Hapla size_t len; 1363ecce1506SVaclav Hapla PetscErrorCode ierr; 1364ecce1506SVaclav Hapla 1365ecce1506SVaclav Hapla PetscFunctionBegin; 1366c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1367c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 13684302210dSVaclav Hapla PetscValidBoolPointer(has,3); 136989e0ef10SVaclav Hapla ierr = PetscStrlen(obj->name, &len);CHKERRQ(ierr); 137089e0ef10SVaclav Hapla if (!len) SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 137189e0ef10SVaclav Hapla ierr = PetscViewerHDF5HasDataset(viewer, obj->name, has);CHKERRQ(ierr); 1372ecce1506SVaclav Hapla PetscFunctionReturn(0); 1373ecce1506SVaclav Hapla } 1374ecce1506SVaclav Hapla 1375df863907SAlex Fikl /*@C 1376ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1377e7bdbf8eSMatthew G. Knepley 1378343a20b2SVaclav Hapla Collective 1379343a20b2SVaclav Hapla 1380e7bdbf8eSMatthew G. Knepley Input Parameters: 1381e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 13824302210dSVaclav Hapla . parent - The parent dataset/group name 1383e7bdbf8eSMatthew G. Knepley - name - The attribute name 1384e7bdbf8eSMatthew G. Knepley 1385e7bdbf8eSMatthew G. Knepley Output Parameter: 1386e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1387e7bdbf8eSMatthew G. Knepley 1388e7bdbf8eSMatthew G. Knepley Level: advanced 1389e7bdbf8eSMatthew G. Knepley 13904302210dSVaclav Hapla Notes: 1391343a20b2SVaclav Hapla If parent starts with '/', it is taken as an absolute path overriding currently pushed group, else parent is relative to the current pushed group. NULL means the current pushed group. 13924302210dSVaclav Hapla 1393577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasObjectAttribute(), PetscViewerHDF5WriteAttribute(), PetscViewerHDF5ReadAttribute(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1394e7bdbf8eSMatthew G. Knepley @*/ 13954302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1396e7bdbf8eSMatthew G. Knepley { 13974302210dSVaclav Hapla char *parentAbsPath; 1398e7bdbf8eSMatthew G. Knepley PetscErrorCode ierr; 1399e7bdbf8eSMatthew G. Knepley 1400e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 14015cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14024302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1403c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 14044302210dSVaclav Hapla PetscValidBoolPointer(has,4); 14054302210dSVaclav Hapla ierr = PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath);CHKERRQ(ierr); 14064302210dSVaclav Hapla ierr = PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL);CHKERRQ(ierr); 14074302210dSVaclav Hapla if (*has) {ierr = PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has);CHKERRQ(ierr);} 14084302210dSVaclav Hapla ierr = PetscFree(parentAbsPath);CHKERRQ(ierr); 140906db490cSVaclav Hapla PetscFunctionReturn(0); 141006db490cSVaclav Hapla } 141106db490cSVaclav Hapla 1412577e0e04SVaclav Hapla /*@C 1413577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1414577e0e04SVaclav Hapla 1415343a20b2SVaclav Hapla Collective 1416343a20b2SVaclav Hapla 1417577e0e04SVaclav Hapla Input Parameters: 1418577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1419577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1420577e0e04SVaclav Hapla - name - The attribute name 1421577e0e04SVaclav Hapla 1422577e0e04SVaclav Hapla Output Parameter: 1423577e0e04SVaclav Hapla . has - Flag for attribute existence 1424577e0e04SVaclav Hapla 1425577e0e04SVaclav Hapla Notes: 1426577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1427577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1428577e0e04SVaclav Hapla 1429577e0e04SVaclav Hapla Level: advanced 1430577e0e04SVaclav Hapla 1431577e0e04SVaclav Hapla .seealso: PetscViewerHDF5Open(), PetscViewerHDF5HasAttribute(), PetscViewerHDF5WriteObjectAttribute(), PetscViewerHDF5ReadObjectAttribute(), PetscViewerHDF5HasObject(), PetscViewerHDF5PushGroup(),PetscViewerHDF5PopGroup(),PetscViewerHDF5GetGroup() 1432577e0e04SVaclav Hapla @*/ 1433577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1434577e0e04SVaclav Hapla { 1435577e0e04SVaclav Hapla PetscErrorCode ierr; 1436577e0e04SVaclav Hapla 1437577e0e04SVaclav Hapla PetscFunctionBegin; 1438577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1439577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1440577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 14414302210dSVaclav Hapla PetscValidBoolPointer(has,4); 1442577e0e04SVaclav Hapla ierr = PetscViewerHDF5CheckNamedObject_Internal(viewer, obj);CHKERRQ(ierr); 1443577e0e04SVaclav Hapla ierr = PetscViewerHDF5HasAttribute(viewer, obj->name, name, has);CHKERRQ(ierr); 1444577e0e04SVaclav Hapla PetscFunctionReturn(0); 1445577e0e04SVaclav Hapla } 1446577e0e04SVaclav Hapla 144706db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 144806db490cSVaclav Hapla { 144906db490cSVaclav Hapla hid_t h5; 145006db490cSVaclav Hapla htri_t hhas; 145106db490cSVaclav Hapla PetscErrorCode ierr; 145206db490cSVaclav Hapla 145306db490cSVaclav Hapla PetscFunctionBegin; 145406db490cSVaclav Hapla ierr = PetscViewerHDF5GetFileId(viewer, &h5);CHKERRQ(ierr); 14552f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 14565cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1457e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1458e7bdbf8eSMatthew G. Knepley } 1459e7bdbf8eSMatthew G. Knepley 1460a75e6a4aSMatthew G. Knepley /* 1461a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1462a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1463a75e6a4aSMatthew G. Knepley */ 1464d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1465a75e6a4aSMatthew G. Knepley 1466a75e6a4aSMatthew G. Knepley /*@C 1467a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1468a75e6a4aSMatthew G. Knepley 1469d083f849SBarry Smith Collective 1470a75e6a4aSMatthew G. Knepley 1471a75e6a4aSMatthew G. Knepley Input Parameter: 1472a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1473a75e6a4aSMatthew G. Knepley 1474a75e6a4aSMatthew G. Knepley Level: intermediate 1475a75e6a4aSMatthew G. Knepley 1476a75e6a4aSMatthew G. Knepley Options Database Keys: 1477a75e6a4aSMatthew G. Knepley . -viewer_hdf5_filename <name> 1478a75e6a4aSMatthew G. Knepley 1479a75e6a4aSMatthew G. Knepley Environmental variables: 1480a75e6a4aSMatthew G. Knepley . PETSC_VIEWER_HDF5_FILENAME 1481a75e6a4aSMatthew G. Knepley 1482a75e6a4aSMatthew G. Knepley Notes: 1483a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1484a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1485a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1486a75e6a4aSMatthew G. Knepley 1487a75e6a4aSMatthew G. Knepley .seealso: PetscViewerHDF5Open(), PetscViewerCreate(), PetscViewerDestroy() 1488a75e6a4aSMatthew G. Knepley @*/ 1489a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1490a75e6a4aSMatthew G. Knepley { 1491a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1492a75e6a4aSMatthew G. Knepley PetscBool flg; 1493a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1494a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1495a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1496a75e6a4aSMatthew G. Knepley 1497a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1498908793a3SLisandro Dalcin ierr = PetscCommDuplicate(comm,&ncomm,NULL);if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1499a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1500908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1501908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1502a75e6a4aSMatthew G. Knepley } 150347435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1504908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1505a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1506a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 15072cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1508a75e6a4aSMatthew G. Knepley if (!flg) { 1509a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 15102cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1511a75e6a4aSMatthew G. Knepley } 1512a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 15132cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1514a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 15152cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 151647435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1517908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1518a75e6a4aSMatthew G. Knepley } 1519a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 15202cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1521a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1522a75e6a4aSMatthew G. Knepley } 1523