16b9fdc0aSVaclav Hapla #include <petsc/private/viewerhdf5impl.h> /*I "petscviewerhdf5.h" I*/ 25c6c1daeSBarry Smith 3bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer, const char[], PetscBool, PetscBool*, H5O_type_t*); 406db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer, const char[], const char[], PetscBool*); 506db490cSVaclav Hapla 64302210dSVaclav Hapla static PetscErrorCode PetscViewerHDF5GetAbsolutePath_Internal(PetscViewer viewer, const char path[], char *abspath[]) 76c132bc1SVaclav Hapla { 84302210dSVaclav Hapla PetscBool relative = PETSC_FALSE; 96c132bc1SVaclav Hapla const char *group; 106c132bc1SVaclav Hapla char buf[PETSC_MAX_PATH_LEN] = ""; 116c132bc1SVaclav Hapla 126c132bc1SVaclav Hapla PetscFunctionBegin; 139566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetGroup(viewer, &group)); 149566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5PathIsRelative(path, PETSC_TRUE, &relative)); 154302210dSVaclav Hapla if (relative) { 169566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(buf, group)); 179566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf, "/")); 189566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf, path)); 199566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(buf, abspath)); 204302210dSVaclav Hapla } else { 219566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(path, abspath)); 224302210dSVaclav Hapla } 236c132bc1SVaclav Hapla PetscFunctionReturn(0); 246c132bc1SVaclav Hapla } 256c132bc1SVaclav Hapla 26577e0e04SVaclav Hapla static PetscErrorCode PetscViewerHDF5CheckNamedObject_Internal(PetscViewer viewer, PetscObject obj) 27577e0e04SVaclav Hapla { 28577e0e04SVaclav Hapla PetscBool has; 29577e0e04SVaclav Hapla 30577e0e04SVaclav Hapla PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasObject(viewer, obj, &has)); 32577e0e04SVaclav Hapla if (!has) { 3389e0ef10SVaclav Hapla const char *group; 349566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetGroup(viewer, &group)); 3598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Object (dataset) \"%s\" not stored in group %s", obj->name, group ? group : "/"); 36577e0e04SVaclav Hapla } 37577e0e04SVaclav Hapla PetscFunctionReturn(0); 38577e0e04SVaclav Hapla } 39577e0e04SVaclav Hapla 404416b707SBarry Smith static PetscErrorCode PetscViewerSetFromOptions_HDF5(PetscOptionItems *PetscOptionsObject,PetscViewer v) 4182ea9b62SBarry Smith { 42ee8b9732SVaclav Hapla PetscBool flg = PETSC_FALSE, set; 4382ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 4482ea9b62SBarry Smith 4582ea9b62SBarry Smith PetscFunctionBegin; 46d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"HDF5 PetscViewer Options"); 479566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_base_dimension2","1d Vectors get 2 dimensions in HDF5","PetscViewerHDF5SetBaseDimension2",hdf5->basedimension2,&hdf5->basedimension2,NULL)); 489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_sp_output","Force data to be written in single precision","PetscViewerHDF5SetSPOutput",hdf5->spoutput,&hdf5->spoutput,NULL)); 499566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_collective","Enable collective transfer mode","PetscViewerHDF5SetCollective",flg,&flg,&set)); 509566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetCollective(v,flg)); 5119a20e4cSMatthew G. Knepley flg = PETSC_FALSE; 529566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-viewer_hdf5_default_timestepping","Set default timestepping state","PetscViewerHDF5SetDefaultTimestepping",flg,&flg,&set)); 539566063dSJacob Faibussowitsch if (set) PetscCall(PetscViewerHDF5SetDefaultTimestepping(v,flg)); 54d0609cedSBarry Smith PetscOptionsHeadEnd(); 5582ea9b62SBarry Smith PetscFunctionReturn(0); 5682ea9b62SBarry Smith } 5782ea9b62SBarry Smith 581b793a25SVaclav Hapla static PetscErrorCode PetscViewerView_HDF5(PetscViewer v,PetscViewer viewer) 591b793a25SVaclav Hapla { 601b793a25SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)v->data; 6103fa8834SVaclav Hapla PetscBool flg; 621b793a25SVaclav Hapla 631b793a25SVaclav Hapla PetscFunctionBegin; 641b793a25SVaclav Hapla if (hdf5->filename) { 659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Filename: %s\n",hdf5->filename)); 661b793a25SVaclav Hapla } 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Vectors with blocksize 1 saved as 2D datasets: %s\n",PetscBools[hdf5->basedimension2])); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Enforce single precision storage: %s\n",PetscBools[hdf5->spoutput])); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetCollective(v,&flg)); 709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"MPI-IO transfer mode: %s\n",flg ? "collective" : "independent")); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Default timestepping: %s\n",PetscBools[hdf5->defTimestepping])); 721b793a25SVaclav Hapla PetscFunctionReturn(0); 731b793a25SVaclav Hapla } 741b793a25SVaclav Hapla 755c6c1daeSBarry Smith static PetscErrorCode PetscViewerFileClose_HDF5(PetscViewer viewer) 765c6c1daeSBarry Smith { 775c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 785c6c1daeSBarry Smith 795c6c1daeSBarry Smith PetscFunctionBegin; 809566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->filename)); 81729ab6d9SBarry Smith if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 825c6c1daeSBarry Smith PetscFunctionReturn(0); 835c6c1daeSBarry Smith } 845c6c1daeSBarry Smith 856226335aSVaclav Hapla static PetscErrorCode PetscViewerFlush_HDF5(PetscViewer viewer) 866226335aSVaclav Hapla { 876226335aSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*)viewer->data; 886226335aSVaclav Hapla 896226335aSVaclav Hapla PetscFunctionBegin; 906226335aSVaclav Hapla if (hdf5->file_id) PetscStackCallHDF5(H5Fflush,(hdf5->file_id, H5F_SCOPE_LOCAL)); 916226335aSVaclav Hapla PetscFunctionReturn(0); 926226335aSVaclav Hapla } 936226335aSVaclav Hapla 949b1403beSVaclav Hapla static PetscErrorCode PetscViewerDestroy_HDF5(PetscViewer viewer) 955c6c1daeSBarry Smith { 965c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 975c6c1daeSBarry Smith 985c6c1daeSBarry Smith PetscFunctionBegin; 999c5072fbSVaclav Hapla PetscStackCallHDF5(H5Pclose,(hdf5->dxpl_id)); 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerFileClose_HDF5(viewer)); 1015c6c1daeSBarry Smith while (hdf5->groups) { 1027d964842SVaclav Hapla PetscViewerHDF5GroupList *tmp = hdf5->groups->next; 1035c6c1daeSBarry Smith 1049566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups->name)); 1059566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5->groups)); 1065c6c1daeSBarry Smith hdf5->groups = tmp; 1075c6c1daeSBarry Smith } 1089566063dSJacob Faibussowitsch PetscCall(PetscFree(hdf5)); 1099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL)); 1109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL)); 1119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL)); 1129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL)); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL)); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetDefaultTimestepping_C",NULL)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetDefaultTimestepping_C",NULL)); 1165c6c1daeSBarry Smith PetscFunctionReturn(0); 1175c6c1daeSBarry Smith } 1185c6c1daeSBarry Smith 1199b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1205c6c1daeSBarry Smith { 1215c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1225c6c1daeSBarry Smith 1235c6c1daeSBarry Smith PetscFunctionBegin; 1245c6c1daeSBarry Smith hdf5->btype = type; 1255c6c1daeSBarry Smith PetscFunctionReturn(0); 1265c6c1daeSBarry Smith } 1275c6c1daeSBarry Smith 1280b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1290b62783dSVaclav Hapla { 1300b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1310b62783dSVaclav Hapla 1320b62783dSVaclav Hapla PetscFunctionBegin; 1330b62783dSVaclav Hapla *type = hdf5->btype; 1340b62783dSVaclav Hapla PetscFunctionReturn(0); 1350b62783dSVaclav Hapla } 1360b62783dSVaclav Hapla 1379b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 13882ea9b62SBarry Smith { 13982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 14082ea9b62SBarry Smith 14182ea9b62SBarry Smith PetscFunctionBegin; 14282ea9b62SBarry Smith hdf5->basedimension2 = flg; 14382ea9b62SBarry Smith PetscFunctionReturn(0); 14482ea9b62SBarry Smith } 14582ea9b62SBarry Smith 146df863907SAlex Fikl /*@ 14782ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 14882ea9b62SBarry Smith dimension of 2. 14982ea9b62SBarry Smith 15082ea9b62SBarry Smith Logically Collective on PetscViewer 15182ea9b62SBarry Smith 15282ea9b62SBarry Smith Input Parameters: 15382ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 15482ea9b62SBarry 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 15582ea9b62SBarry Smith 15682ea9b62SBarry Smith Options Database: 15782ea9b62SBarry 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 15882ea9b62SBarry Smith 15995452b02SPatrick Sanan Notes: 16095452b02SPatrick 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 16182ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 16282ea9b62SBarry Smith 16382ea9b62SBarry Smith Level: intermediate 16482ea9b62SBarry Smith 165db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 16682ea9b62SBarry Smith 16782ea9b62SBarry Smith @*/ 16882ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 16982ea9b62SBarry Smith { 17082ea9b62SBarry Smith PetscFunctionBegin; 17182ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 172cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg)); 17382ea9b62SBarry Smith PetscFunctionReturn(0); 17482ea9b62SBarry Smith } 17582ea9b62SBarry Smith 176df863907SAlex Fikl /*@ 17782ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 17882ea9b62SBarry Smith dimension of 2. 17982ea9b62SBarry Smith 18082ea9b62SBarry Smith Logically Collective on PetscViewer 18182ea9b62SBarry Smith 18282ea9b62SBarry Smith Input Parameter: 18382ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18482ea9b62SBarry Smith 18582ea9b62SBarry Smith Output Parameter: 18682ea9b62SBarry 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 18782ea9b62SBarry Smith 18895452b02SPatrick Sanan Notes: 18995452b02SPatrick 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 19082ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19182ea9b62SBarry Smith 19282ea9b62SBarry Smith Level: intermediate 19382ea9b62SBarry Smith 194db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 19582ea9b62SBarry Smith 19682ea9b62SBarry Smith @*/ 19782ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 19882ea9b62SBarry Smith { 19982ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 20082ea9b62SBarry Smith 20182ea9b62SBarry Smith PetscFunctionBegin; 20282ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 20382ea9b62SBarry Smith *flg = hdf5->basedimension2; 20482ea9b62SBarry Smith PetscFunctionReturn(0); 20582ea9b62SBarry Smith } 20682ea9b62SBarry Smith 2079b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2089a0502c6SHåkon Strandenes { 2099a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2109a0502c6SHåkon Strandenes 2119a0502c6SHåkon Strandenes PetscFunctionBegin; 2129a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2139a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2149a0502c6SHåkon Strandenes } 2159a0502c6SHåkon Strandenes 216df863907SAlex Fikl /*@ 2179a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2189a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2199a0502c6SHåkon Strandenes 2209a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2219a0502c6SHåkon Strandenes 2229a0502c6SHåkon Strandenes Input Parameters: 2239a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2249a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2259a0502c6SHåkon Strandenes 2269a0502c6SHåkon Strandenes Options Database: 2279a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2289a0502c6SHåkon Strandenes 22995452b02SPatrick Sanan Notes: 23095452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2319a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2329a0502c6SHåkon Strandenes 2339a0502c6SHåkon Strandenes Level: intermediate 2349a0502c6SHåkon Strandenes 235db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 236db781477SPatrick Sanan `PetscReal` 2379a0502c6SHåkon Strandenes 2389a0502c6SHåkon Strandenes @*/ 2399a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2409a0502c6SHåkon Strandenes { 2419a0502c6SHåkon Strandenes PetscFunctionBegin; 2429a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 243cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg)); 2449a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2459a0502c6SHåkon Strandenes } 2469a0502c6SHåkon Strandenes 247df863907SAlex Fikl /*@ 2489a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2499a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2509a0502c6SHåkon Strandenes 2519a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2529a0502c6SHåkon Strandenes 2539a0502c6SHåkon Strandenes Input Parameter: 2549a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2559a0502c6SHåkon Strandenes 2569a0502c6SHåkon Strandenes Output Parameter: 2579a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2589a0502c6SHåkon Strandenes 25995452b02SPatrick Sanan Notes: 26095452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2619a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2629a0502c6SHåkon Strandenes 2639a0502c6SHåkon Strandenes Level: intermediate 2649a0502c6SHåkon Strandenes 265db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 266db781477SPatrick Sanan `PetscReal` 2679a0502c6SHåkon Strandenes 2689a0502c6SHåkon Strandenes @*/ 2699a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2709a0502c6SHåkon Strandenes { 2719a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2729a0502c6SHåkon Strandenes 2739a0502c6SHåkon Strandenes PetscFunctionBegin; 2749a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2759a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2769a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2779a0502c6SHåkon Strandenes } 2789a0502c6SHåkon Strandenes 279ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 280ee8b9732SVaclav Hapla { 281ee8b9732SVaclav Hapla PetscFunctionBegin; 282ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 283ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 284c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 285ee8b9732SVaclav Hapla { 286ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 287ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 288ee8b9732SVaclav Hapla } 289ee8b9732SVaclav Hapla #else 2909566063dSJacob Faibussowitsch if (flg) PetscCall(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")); 291ee8b9732SVaclav Hapla #endif 292ee8b9732SVaclav Hapla PetscFunctionReturn(0); 293ee8b9732SVaclav Hapla } 294ee8b9732SVaclav Hapla 295ee8b9732SVaclav Hapla /*@ 296ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 297ee8b9732SVaclav Hapla 298ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 299ee8b9732SVaclav Hapla 300ee8b9732SVaclav Hapla Input Parameters: 301ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 302ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 303ee8b9732SVaclav Hapla 304ee8b9732SVaclav Hapla Options Database: 305ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 306ee8b9732SVaclav Hapla 307ee8b9732SVaclav Hapla Notes: 308ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 30953effed7SBarry 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. 310ee8b9732SVaclav Hapla 311ee8b9732SVaclav Hapla Developer notes: 312ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 313ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 314ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 315ee8b9732SVaclav Hapla 316ee8b9732SVaclav Hapla Level: intermediate 317ee8b9732SVaclav Hapla 318db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 319ee8b9732SVaclav Hapla 320ee8b9732SVaclav Hapla @*/ 321ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 322ee8b9732SVaclav Hapla { 323ee8b9732SVaclav Hapla PetscFunctionBegin; 324ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 325ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 326cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg)); 327ee8b9732SVaclav Hapla PetscFunctionReturn(0); 328ee8b9732SVaclav Hapla } 329ee8b9732SVaclav Hapla 330ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 331ee8b9732SVaclav Hapla { 332c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 333ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 334ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 335c796909eSBarry Smith #endif 336ee8b9732SVaclav Hapla 337ee8b9732SVaclav Hapla PetscFunctionBegin; 338c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 339c796909eSBarry Smith *flg = PETSC_FALSE; 340c796909eSBarry Smith #else 341ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 342ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 343c796909eSBarry Smith #endif 344ee8b9732SVaclav Hapla PetscFunctionReturn(0); 345ee8b9732SVaclav Hapla } 346ee8b9732SVaclav Hapla 347ee8b9732SVaclav Hapla /*@ 348ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 349ee8b9732SVaclav Hapla 350ee8b9732SVaclav Hapla Not Collective 351ee8b9732SVaclav Hapla 352ee8b9732SVaclav Hapla Input Parameters: 353ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 354ee8b9732SVaclav Hapla 355ee8b9732SVaclav Hapla Output Parameters: 356ee8b9732SVaclav Hapla . flg - the flag 357ee8b9732SVaclav Hapla 358ee8b9732SVaclav Hapla Level: intermediate 359ee8b9732SVaclav Hapla 360ee8b9732SVaclav Hapla Notes: 361c796909eSBarry 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. 362ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 363ee8b9732SVaclav Hapla 364db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 365ee8b9732SVaclav Hapla 366ee8b9732SVaclav Hapla @*/ 367ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 368ee8b9732SVaclav Hapla { 369ee8b9732SVaclav Hapla PetscFunctionBegin; 370ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 371534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 372ee8b9732SVaclav Hapla 373cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg)); 374ee8b9732SVaclav Hapla PetscFunctionReturn(0); 375ee8b9732SVaclav Hapla } 376ee8b9732SVaclav Hapla 3779b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3785c6c1daeSBarry Smith { 3795c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3805c6c1daeSBarry Smith hid_t plist_id; 3815c6c1daeSBarry Smith 3825c6c1daeSBarry Smith PetscFunctionBegin; 383f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 3849566063dSJacob Faibussowitsch if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 3859566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &hdf5->filename)); 3865c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 387729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 388c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 389d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 390c796909eSBarry Smith #endif 3915c6c1daeSBarry Smith /* Create or open the file collectively */ 3925c6c1daeSBarry Smith switch (hdf5->btype) { 3935c6c1daeSBarry Smith case FILE_MODE_READ: 3948a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 3958a2871f6SBarry Smith PetscMPIInt rank; 3968a2871f6SBarry Smith PetscBool flg; 3978a2871f6SBarry Smith 3989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 3998a2871f6SBarry Smith if (rank == 0) { 4009566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4015f80ce2aSJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"File %s requested for reading does not exist",hdf5->filename); 4028a2871f6SBarry Smith } 4039566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 4048a2871f6SBarry Smith } 405729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4065c6c1daeSBarry Smith break; 4075c6c1daeSBarry Smith case FILE_MODE_APPEND: 4087e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4097e4fd573SVaclav Hapla { 4107e4fd573SVaclav Hapla PetscBool flg; 4119566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4127e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4137e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4145c6c1daeSBarry Smith break; 4157e4fd573SVaclav Hapla } 4165c6c1daeSBarry Smith case FILE_MODE_WRITE: 417729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4185c6c1daeSBarry Smith break; 4197e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4207e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4215c6c1daeSBarry Smith default: 42298921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4235c6c1daeSBarry Smith } 4245f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->file_id >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 425729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4265c6c1daeSBarry Smith PetscFunctionReturn(0); 4275c6c1daeSBarry Smith } 4285c6c1daeSBarry Smith 429d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 430d1232d7fSToby Isaac { 431d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 432d1232d7fSToby Isaac 433d1232d7fSToby Isaac PetscFunctionBegin; 434d1232d7fSToby Isaac *name = vhdf5->filename; 435d1232d7fSToby Isaac PetscFunctionReturn(0); 436d1232d7fSToby Isaac } 437d1232d7fSToby Isaac 438b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 439b723ab35SVaclav Hapla { 4405dc64a97SVaclav Hapla /* 441b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4425dc64a97SVaclav Hapla */ 443b723ab35SVaclav Hapla 444b723ab35SVaclav Hapla PetscFunctionBegin; 445b723ab35SVaclav Hapla PetscFunctionReturn(0); 446b723ab35SVaclav Hapla } 447b723ab35SVaclav Hapla 44819a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 44919a20e4cSMatthew G. Knepley { 45019a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 45119a20e4cSMatthew G. Knepley 45219a20e4cSMatthew G. Knepley PetscFunctionBegin; 45319a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 45419a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 45519a20e4cSMatthew G. Knepley } 45619a20e4cSMatthew G. Knepley 45719a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 45819a20e4cSMatthew G. Knepley { 45919a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 46019a20e4cSMatthew G. Knepley 46119a20e4cSMatthew G. Knepley PetscFunctionBegin; 46219a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 46319a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 46419a20e4cSMatthew G. Knepley } 46519a20e4cSMatthew G. Knepley 46619a20e4cSMatthew G. Knepley /*@ 46719a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 46819a20e4cSMatthew G. Knepley 46919a20e4cSMatthew G. Knepley Logically Collective on PetscViewer 47019a20e4cSMatthew G. Knepley 47119a20e4cSMatthew G. Knepley Input Parameters: 47219a20e4cSMatthew G. Knepley + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 47319a20e4cSMatthew G. Knepley - flg - if PETSC_TRUE we will assume that timestepping is on 47419a20e4cSMatthew G. Knepley 47519a20e4cSMatthew G. Knepley Options Database: 47619a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 47719a20e4cSMatthew G. Knepley 47819a20e4cSMatthew G. Knepley Notes: 47919a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 48019a20e4cSMatthew G. Knepley 48119a20e4cSMatthew G. Knepley Level: intermediate 48219a20e4cSMatthew G. Knepley 483db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 48419a20e4cSMatthew G. Knepley @*/ 48519a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 48619a20e4cSMatthew G. Knepley { 48719a20e4cSMatthew G. Knepley PetscFunctionBegin; 48819a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 489cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg)); 49019a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 49119a20e4cSMatthew G. Knepley } 49219a20e4cSMatthew G. Knepley 49319a20e4cSMatthew G. Knepley /*@ 49419a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 49519a20e4cSMatthew G. Knepley 49619a20e4cSMatthew G. Knepley Not collective 49719a20e4cSMatthew G. Knepley 49819a20e4cSMatthew G. Knepley Input Parameter: 49919a20e4cSMatthew G. Knepley . viewer - the PetscViewer 50019a20e4cSMatthew G. Knepley 50119a20e4cSMatthew G. Knepley Output Parameter: 50219a20e4cSMatthew G. Knepley . flg - if PETSC_TRUE we will assume that timestepping is on 50319a20e4cSMatthew G. Knepley 50419a20e4cSMatthew G. Knepley Notes: 50519a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 50619a20e4cSMatthew G. Knepley 50719a20e4cSMatthew G. Knepley Level: intermediate 50819a20e4cSMatthew G. Knepley 509db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 51019a20e4cSMatthew G. Knepley @*/ 51119a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 51219a20e4cSMatthew G. Knepley { 51319a20e4cSMatthew G. Knepley PetscFunctionBegin; 51419a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 515cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg)); 51619a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 51719a20e4cSMatthew G. Knepley } 51819a20e4cSMatthew G. Knepley 5198556b5ebSBarry Smith /*MC 5208556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5218556b5ebSBarry Smith 522db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 523db781477SPatrick Sanan `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 524db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 525db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 5268556b5ebSBarry Smith 5271b266c99SBarry Smith Level: beginner 5288556b5ebSBarry Smith M*/ 529d1232d7fSToby Isaac 5308cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 5315c6c1daeSBarry Smith { 5325c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5335c6c1daeSBarry Smith 5345c6c1daeSBarry Smith PetscFunctionBegin; 53599335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 53699335c34SVaclav Hapla { 53799335c34SVaclav Hapla PetscMPIInt size; 5389566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 5395f80ce2aSJacob Faibussowitsch PetscCheck(size <= 1,PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "Cannot use parallel HDF5 viewer since the given HDF5 does not support parallel I/O (H5_HAVE_PARALLEL is unset)"); 54099335c34SVaclav Hapla } 54199335c34SVaclav Hapla #endif 54299335c34SVaclav Hapla 5439566063dSJacob Faibussowitsch PetscCall(PetscNewLog(v,&hdf5)); 5445c6c1daeSBarry Smith 5455c6c1daeSBarry Smith v->data = (void*) hdf5; 5465c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 54782ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 548b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5491b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5506226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5517e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 552908793a3SLisandro Dalcin hdf5->filename = NULL; 5535c6c1daeSBarry Smith hdf5->timestep = -1; 5540298fd71SBarry Smith hdf5->groups = NULL; 5555c6c1daeSBarry Smith 5569c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 5579c5072fbSVaclav Hapla 5589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5)); 5599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5)); 5609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5)); 5619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5)); 5629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5)); 5639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5)); 5649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5)); 5659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5)); 5669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5)); 5679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5)); 5685c6c1daeSBarry Smith PetscFunctionReturn(0); 5695c6c1daeSBarry Smith } 5705c6c1daeSBarry Smith 5715c6c1daeSBarry Smith /*@C 5725c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 5735c6c1daeSBarry Smith 574d083f849SBarry Smith Collective 5755c6c1daeSBarry Smith 5765c6c1daeSBarry Smith Input Parameters: 5775c6c1daeSBarry Smith + comm - MPI communicator 5785c6c1daeSBarry Smith . name - name of file 5795c6c1daeSBarry Smith - type - type of file 5805c6c1daeSBarry Smith 5815c6c1daeSBarry Smith Output Parameter: 5825c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5835c6c1daeSBarry Smith 58482ea9b62SBarry Smith Options Database: 585a2b725a8SWilliam 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 586a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 58782ea9b62SBarry Smith 5885c6c1daeSBarry Smith Level: beginner 5895c6c1daeSBarry Smith 5907e4fd573SVaclav Hapla Notes: 5917e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5927e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5937e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5947e4fd573SVaclav 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] 5957e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5967e4fd573SVaclav Hapla 5977e4fd573SVaclav 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. 5987e4fd573SVaclav Hapla 5995c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 6005c6c1daeSBarry Smith 601db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 602db781477SPatrick Sanan `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 603db781477SPatrick Sanan `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 6045c6c1daeSBarry Smith @*/ 6055c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 6065c6c1daeSBarry Smith { 6075c6c1daeSBarry Smith PetscFunctionBegin; 6089566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, hdf5v)); 6099566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 6109566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 6119566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*hdf5v, name)); 6129566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*hdf5v)); 6135c6c1daeSBarry Smith PetscFunctionReturn(0); 6145c6c1daeSBarry Smith } 6155c6c1daeSBarry Smith 6165c6c1daeSBarry Smith /*@C 6175c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6185c6c1daeSBarry Smith 6195c6c1daeSBarry Smith Not collective 6205c6c1daeSBarry Smith 6215c6c1daeSBarry Smith Input Parameter: 6225c6c1daeSBarry Smith . viewer - the PetscViewer 6235c6c1daeSBarry Smith 6245c6c1daeSBarry Smith Output Parameter: 6255c6c1daeSBarry Smith . file_id - The file id 6265c6c1daeSBarry Smith 6275c6c1daeSBarry Smith Level: intermediate 6285c6c1daeSBarry Smith 629db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()` 6305c6c1daeSBarry Smith @*/ 6315c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 6325c6c1daeSBarry Smith { 6335c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6345c6c1daeSBarry Smith 6355c6c1daeSBarry Smith PetscFunctionBegin; 6365c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6375c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6385c6c1daeSBarry Smith PetscFunctionReturn(0); 6395c6c1daeSBarry Smith } 6405c6c1daeSBarry Smith 6415c6c1daeSBarry Smith /*@C 6425c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith Not collective 6455c6c1daeSBarry Smith 6465c6c1daeSBarry Smith Input Parameters: 6475c6c1daeSBarry Smith + viewer - the PetscViewer 6485c6c1daeSBarry Smith - name - The group name 6495c6c1daeSBarry Smith 6505c6c1daeSBarry Smith Level: intermediate 6515c6c1daeSBarry Smith 6522e361344SVaclav Hapla Notes: 6532e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 6542e361344SVaclav 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. 6552e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 6562e361344SVaclav Hapla - "." means the current group is pushed again. 6572e361344SVaclav Hapla 6582e361344SVaclav Hapla Example: 6592e361344SVaclav Hapla Suppose the current group is "/a". 6602e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6612e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 6622e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 6632e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 6642e361344SVaclav Hapla 6652e361344SVaclav Hapla Developer Notes: 6662e361344SVaclav Hapla The root group "/" is internally stored as NULL. 667820fc2d1SVaclav Hapla 668*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 6695c6c1daeSBarry Smith @*/ 670be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 6715c6c1daeSBarry Smith { 6725c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6737d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6742e361344SVaclav Hapla size_t i,len; 6752e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6762e361344SVaclav Hapla const char *gname; 6775c6c1daeSBarry Smith 6785c6c1daeSBarry Smith PetscFunctionBegin; 6795c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 680820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 6819566063dSJacob Faibussowitsch PetscCall(PetscStrlen(name, &len)); 6822e361344SVaclav Hapla gname = NULL; 6832e361344SVaclav Hapla if (len) { 6842e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6852e361344SVaclav Hapla /* use current name */ 6862e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6872e361344SVaclav Hapla } else if (name[0] == '/') { 6882e361344SVaclav Hapla /* absolute */ 6892e361344SVaclav Hapla for (i=1; i<len; i++) { 6902e361344SVaclav Hapla if (name[i] != '/') { 6912e361344SVaclav Hapla gname = name; 6922e361344SVaclav Hapla break; 6932e361344SVaclav Hapla } 6942e361344SVaclav Hapla } 6952e361344SVaclav Hapla } else { 6962e361344SVaclav Hapla /* relative */ 6972e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 6989566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 6992e361344SVaclav Hapla gname = buf; 7002e361344SVaclav Hapla } 7012e361344SVaclav Hapla } 7029566063dSJacob Faibussowitsch PetscCall(PetscNew(&groupNode)); 7039566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(gname, (char**) &groupNode->name)); 7045c6c1daeSBarry Smith groupNode->next = hdf5->groups; 7055c6c1daeSBarry Smith hdf5->groups = groupNode; 7065c6c1daeSBarry Smith PetscFunctionReturn(0); 7075c6c1daeSBarry Smith } 7085c6c1daeSBarry Smith 7093ef9c667SSatish Balay /*@ 7105c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 7115c6c1daeSBarry Smith 7125c6c1daeSBarry Smith Not collective 7135c6c1daeSBarry Smith 7145c6c1daeSBarry Smith Input Parameter: 7155c6c1daeSBarry Smith . viewer - the PetscViewer 7165c6c1daeSBarry Smith 7175c6c1daeSBarry Smith Level: intermediate 7185c6c1daeSBarry Smith 719*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 7205c6c1daeSBarry Smith @*/ 7215c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 7225c6c1daeSBarry Smith { 7235c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7247d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7255c6c1daeSBarry Smith 7265c6c1daeSBarry Smith PetscFunctionBegin; 7275c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7285f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->groups,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7295c6c1daeSBarry Smith groupNode = hdf5->groups; 7305c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7319566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode->name)); 7329566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode)); 7335c6c1daeSBarry Smith PetscFunctionReturn(0); 7345c6c1daeSBarry Smith } 7355c6c1daeSBarry Smith 7365c6c1daeSBarry Smith /*@C 737874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 738874270d9SVaclav Hapla If none has been assigned, returns NULL. 7395c6c1daeSBarry Smith 7405c6c1daeSBarry Smith Not collective 7415c6c1daeSBarry Smith 7425c6c1daeSBarry Smith Input Parameter: 7435c6c1daeSBarry Smith . viewer - the PetscViewer 7445c6c1daeSBarry Smith 7455c6c1daeSBarry Smith Output Parameter: 7465c6c1daeSBarry Smith . name - The group name 7475c6c1daeSBarry Smith 7485c6c1daeSBarry Smith Level: intermediate 7495c6c1daeSBarry Smith 750*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()` 7515c6c1daeSBarry Smith @*/ 752be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 7535c6c1daeSBarry Smith { 7545c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 7555c6c1daeSBarry Smith 7565c6c1daeSBarry Smith PetscFunctionBegin; 7575c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 758c959eef4SJed Brown PetscValidPointer(name,2); 759a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 7600298fd71SBarry Smith else *name = NULL; 7615c6c1daeSBarry Smith PetscFunctionReturn(0); 7625c6c1daeSBarry Smith } 7635c6c1daeSBarry Smith 7648c557b5aSVaclav Hapla /*@ 765874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 766874270d9SVaclav Hapla and return this group's ID and file ID. 767874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 768874270d9SVaclav Hapla 769874270d9SVaclav Hapla Not collective 770874270d9SVaclav Hapla 771874270d9SVaclav Hapla Input Parameter: 772874270d9SVaclav Hapla . viewer - the PetscViewer 773874270d9SVaclav Hapla 774d8d19677SJose E. Roman Output Parameters: 775874270d9SVaclav Hapla + fileId - The HDF5 file ID 776874270d9SVaclav Hapla - groupId - The HDF5 group ID 777874270d9SVaclav Hapla 778e74616beSVaclav Hapla Notes: 779e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 780e74616beSVaclav Hapla 781874270d9SVaclav Hapla Level: intermediate 782874270d9SVaclav Hapla 783*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 784874270d9SVaclav Hapla @*/ 78554dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 78654dbf706SVaclav Hapla { 78790e03537SVaclav Hapla hid_t file_id; 78890e03537SVaclav Hapla H5O_type_t type; 78976d59af2SVaclav Hapla const char *groupName = NULL, *fileName = NULL; 79076d59af2SVaclav Hapla PetscBool writable, has; 79154dbf706SVaclav Hapla 79254dbf706SVaclav Hapla PetscFunctionBegin; 7939566063dSJacob Faibussowitsch PetscCall(PetscViewerWritable(viewer, &writable)); 7949566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 7959566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetName(viewer, &fileName)); 7969566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName)); 7979566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 79876d59af2SVaclav Hapla if (!has) { 7995f80ce2aSJacob Faibussowitsch PetscCheck(writable,PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Group %s does not exist and file %s is not open for writing", groupName, fileName); 80098921bdaSJacob Faibussowitsch else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 80176d59af2SVaclav Hapla } 8025f80ce2aSJacob Faibussowitsch PetscCheck(type == H5O_TYPE_GROUP,PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Path %s in file %s resolves to something which is not a group", groupName, fileName); 80390e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 80454dbf706SVaclav Hapla *fileId = file_id; 80554dbf706SVaclav Hapla PetscFunctionReturn(0); 80654dbf706SVaclav Hapla } 80754dbf706SVaclav Hapla 8083ef9c667SSatish Balay /*@ 809d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8105c6c1daeSBarry Smith 8115c6c1daeSBarry Smith Not collective 8125c6c1daeSBarry Smith 8135c6c1daeSBarry Smith Input Parameter: 8145c6c1daeSBarry Smith . viewer - the PetscViewer 8155c6c1daeSBarry Smith 8165c6c1daeSBarry Smith Level: intermediate 8175c6c1daeSBarry Smith 818d7dd068bSVaclav Hapla Notes: 819d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 820d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 821d7dd068bSVaclav 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()]. 822d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 823d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 824d7dd068bSVaclav Hapla 825d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 826d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 827d7dd068bSVaclav Hapla 828d7dd068bSVaclav Hapla Developer notes: 829d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 830d7dd068bSVaclav Hapla 831db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 832d7dd068bSVaclav Hapla @*/ 833d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 834d7dd068bSVaclav Hapla { 835d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 836d7dd068bSVaclav Hapla 837d7dd068bSVaclav Hapla PetscFunctionBegin; 838d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 8395f80ce2aSJacob Faibussowitsch PetscCheck(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 840d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 841d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 842d7dd068bSVaclav Hapla PetscFunctionReturn(0); 843d7dd068bSVaclav Hapla } 844d7dd068bSVaclav Hapla 845d7dd068bSVaclav Hapla /*@ 846d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 847d7dd068bSVaclav Hapla 848d7dd068bSVaclav Hapla Not collective 849d7dd068bSVaclav Hapla 850d7dd068bSVaclav Hapla Input Parameter: 851d7dd068bSVaclav Hapla . viewer - the PetscViewer 852d7dd068bSVaclav Hapla 853d7dd068bSVaclav Hapla Level: intermediate 854d7dd068bSVaclav Hapla 855d7dd068bSVaclav Hapla Notes: 856d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 857d7dd068bSVaclav Hapla 858db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 859d7dd068bSVaclav Hapla @*/ 860d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 861d7dd068bSVaclav Hapla { 862d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 863d7dd068bSVaclav Hapla 864d7dd068bSVaclav Hapla PetscFunctionBegin; 865d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 8665f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 867d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 868d7dd068bSVaclav Hapla PetscFunctionReturn(0); 869d7dd068bSVaclav Hapla } 870d7dd068bSVaclav Hapla 871d7dd068bSVaclav Hapla /*@ 872d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 873d7dd068bSVaclav Hapla 874d7dd068bSVaclav Hapla Not collective 875d7dd068bSVaclav Hapla 876d7dd068bSVaclav Hapla Input Parameter: 877d7dd068bSVaclav Hapla . viewer - the PetscViewer 878d7dd068bSVaclav Hapla 879d7dd068bSVaclav Hapla Output Parameter: 880d7dd068bSVaclav Hapla . flg - is timestepping active? 881d7dd068bSVaclav Hapla 882d7dd068bSVaclav Hapla Level: intermediate 883d7dd068bSVaclav Hapla 884d7dd068bSVaclav Hapla Notes: 885d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 886d7dd068bSVaclav Hapla 887db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 888d7dd068bSVaclav Hapla @*/ 889d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 890d7dd068bSVaclav Hapla { 891d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 892d7dd068bSVaclav Hapla 893d7dd068bSVaclav Hapla PetscFunctionBegin; 894d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 895d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 896d7dd068bSVaclav Hapla PetscFunctionReturn(0); 897d7dd068bSVaclav Hapla } 898d7dd068bSVaclav Hapla 899d7dd068bSVaclav Hapla /*@ 900d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 901d7dd068bSVaclav Hapla 902d7dd068bSVaclav Hapla Not collective 903d7dd068bSVaclav Hapla 904d7dd068bSVaclav Hapla Input Parameter: 905d7dd068bSVaclav Hapla . viewer - the PetscViewer 906d7dd068bSVaclav Hapla 907d7dd068bSVaclav Hapla Level: intermediate 908d7dd068bSVaclav Hapla 909d7dd068bSVaclav Hapla Notes: 910d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 911d7dd068bSVaclav Hapla 912db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 9135c6c1daeSBarry Smith @*/ 9145c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 9155c6c1daeSBarry Smith { 9165c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9175c6c1daeSBarry Smith 9185c6c1daeSBarry Smith PetscFunctionBegin; 9195c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9205f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9215c6c1daeSBarry Smith ++hdf5->timestep; 9225c6c1daeSBarry Smith PetscFunctionReturn(0); 9235c6c1daeSBarry Smith } 9245c6c1daeSBarry Smith 9253ef9c667SSatish Balay /*@ 926d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9275c6c1daeSBarry Smith 928d7dd068bSVaclav Hapla Logically collective 9295c6c1daeSBarry Smith 9305c6c1daeSBarry Smith Input Parameters: 9315c6c1daeSBarry Smith + viewer - the PetscViewer 932d7dd068bSVaclav Hapla - timestep - The timestep 9335c6c1daeSBarry Smith 9345c6c1daeSBarry Smith Level: intermediate 9355c6c1daeSBarry Smith 936d7dd068bSVaclav Hapla Notes: 937d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 938d7dd068bSVaclav Hapla 939db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 9405c6c1daeSBarry Smith @*/ 9415c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 9425c6c1daeSBarry Smith { 9435c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9445c6c1daeSBarry Smith 9455c6c1daeSBarry Smith PetscFunctionBegin; 9465c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 947d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 9485f80ce2aSJacob Faibussowitsch PetscCheck(timestep >= 0,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 9495f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9505c6c1daeSBarry Smith hdf5->timestep = timestep; 9515c6c1daeSBarry Smith PetscFunctionReturn(0); 9525c6c1daeSBarry Smith } 9535c6c1daeSBarry Smith 9543ef9c667SSatish Balay /*@ 9555c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 9565c6c1daeSBarry Smith 9575c6c1daeSBarry Smith Not collective 9585c6c1daeSBarry Smith 9595c6c1daeSBarry Smith Input Parameter: 9605c6c1daeSBarry Smith . viewer - the PetscViewer 9615c6c1daeSBarry Smith 9625c6c1daeSBarry Smith Output Parameter: 963d7dd068bSVaclav Hapla . timestep - The timestep 9645c6c1daeSBarry Smith 9655c6c1daeSBarry Smith Level: intermediate 9665c6c1daeSBarry Smith 967d7dd068bSVaclav Hapla Notes: 968d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 969d7dd068bSVaclav Hapla 970db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 9715c6c1daeSBarry Smith @*/ 9725c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 9735c6c1daeSBarry Smith { 9745c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9755c6c1daeSBarry Smith 9765c6c1daeSBarry Smith PetscFunctionBegin; 9775c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 978d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 9795f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9805c6c1daeSBarry Smith *timestep = hdf5->timestep; 9815c6c1daeSBarry Smith PetscFunctionReturn(0); 9825c6c1daeSBarry Smith } 9835c6c1daeSBarry Smith 98436bce6e8SMatthew G. Knepley /*@C 98536bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 98636bce6e8SMatthew G. Knepley 98736bce6e8SMatthew G. Knepley Not collective 98836bce6e8SMatthew G. Knepley 98936bce6e8SMatthew G. Knepley Input Parameter: 99036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 99136bce6e8SMatthew G. Knepley 99236bce6e8SMatthew G. Knepley Output Parameter: 99336bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 99436bce6e8SMatthew G. Knepley 99536bce6e8SMatthew G. Knepley Level: advanced 99636bce6e8SMatthew G. Knepley 997db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 99836bce6e8SMatthew G. Knepley @*/ 99936bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 100036bce6e8SMatthew G. Knepley { 100136bce6e8SMatthew G. Knepley PetscFunctionBegin; 100236bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 100336bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 100436bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 100536bce6e8SMatthew G. Knepley #else 100636bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 100736bce6e8SMatthew G. Knepley #endif 100836bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 100936bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 101036bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1011c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1012de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 101336bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 101436bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 101536bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10167e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 101736bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 101836bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 101936bce6e8SMatthew G. Knepley } 102036bce6e8SMatthew G. Knepley 102136bce6e8SMatthew G. Knepley /*@C 102236bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 102336bce6e8SMatthew G. Knepley 102436bce6e8SMatthew G. Knepley Not collective 102536bce6e8SMatthew G. Knepley 102636bce6e8SMatthew G. Knepley Input Parameter: 102736bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 102836bce6e8SMatthew G. Knepley 102936bce6e8SMatthew G. Knepley Output Parameter: 103036bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 103136bce6e8SMatthew G. Knepley 103236bce6e8SMatthew G. Knepley Level: advanced 103336bce6e8SMatthew G. Knepley 1034db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 103536bce6e8SMatthew G. Knepley @*/ 103636bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 103736bce6e8SMatthew G. Knepley { 103836bce6e8SMatthew G. Knepley PetscFunctionBegin; 103936bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 104036bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 104136bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 104236bce6e8SMatthew G. Knepley #else 104336bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 104436bce6e8SMatthew G. Knepley #endif 104536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 104636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 104736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 104836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 104936bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 105036bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 10517e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 105236bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 105336bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 105436bce6e8SMatthew G. Knepley } 105536bce6e8SMatthew G. Knepley 1056df863907SAlex Fikl /*@C 1057b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 105836bce6e8SMatthew G. Knepley 1059343a20b2SVaclav Hapla Collective 1060343a20b2SVaclav Hapla 106136bce6e8SMatthew G. Knepley Input Parameters: 106236bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 10634302210dSVaclav Hapla . parent - The parent dataset/group name 106436bce6e8SMatthew G. Knepley . name - The attribute name 106536bce6e8SMatthew G. Knepley . datatype - The attribute type 106636bce6e8SMatthew G. Knepley - value - The attribute value 106736bce6e8SMatthew G. Knepley 106836bce6e8SMatthew G. Knepley Level: advanced 106936bce6e8SMatthew G. Knepley 10704302210dSVaclav Hapla Notes: 1071343a20b2SVaclav 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. 10724302210dSVaclav Hapla 1073*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 107436bce6e8SMatthew G. Knepley @*/ 10754302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 107636bce6e8SMatthew G. Knepley { 10774302210dSVaclav Hapla char *parentAbsPath; 107860568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1079080f144cSVaclav Hapla PetscBool has; 108036bce6e8SMatthew G. Knepley 108136bce6e8SMatthew G. Knepley PetscFunctionBegin; 10825cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 10834302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1084c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 10854302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1086b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 10879566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 10889566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 10899566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 10909566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 10917e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 10927e97c476SMatthew G. Knepley size_t len; 10939566063dSJacob Faibussowitsch PetscCall(PetscStrlen((const char *) value, &len)); 1094729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 10957e97c476SMatthew G. Knepley } 10969566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1097729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 10984302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1099080f144cSVaclav Hapla if (has) { 1100080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1101080f144cSVaclav Hapla } else { 110260568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1103080f144cSVaclav Hapla } 1104729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1105729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1106729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 110760568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1108729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 11099566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 111036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 111136bce6e8SMatthew G. Knepley } 111236bce6e8SMatthew G. Knepley 1113df863907SAlex Fikl /*@C 1114577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1115577e0e04SVaclav Hapla 1116343a20b2SVaclav Hapla Collective 1117343a20b2SVaclav Hapla 1118577e0e04SVaclav Hapla Input Parameters: 1119577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1120577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1121577e0e04SVaclav Hapla . name - The attribute name 1122577e0e04SVaclav Hapla . datatype - The attribute type 1123577e0e04SVaclav Hapla - value - The attribute value 1124577e0e04SVaclav Hapla 1125577e0e04SVaclav Hapla Notes: 1126343a20b2SVaclav 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). 1127577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1128577e0e04SVaclav Hapla 1129577e0e04SVaclav Hapla Level: advanced 1130577e0e04SVaclav Hapla 1131*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1132577e0e04SVaclav Hapla @*/ 1133577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1134577e0e04SVaclav Hapla { 1135577e0e04SVaclav Hapla PetscFunctionBegin; 1136577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1137577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1138577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1139b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 11409566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 11419566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 1142577e0e04SVaclav Hapla PetscFunctionReturn(0); 1143577e0e04SVaclav Hapla } 1144577e0e04SVaclav Hapla 1145577e0e04SVaclav Hapla /*@C 1146b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1147ad0c61feSMatthew G. Knepley 1148343a20b2SVaclav Hapla Collective 1149343a20b2SVaclav Hapla 1150ad0c61feSMatthew G. Knepley Input Parameters: 1151ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 11524302210dSVaclav Hapla . parent - The parent dataset/group name 1153ad0c61feSMatthew G. Knepley . name - The attribute name 1154a2d6be1bSVaclav Hapla . datatype - The attribute type 1155a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1156ad0c61feSMatthew G. Knepley 1157ad0c61feSMatthew G. Knepley Output Parameter: 1158a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1159ad0c61feSMatthew G. Knepley 1160a2d6be1bSVaclav Hapla Notes: 1161a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1162a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1163a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1164a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 11659566063dSJacob Faibussowitsch $ PetscCall(PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 1166a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1167a2d6be1bSVaclav Hapla 1168a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1169708d4cb3SBarry Smith 1170343a20b2SVaclav 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. 1171ad0c61feSMatthew G. Knepley 1172343a20b2SVaclav Hapla Level: advanced 11734302210dSVaclav Hapla 1174*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1175ad0c61feSMatthew G. Knepley @*/ 1176d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1177ad0c61feSMatthew G. Knepley { 11784302210dSVaclav Hapla char *parentAbsPath; 1179a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1180969748fdSVaclav Hapla PetscBool has; 1181ad0c61feSMatthew G. Knepley 1182ad0c61feSMatthew G. Knepley PetscFunctionBegin; 11835cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11844302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1185c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1186a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1187a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 11889566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 11899566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 11909566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 11919566063dSJacob Faibussowitsch if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1192a2d6be1bSVaclav Hapla if (!has) { 1193a2d6be1bSVaclav Hapla if (defaultValue) { 1194a2d6be1bSVaclav Hapla if (defaultValue != value) { 1195a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 11969566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(*(char**)defaultValue, (char**)value)); 1197a2d6be1bSVaclav Hapla } else { 1198a2d6be1bSVaclav Hapla size_t len; 1199ece88022SPierre Jolivet PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(dtype)); 12009566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(value, defaultValue, len)); 1201a2d6be1bSVaclav Hapla } 1202a2d6be1bSVaclav Hapla } 12039566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1204a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 120598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1206a2d6be1bSVaclav Hapla } 12079566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 12084302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 120960568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1210f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1211f0b58479SMatthew G. Knepley size_t len; 1212a2d6be1bSVaclav Hapla hid_t atype; 1213e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1214ece88022SPierre Jolivet PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(atype)); 12159566063dSJacob Faibussowitsch PetscCall(PetscMalloc((len+1) * sizeof(char), value)); 1216708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1217708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1218708d4cb3SBarry Smith } else { 121970efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1220708d4cb3SBarry Smith } 1221729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1222e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1223e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 12249566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1225ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1226ad0c61feSMatthew G. Knepley } 1227ad0c61feSMatthew G. Knepley 1228577e0e04SVaclav Hapla /*@C 1229577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1230577e0e04SVaclav Hapla 1231343a20b2SVaclav Hapla Collective 1232343a20b2SVaclav Hapla 1233577e0e04SVaclav Hapla Input Parameters: 1234577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1235577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1236577e0e04SVaclav Hapla . name - The attribute name 1237577e0e04SVaclav Hapla - datatype - The attribute type 1238577e0e04SVaclav Hapla 1239577e0e04SVaclav Hapla Output Parameter: 1240577e0e04SVaclav Hapla . value - The attribute value 1241577e0e04SVaclav Hapla 1242577e0e04SVaclav Hapla Notes: 1243577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1244577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1245577e0e04SVaclav Hapla 1246577e0e04SVaclav Hapla Level: advanced 1247577e0e04SVaclav Hapla 1248*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1249577e0e04SVaclav Hapla @*/ 1250a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1251577e0e04SVaclav Hapla { 1252577e0e04SVaclav Hapla PetscFunctionBegin; 1253577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1254577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1255577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1256064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 12579566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 12589566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 1259577e0e04SVaclav Hapla PetscFunctionReturn(0); 1260577e0e04SVaclav Hapla } 1261577e0e04SVaclav Hapla 12629fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1263a75c3fd4SVaclav Hapla { 1264a75c3fd4SVaclav Hapla htri_t exists; 1265a75c3fd4SVaclav Hapla hid_t group; 1266a75c3fd4SVaclav Hapla 1267a75c3fd4SVaclav Hapla PetscFunctionBegin; 1268a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1269a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1270a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1271a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1272a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1273a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1274a75c3fd4SVaclav Hapla } 1275a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1276a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1277a75c3fd4SVaclav Hapla } 1278a75c3fd4SVaclav Hapla 1279bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 12805cdeb108SMatthew G. Knepley { 128190e03537SVaclav Hapla const char rootGroupName[] = "/"; 12825cdeb108SMatthew G. Knepley hid_t h5; 1283e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 128415b861d2SVaclav Hapla PetscInt i; 128515b861d2SVaclav Hapla int n; 128685688ae2SVaclav Hapla char **hierarchy; 128785688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 12885cdeb108SMatthew G. Knepley 12895cdeb108SMatthew G. Knepley PetscFunctionBegin; 12905cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 129190e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 129290e03537SVaclav Hapla else name = rootGroupName; 1293ccd66a83SVaclav Hapla if (has) { 1294064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 12955cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1296ccd66a83SVaclav Hapla } 1297ccd66a83SVaclav Hapla if (otype) { 1298064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 129956cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1300ccd66a83SVaclav Hapla } 13019566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 130285688ae2SVaclav Hapla 130385688ae2SVaclav Hapla /* 130485688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 130585688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 130685688ae2SVaclav Hapla 1) whether it's a valid link 130785688ae2SVaclav Hapla 2) whether this link resolves to an object 130885688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 130985688ae2SVaclav Hapla */ 13109566063dSJacob Faibussowitsch PetscCall(PetscStrToArray(name,'/',&n,&hierarchy)); 131185688ae2SVaclav Hapla if (!n) { 131285688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1313ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1314ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 13159566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n,hierarchy)); 131685688ae2SVaclav Hapla PetscFunctionReturn(0); 131785688ae2SVaclav Hapla } 131885688ae2SVaclav Hapla for (i=0; i<n; i++) { 13199566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf,"/")); 13209566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf,hierarchy[i])); 13219566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1322a75c3fd4SVaclav Hapla if (!exists) break; 132385688ae2SVaclav Hapla } 13249566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n,hierarchy)); 132585688ae2SVaclav Hapla 132685688ae2SVaclav Hapla /* If the object exists, get its type */ 1327ccd66a83SVaclav Hapla if (exists && otype) { 13285cdeb108SMatthew G. Knepley H5O_info_t info; 13295cdeb108SMatthew G. Knepley 133076276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 133104633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 133256cc0592SVaclav Hapla *otype = info.type; 13335cdeb108SMatthew G. Knepley } 1334ccd66a83SVaclav Hapla if (has) *has = exists; 13355cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 13365cdeb108SMatthew G. Knepley } 13375cdeb108SMatthew G. Knepley 13384302210dSVaclav Hapla /*@C 13390a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13400a9f38efSVaclav Hapla 1341343a20b2SVaclav Hapla Collective 1342343a20b2SVaclav Hapla 13430a9f38efSVaclav Hapla Input Parameters: 13444302210dSVaclav Hapla + viewer - The HDF5 viewer 13454302210dSVaclav Hapla - path - The group path 13460a9f38efSVaclav Hapla 13470a9f38efSVaclav Hapla Output Parameter: 13480a9f38efSVaclav Hapla . has - Flag for group existence 13490a9f38efSVaclav Hapla 13500a9f38efSVaclav Hapla Level: advanced 13510a9f38efSVaclav Hapla 13524302210dSVaclav Hapla Notes: 1353785c443eSVaclav 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. 1354785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1355343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 13564302210dSVaclav Hapla 1357db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 13580a9f38efSVaclav Hapla @*/ 13594302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 13600a9f38efSVaclav Hapla { 13610a9f38efSVaclav Hapla H5O_type_t type; 13624302210dSVaclav Hapla char *abspath; 13630a9f38efSVaclav Hapla 13640a9f38efSVaclav Hapla PetscFunctionBegin; 13650a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 13664302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 13674302210dSVaclav Hapla PetscValidBoolPointer(has,3); 13689566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath)); 13699566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 13704302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 13719566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 13720a9f38efSVaclav Hapla PetscFunctionReturn(0); 13730a9f38efSVaclav Hapla } 13740a9f38efSVaclav Hapla 137589e0ef10SVaclav Hapla /*@C 137689e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 137789e0ef10SVaclav Hapla 1378343a20b2SVaclav Hapla Collective 1379343a20b2SVaclav Hapla 138089e0ef10SVaclav Hapla Input Parameters: 138189e0ef10SVaclav Hapla + viewer - The HDF5 viewer 138289e0ef10SVaclav Hapla - path - The dataset path 138389e0ef10SVaclav Hapla 138489e0ef10SVaclav Hapla Output Parameter: 138589e0ef10SVaclav Hapla . has - Flag whether dataset exists 138689e0ef10SVaclav Hapla 138789e0ef10SVaclav Hapla Level: advanced 138889e0ef10SVaclav Hapla 138989e0ef10SVaclav Hapla Notes: 1390343a20b2SVaclav 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. 139189e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1392343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 139389e0ef10SVaclav Hapla 1394db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 139589e0ef10SVaclav Hapla @*/ 139689e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 139789e0ef10SVaclav Hapla { 139889e0ef10SVaclav Hapla H5O_type_t type; 139989e0ef10SVaclav Hapla char *abspath; 140089e0ef10SVaclav Hapla 140189e0ef10SVaclav Hapla PetscFunctionBegin; 140289e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 140389e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 140489e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 14059566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath)); 14069566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 140789e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 14089566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 140989e0ef10SVaclav Hapla PetscFunctionReturn(0); 141089e0ef10SVaclav Hapla } 141189e0ef10SVaclav Hapla 14120a9f38efSVaclav Hapla /*@ 1413e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1414ecce1506SVaclav Hapla 1415343a20b2SVaclav Hapla Collective 1416343a20b2SVaclav Hapla 1417ecce1506SVaclav Hapla Input Parameters: 1418ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1419ecce1506SVaclav Hapla - obj - The named object 1420ecce1506SVaclav Hapla 1421ecce1506SVaclav Hapla Output Parameter: 142289e0ef10SVaclav Hapla . has - Flag for dataset existence 1423ecce1506SVaclav Hapla 1424e3f143f7SVaclav Hapla Notes: 142589e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1426343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1427e3f143f7SVaclav Hapla 1428ecce1506SVaclav Hapla Level: advanced 1429ecce1506SVaclav Hapla 1430db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1431ecce1506SVaclav Hapla @*/ 1432ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1433ecce1506SVaclav Hapla { 143489e0ef10SVaclav Hapla size_t len; 1435ecce1506SVaclav Hapla 1436ecce1506SVaclav Hapla PetscFunctionBegin; 1437c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1438c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 14394302210dSVaclav Hapla PetscValidBoolPointer(has,3); 14409566063dSJacob Faibussowitsch PetscCall(PetscStrlen(obj->name, &len)); 14415f80ce2aSJacob Faibussowitsch PetscCheck(len,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 14429566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 1443ecce1506SVaclav Hapla PetscFunctionReturn(0); 1444ecce1506SVaclav Hapla } 1445ecce1506SVaclav Hapla 1446df863907SAlex Fikl /*@C 1447ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1448e7bdbf8eSMatthew G. Knepley 1449343a20b2SVaclav Hapla Collective 1450343a20b2SVaclav Hapla 1451e7bdbf8eSMatthew G. Knepley Input Parameters: 1452e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 14534302210dSVaclav Hapla . parent - The parent dataset/group name 1454e7bdbf8eSMatthew G. Knepley - name - The attribute name 1455e7bdbf8eSMatthew G. Knepley 1456e7bdbf8eSMatthew G. Knepley Output Parameter: 1457e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1458e7bdbf8eSMatthew G. Knepley 1459e7bdbf8eSMatthew G. Knepley Level: advanced 1460e7bdbf8eSMatthew G. Knepley 14614302210dSVaclav Hapla Notes: 1462343a20b2SVaclav 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. 14634302210dSVaclav Hapla 1464*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1465e7bdbf8eSMatthew G. Knepley @*/ 14664302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1467e7bdbf8eSMatthew G. Knepley { 14684302210dSVaclav Hapla char *parentAbsPath; 1469e7bdbf8eSMatthew G. Knepley 1470e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 14715cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14724302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1473c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 14744302210dSVaclav Hapla PetscValidBoolPointer(has,4); 14759566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 14769566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 14779566063dSJacob Faibussowitsch if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 14789566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 147906db490cSVaclav Hapla PetscFunctionReturn(0); 148006db490cSVaclav Hapla } 148106db490cSVaclav Hapla 1482577e0e04SVaclav Hapla /*@C 1483577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1484577e0e04SVaclav Hapla 1485343a20b2SVaclav Hapla Collective 1486343a20b2SVaclav Hapla 1487577e0e04SVaclav Hapla Input Parameters: 1488577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1489577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1490577e0e04SVaclav Hapla - name - The attribute name 1491577e0e04SVaclav Hapla 1492577e0e04SVaclav Hapla Output Parameter: 1493577e0e04SVaclav Hapla . has - Flag for attribute existence 1494577e0e04SVaclav Hapla 1495577e0e04SVaclav Hapla Notes: 1496577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1497577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1498577e0e04SVaclav Hapla 1499577e0e04SVaclav Hapla Level: advanced 1500577e0e04SVaclav Hapla 1501*c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1502577e0e04SVaclav Hapla @*/ 1503577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1504577e0e04SVaclav Hapla { 1505577e0e04SVaclav Hapla PetscFunctionBegin; 1506577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1507577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1508577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 15094302210dSVaclav Hapla PetscValidBoolPointer(has,4); 15109566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 15119566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 1512577e0e04SVaclav Hapla PetscFunctionReturn(0); 1513577e0e04SVaclav Hapla } 1514577e0e04SVaclav Hapla 151506db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 151606db490cSVaclav Hapla { 151706db490cSVaclav Hapla hid_t h5; 151806db490cSVaclav Hapla htri_t hhas; 151906db490cSVaclav Hapla 152006db490cSVaclav Hapla PetscFunctionBegin; 15219566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 15222f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 15235cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1524e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1525e7bdbf8eSMatthew G. Knepley } 1526e7bdbf8eSMatthew G. Knepley 1527a75e6a4aSMatthew G. Knepley /* 1528a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1529a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1530a75e6a4aSMatthew G. Knepley */ 1531d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1532a75e6a4aSMatthew G. Knepley 1533a75e6a4aSMatthew G. Knepley /*@C 1534a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1535a75e6a4aSMatthew G. Knepley 1536d083f849SBarry Smith Collective 1537a75e6a4aSMatthew G. Knepley 1538a75e6a4aSMatthew G. Knepley Input Parameter: 1539a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1540a75e6a4aSMatthew G. Knepley 1541a75e6a4aSMatthew G. Knepley Level: intermediate 1542a75e6a4aSMatthew G. Knepley 1543a75e6a4aSMatthew G. Knepley Options Database Keys: 154410699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1545a75e6a4aSMatthew G. Knepley 1546a75e6a4aSMatthew G. Knepley Environmental variables: 154710699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1548a75e6a4aSMatthew G. Knepley 1549a75e6a4aSMatthew G. Knepley Notes: 1550a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1551a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1552a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1553a75e6a4aSMatthew G. Knepley 1554db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1555a75e6a4aSMatthew G. Knepley @*/ 1556a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1557a75e6a4aSMatthew G. Knepley { 1558a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1559a75e6a4aSMatthew G. Knepley PetscBool flg; 1560a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1561a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1562a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1563a75e6a4aSMatthew G. Knepley 1564a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1565908793a3SLisandro 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);} 1566a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1567908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1568908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1569a75e6a4aSMatthew G. Knepley } 157047435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1571908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1572a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1573a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 15742cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1575a75e6a4aSMatthew G. Knepley if (!flg) { 1576a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 15772cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1578a75e6a4aSMatthew G. Knepley } 1579a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 15802cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1581a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 15822cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 158347435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1584908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1585a75e6a4aSMatthew G. Knepley } 1586a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 15872cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1588a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1589a75e6a4aSMatthew G. Knepley } 1590