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)); 1122e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL)); 1139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetBaseDimension2_C",NULL)); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetSPOutput_C",NULL)); 1152e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetCollective_C",NULL)); 1162e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetCollective_C",NULL)); 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5GetDefaultTimestepping_C",NULL)); 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerHDF5SetDefaultTimestepping_C",NULL)); 1195c6c1daeSBarry Smith PetscFunctionReturn(0); 1205c6c1daeSBarry Smith } 1215c6c1daeSBarry Smith 1229b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetMode_HDF5(PetscViewer viewer, PetscFileMode type) 1235c6c1daeSBarry Smith { 1245c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1255c6c1daeSBarry Smith 1265c6c1daeSBarry Smith PetscFunctionBegin; 1275c6c1daeSBarry Smith hdf5->btype = type; 1285c6c1daeSBarry Smith PetscFunctionReturn(0); 1295c6c1daeSBarry Smith } 1305c6c1daeSBarry Smith 1310b62783dSVaclav Hapla static PetscErrorCode PetscViewerFileGetMode_HDF5(PetscViewer viewer, PetscFileMode *type) 1320b62783dSVaclav Hapla { 1330b62783dSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 1340b62783dSVaclav Hapla 1350b62783dSVaclav Hapla PetscFunctionBegin; 1360b62783dSVaclav Hapla *type = hdf5->btype; 1370b62783dSVaclav Hapla PetscFunctionReturn(0); 1380b62783dSVaclav Hapla } 1390b62783dSVaclav Hapla 1409b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetBaseDimension2_HDF5(PetscViewer viewer, PetscBool flg) 14182ea9b62SBarry Smith { 14282ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 14382ea9b62SBarry Smith 14482ea9b62SBarry Smith PetscFunctionBegin; 14582ea9b62SBarry Smith hdf5->basedimension2 = flg; 14682ea9b62SBarry Smith PetscFunctionReturn(0); 14782ea9b62SBarry Smith } 14882ea9b62SBarry Smith 149df863907SAlex Fikl /*@ 15082ea9b62SBarry Smith PetscViewerHDF5SetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 15182ea9b62SBarry Smith dimension of 2. 15282ea9b62SBarry Smith 15382ea9b62SBarry Smith Logically Collective on PetscViewer 15482ea9b62SBarry Smith 15582ea9b62SBarry Smith Input Parameters: 15682ea9b62SBarry Smith + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 15782ea9b62SBarry 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 15882ea9b62SBarry Smith 15982ea9b62SBarry Smith Options Database: 16082ea9b62SBarry 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 16182ea9b62SBarry Smith 16295452b02SPatrick Sanan Notes: 16395452b02SPatrick 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 16482ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 16582ea9b62SBarry Smith 16682ea9b62SBarry Smith Level: intermediate 16782ea9b62SBarry Smith 168db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 16982ea9b62SBarry Smith 17082ea9b62SBarry Smith @*/ 17182ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5SetBaseDimension2(PetscViewer viewer,PetscBool flg) 17282ea9b62SBarry Smith { 17382ea9b62SBarry Smith PetscFunctionBegin; 17482ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 175cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetBaseDimension2_C",(PetscViewer,PetscBool),(viewer,flg)); 17682ea9b62SBarry Smith PetscFunctionReturn(0); 17782ea9b62SBarry Smith } 17882ea9b62SBarry Smith 179df863907SAlex Fikl /*@ 18082ea9b62SBarry Smith PetscViewerHDF5GetBaseDimension2 - Vectors of 1 dimension (i.e. bs/dof is 1) will be saved in the HDF5 file with a 18182ea9b62SBarry Smith dimension of 2. 18282ea9b62SBarry Smith 18382ea9b62SBarry Smith Logically Collective on PetscViewer 18482ea9b62SBarry Smith 18582ea9b62SBarry Smith Input Parameter: 18682ea9b62SBarry Smith . viewer - the PetscViewer, must be of type HDF5 18782ea9b62SBarry Smith 18882ea9b62SBarry Smith Output Parameter: 18982ea9b62SBarry 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 19082ea9b62SBarry Smith 19195452b02SPatrick Sanan Notes: 19295452b02SPatrick 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 19382ea9b62SBarry Smith of one when the dimension is lower. Others think the option is crazy. 19482ea9b62SBarry Smith 19582ea9b62SBarry Smith Level: intermediate 19682ea9b62SBarry Smith 197db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()` 19882ea9b62SBarry Smith 19982ea9b62SBarry Smith @*/ 20082ea9b62SBarry Smith PetscErrorCode PetscViewerHDF5GetBaseDimension2(PetscViewer viewer,PetscBool *flg) 20182ea9b62SBarry Smith { 20282ea9b62SBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 20382ea9b62SBarry Smith 20482ea9b62SBarry Smith PetscFunctionBegin; 20582ea9b62SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 20682ea9b62SBarry Smith *flg = hdf5->basedimension2; 20782ea9b62SBarry Smith PetscFunctionReturn(0); 20882ea9b62SBarry Smith } 20982ea9b62SBarry Smith 2109b1403beSVaclav Hapla static PetscErrorCode PetscViewerHDF5SetSPOutput_HDF5(PetscViewer viewer, PetscBool flg) 2119a0502c6SHåkon Strandenes { 2129a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2139a0502c6SHåkon Strandenes 2149a0502c6SHåkon Strandenes PetscFunctionBegin; 2159a0502c6SHåkon Strandenes hdf5->spoutput = flg; 2169a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2179a0502c6SHåkon Strandenes } 2189a0502c6SHåkon Strandenes 219df863907SAlex Fikl /*@ 2209a0502c6SHåkon Strandenes PetscViewerHDF5SetSPOutput - Data is written to disk in single precision even if PETSc is 2219a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2229a0502c6SHåkon Strandenes 2239a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2249a0502c6SHåkon Strandenes 2259a0502c6SHåkon Strandenes Input Parameters: 2269a0502c6SHåkon Strandenes + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 2279a0502c6SHåkon Strandenes - flg - if PETSC_TRUE the data will be written to disk with single precision 2289a0502c6SHåkon Strandenes 2299a0502c6SHåkon Strandenes Options Database: 2309a0502c6SHåkon Strandenes . -viewer_hdf5_sp_output - turns on (true) or off (false) output in single precision 2319a0502c6SHåkon Strandenes 23295452b02SPatrick Sanan Notes: 23395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2349a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2359a0502c6SHåkon Strandenes 2369a0502c6SHåkon Strandenes Level: intermediate 2379a0502c6SHåkon Strandenes 238db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 239db781477SPatrick Sanan `PetscReal` 2409a0502c6SHåkon Strandenes 2419a0502c6SHåkon Strandenes @*/ 2429a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5SetSPOutput(PetscViewer viewer,PetscBool flg) 2439a0502c6SHåkon Strandenes { 2449a0502c6SHåkon Strandenes PetscFunctionBegin; 2459a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 246cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetSPOutput_C",(PetscViewer,PetscBool),(viewer,flg)); 2479a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2489a0502c6SHåkon Strandenes } 2499a0502c6SHåkon Strandenes 250df863907SAlex Fikl /*@ 2519a0502c6SHåkon Strandenes PetscViewerHDF5GetSPOutput - Data is written to disk in single precision even if PETSc is 2529a0502c6SHåkon Strandenes compiled with double precision PetscReal. 2539a0502c6SHåkon Strandenes 2549a0502c6SHåkon Strandenes Logically Collective on PetscViewer 2559a0502c6SHåkon Strandenes 2569a0502c6SHåkon Strandenes Input Parameter: 2579a0502c6SHåkon Strandenes . viewer - the PetscViewer, must be of type HDF5 2589a0502c6SHåkon Strandenes 2599a0502c6SHåkon Strandenes Output Parameter: 2609a0502c6SHåkon Strandenes . flg - if PETSC_TRUE the data will be written to disk with single precision 2619a0502c6SHåkon Strandenes 26295452b02SPatrick Sanan Notes: 26395452b02SPatrick Sanan Setting this option does not make any difference if PETSc is compiled with single precision 2649a0502c6SHåkon Strandenes in the first place. It does not affect reading datasets (HDF5 handle this internally). 2659a0502c6SHåkon Strandenes 2669a0502c6SHåkon Strandenes Level: intermediate 2679a0502c6SHåkon Strandenes 268db781477SPatrick Sanan .seealso: `PetscViewerFileSetMode()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerBinaryOpen()`, 269db781477SPatrick Sanan `PetscReal` 2709a0502c6SHåkon Strandenes 2719a0502c6SHåkon Strandenes @*/ 2729a0502c6SHåkon Strandenes PetscErrorCode PetscViewerHDF5GetSPOutput(PetscViewer viewer,PetscBool *flg) 2739a0502c6SHåkon Strandenes { 2749a0502c6SHåkon Strandenes PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 2759a0502c6SHåkon Strandenes 2769a0502c6SHåkon Strandenes PetscFunctionBegin; 2779a0502c6SHåkon Strandenes PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 2789a0502c6SHåkon Strandenes *flg = hdf5->spoutput; 2799a0502c6SHåkon Strandenes PetscFunctionReturn(0); 2809a0502c6SHåkon Strandenes } 2819a0502c6SHåkon Strandenes 282ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5SetCollective_HDF5(PetscViewer viewer, PetscBool flg) 283ee8b9732SVaclav Hapla { 284ee8b9732SVaclav Hapla PetscFunctionBegin; 285ee8b9732SVaclav Hapla /* H5FD_MPIO_COLLECTIVE is wrong in hdf5 1.10.2, and is the same as H5FD_MPIO_INDEPENDENT in earlier versions 286ee8b9732SVaclav Hapla - see e.g. https://gitlab.cosma.dur.ac.uk/swift/swiftsim/issues/431 */ 287c796909eSBarry Smith #if H5_VERSION_GE(1,10,3) && defined(H5_HAVE_PARALLEL) 288ee8b9732SVaclav Hapla { 289ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 290ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pset_dxpl_mpio,(hdf5->dxpl_id, flg ? H5FD_MPIO_COLLECTIVE : H5FD_MPIO_INDEPENDENT)); 291ee8b9732SVaclav Hapla } 292ee8b9732SVaclav Hapla #else 2939566063dSJacob 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")); 294ee8b9732SVaclav Hapla #endif 295ee8b9732SVaclav Hapla PetscFunctionReturn(0); 296ee8b9732SVaclav Hapla } 297ee8b9732SVaclav Hapla 298ee8b9732SVaclav Hapla /*@ 299ee8b9732SVaclav Hapla PetscViewerHDF5SetCollective - Use collective MPI-IO transfer mode for HDF5 reads and writes. 300ee8b9732SVaclav Hapla 301ee8b9732SVaclav Hapla Logically Collective; flg must contain common value 302ee8b9732SVaclav Hapla 303ee8b9732SVaclav Hapla Input Parameters: 304ee8b9732SVaclav Hapla + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 305ee8b9732SVaclav Hapla - flg - PETSC_TRUE for collective mode; PETSC_FALSE for independent mode (default) 306ee8b9732SVaclav Hapla 307ee8b9732SVaclav Hapla Options Database: 308ee8b9732SVaclav Hapla . -viewer_hdf5_collective - turns on (true) or off (false) collective transfers 309ee8b9732SVaclav Hapla 310ee8b9732SVaclav Hapla Notes: 311ee8b9732SVaclav Hapla Collective mode gives the MPI-IO layer underneath HDF5 a chance to do some additional collective optimizations and hence can perform better. 31253effed7SBarry 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. 313ee8b9732SVaclav Hapla 314ee8b9732SVaclav Hapla Developer notes: 315ee8b9732SVaclav Hapla In the HDF5 layer, PETSC_TRUE / PETSC_FALSE means H5Pset_dxpl_mpio() is called with H5FD_MPIO_COLLECTIVE / H5FD_MPIO_INDEPENDENT, respectively. 316ee8b9732SVaclav Hapla This in turn means use of MPI_File_{read,write}_all / MPI_File_{read,write} in the MPI-IO layer, respectively. 317ee8b9732SVaclav Hapla See HDF5 documentation and MPI-IO documentation for details. 318ee8b9732SVaclav Hapla 319ee8b9732SVaclav Hapla Level: intermediate 320ee8b9732SVaclav Hapla 321db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 322ee8b9732SVaclav Hapla 323ee8b9732SVaclav Hapla @*/ 324ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5SetCollective(PetscViewer viewer,PetscBool flg) 325ee8b9732SVaclav Hapla { 326ee8b9732SVaclav Hapla PetscFunctionBegin; 327ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 328ee8b9732SVaclav Hapla PetscValidLogicalCollectiveBool(viewer,flg,2); 329cac4c232SBarry Smith PetscTryMethod(viewer,"PetscViewerHDF5SetCollective_C",(PetscViewer,PetscBool),(viewer,flg)); 330ee8b9732SVaclav Hapla PetscFunctionReturn(0); 331ee8b9732SVaclav Hapla } 332ee8b9732SVaclav Hapla 333ee8b9732SVaclav Hapla static PetscErrorCode PetscViewerHDF5GetCollective_HDF5(PetscViewer viewer, PetscBool *flg) 334ee8b9732SVaclav Hapla { 335c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 336ee8b9732SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 337ee8b9732SVaclav Hapla H5FD_mpio_xfer_t mode; 338c796909eSBarry Smith #endif 339ee8b9732SVaclav Hapla 340ee8b9732SVaclav Hapla PetscFunctionBegin; 341c796909eSBarry Smith #if !defined(H5_HAVE_PARALLEL) 342c796909eSBarry Smith *flg = PETSC_FALSE; 343c796909eSBarry Smith #else 344ee8b9732SVaclav Hapla PetscStackCallHDF5(H5Pget_dxpl_mpio,(hdf5->dxpl_id, &mode)); 345ee8b9732SVaclav Hapla *flg = (mode == H5FD_MPIO_COLLECTIVE) ? PETSC_TRUE : PETSC_FALSE; 346c796909eSBarry Smith #endif 347ee8b9732SVaclav Hapla PetscFunctionReturn(0); 348ee8b9732SVaclav Hapla } 349ee8b9732SVaclav Hapla 350ee8b9732SVaclav Hapla /*@ 351ee8b9732SVaclav Hapla PetscViewerHDF5GetCollective - Return flag whether collective MPI-IO transfer mode is used for HDF5 reads and writes. 352ee8b9732SVaclav Hapla 353ee8b9732SVaclav Hapla Not Collective 354ee8b9732SVaclav Hapla 355ee8b9732SVaclav Hapla Input Parameters: 356ee8b9732SVaclav Hapla . viewer - the HDF5 PetscViewer 357ee8b9732SVaclav Hapla 358ee8b9732SVaclav Hapla Output Parameters: 359ee8b9732SVaclav Hapla . flg - the flag 360ee8b9732SVaclav Hapla 361ee8b9732SVaclav Hapla Level: intermediate 362ee8b9732SVaclav Hapla 363ee8b9732SVaclav Hapla Notes: 364c796909eSBarry 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. 365ee8b9732SVaclav Hapla For more details, see PetscViewerHDF5SetCollective(). 366ee8b9732SVaclav Hapla 367db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetCollective()`, `PetscViewerCreate()`, `PetscViewerSetType()`, `PetscViewerHDF5Open()` 368ee8b9732SVaclav Hapla 369ee8b9732SVaclav Hapla @*/ 370ee8b9732SVaclav Hapla PetscErrorCode PetscViewerHDF5GetCollective(PetscViewer viewer,PetscBool *flg) 371ee8b9732SVaclav Hapla { 372ee8b9732SVaclav Hapla PetscFunctionBegin; 373ee8b9732SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 374534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 375ee8b9732SVaclav Hapla 376cac4c232SBarry Smith PetscUseMethod(viewer,"PetscViewerHDF5GetCollective_C",(PetscViewer,PetscBool*),(viewer,flg)); 377ee8b9732SVaclav Hapla PetscFunctionReturn(0); 378ee8b9732SVaclav Hapla } 379ee8b9732SVaclav Hapla 3809b1403beSVaclav Hapla static PetscErrorCode PetscViewerFileSetName_HDF5(PetscViewer viewer, const char name[]) 3815c6c1daeSBarry Smith { 3825c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 3835c6c1daeSBarry Smith hid_t plist_id; 3845c6c1daeSBarry Smith 3855c6c1daeSBarry Smith PetscFunctionBegin; 386f683cc58SToby Isaac if (hdf5->file_id) PetscStackCallHDF5(H5Fclose,(hdf5->file_id)); 3879566063dSJacob Faibussowitsch if (hdf5->filename) PetscCall(PetscFree(hdf5->filename)); 3889566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &hdf5->filename)); 3895c6c1daeSBarry Smith /* Set up file access property list with parallel I/O access */ 390729ab6d9SBarry Smith PetscStackCallHDF5Return(plist_id,H5Pcreate,(H5P_FILE_ACCESS)); 391c796909eSBarry Smith #if defined(H5_HAVE_PARALLEL) 392d28bfa55SVaclav Hapla PetscStackCallHDF5(H5Pset_fapl_mpio,(plist_id, PetscObjectComm((PetscObject)viewer), MPI_INFO_NULL)); 393c796909eSBarry Smith #endif 3945c6c1daeSBarry Smith /* Create or open the file collectively */ 3955c6c1daeSBarry Smith switch (hdf5->btype) { 3965c6c1daeSBarry Smith case FILE_MODE_READ: 3978a2871f6SBarry Smith if (PetscDefined(USE_DEBUG)) { 3988a2871f6SBarry Smith PetscMPIInt rank; 3998a2871f6SBarry Smith PetscBool flg; 4008a2871f6SBarry Smith 4019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank)); 4028a2871f6SBarry Smith if (rank == 0) { 4039566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4045f80ce2aSJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_FILE_OPEN,"File %s requested for reading does not exist",hdf5->filename); 4058a2871f6SBarry Smith } 4069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 4078a2871f6SBarry Smith } 408729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDONLY, plist_id)); 4095c6c1daeSBarry Smith break; 4105c6c1daeSBarry Smith case FILE_MODE_APPEND: 4117e4fd573SVaclav Hapla case FILE_MODE_UPDATE: 4127e4fd573SVaclav Hapla { 4137e4fd573SVaclav Hapla PetscBool flg; 4149566063dSJacob Faibussowitsch PetscCall(PetscTestFile(hdf5->filename, 'r', &flg)); 4157e4fd573SVaclav Hapla if (flg) PetscStackCallHDF5Return(hdf5->file_id,H5Fopen,(name, H5F_ACC_RDWR, plist_id)); 4167e4fd573SVaclav Hapla else PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_EXCL, H5P_DEFAULT, plist_id)); 4175c6c1daeSBarry Smith break; 4187e4fd573SVaclav Hapla } 4195c6c1daeSBarry Smith case FILE_MODE_WRITE: 420729ab6d9SBarry Smith PetscStackCallHDF5Return(hdf5->file_id,H5Fcreate,(name, H5F_ACC_TRUNC, H5P_DEFAULT, plist_id)); 4215c6c1daeSBarry Smith break; 4227e4fd573SVaclav Hapla case FILE_MODE_UNDEFINED: 4237e4fd573SVaclav Hapla SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ORDER, "Must call PetscViewerFileSetMode() before PetscViewerFileSetName()"); 4245c6c1daeSBarry Smith default: 42598921bdaSJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP, "Unsupported file mode %s",PetscFileModes[hdf5->btype]); 4265c6c1daeSBarry Smith } 4275f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->file_id >= 0,PETSC_COMM_SELF,PETSC_ERR_LIB, "H5Fcreate failed for %s", name); 428729ab6d9SBarry Smith PetscStackCallHDF5(H5Pclose,(plist_id)); 4295c6c1daeSBarry Smith PetscFunctionReturn(0); 4305c6c1daeSBarry Smith } 4315c6c1daeSBarry Smith 432d1232d7fSToby Isaac static PetscErrorCode PetscViewerFileGetName_HDF5(PetscViewer viewer,const char **name) 433d1232d7fSToby Isaac { 434d1232d7fSToby Isaac PetscViewer_HDF5 *vhdf5 = (PetscViewer_HDF5*)viewer->data; 435d1232d7fSToby Isaac 436d1232d7fSToby Isaac PetscFunctionBegin; 437d1232d7fSToby Isaac *name = vhdf5->filename; 438d1232d7fSToby Isaac PetscFunctionReturn(0); 439d1232d7fSToby Isaac } 440d1232d7fSToby Isaac 441b723ab35SVaclav Hapla static PetscErrorCode PetscViewerSetUp_HDF5(PetscViewer viewer) 442b723ab35SVaclav Hapla { 4435dc64a97SVaclav Hapla /* 444b723ab35SVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 4455dc64a97SVaclav Hapla */ 446b723ab35SVaclav Hapla 447b723ab35SVaclav Hapla PetscFunctionBegin; 448b723ab35SVaclav Hapla PetscFunctionReturn(0); 449b723ab35SVaclav Hapla } 450b723ab35SVaclav Hapla 45119a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5SetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool flg) 45219a20e4cSMatthew G. Knepley { 45319a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 45419a20e4cSMatthew G. Knepley 45519a20e4cSMatthew G. Knepley PetscFunctionBegin; 45619a20e4cSMatthew G. Knepley hdf5->defTimestepping = flg; 45719a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 45819a20e4cSMatthew G. Knepley } 45919a20e4cSMatthew G. Knepley 46019a20e4cSMatthew G. Knepley static PetscErrorCode PetscViewerHDF5GetDefaultTimestepping_HDF5(PetscViewer viewer, PetscBool *flg) 46119a20e4cSMatthew G. Knepley { 46219a20e4cSMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 46319a20e4cSMatthew G. Knepley 46419a20e4cSMatthew G. Knepley PetscFunctionBegin; 46519a20e4cSMatthew G. Knepley *flg = hdf5->defTimestepping; 46619a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 46719a20e4cSMatthew G. Knepley } 46819a20e4cSMatthew G. Knepley 46919a20e4cSMatthew G. Knepley /*@ 47019a20e4cSMatthew G. Knepley PetscViewerHDF5SetDefaultTimestepping - Set the flag for default timestepping 47119a20e4cSMatthew G. Knepley 47219a20e4cSMatthew G. Knepley Logically Collective on PetscViewer 47319a20e4cSMatthew G. Knepley 47419a20e4cSMatthew G. Knepley Input Parameters: 47519a20e4cSMatthew G. Knepley + viewer - the PetscViewer; if it is not hdf5 then this command is ignored 47619a20e4cSMatthew G. Knepley - flg - if PETSC_TRUE we will assume that timestepping is on 47719a20e4cSMatthew G. Knepley 47819a20e4cSMatthew G. Knepley Options Database: 47919a20e4cSMatthew G. Knepley . -viewer_hdf5_default_timestepping - turns on (true) or off (false) default timestepping 48019a20e4cSMatthew G. Knepley 48119a20e4cSMatthew G. Knepley Notes: 48219a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 48319a20e4cSMatthew G. Knepley 48419a20e4cSMatthew G. Knepley Level: intermediate 48519a20e4cSMatthew G. Knepley 486db781477SPatrick Sanan .seealso: `PetscViewerHDF5GetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 48719a20e4cSMatthew G. Knepley @*/ 48819a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5SetDefaultTimestepping(PetscViewer viewer, PetscBool flg) 48919a20e4cSMatthew G. Knepley { 49019a20e4cSMatthew G. Knepley PetscFunctionBegin; 49119a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 492cac4c232SBarry Smith PetscTryMethod(viewer, "PetscViewerHDF5SetDefaultTimestepping_C", (PetscViewer,PetscBool), (viewer,flg)); 49319a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 49419a20e4cSMatthew G. Knepley } 49519a20e4cSMatthew G. Knepley 49619a20e4cSMatthew G. Knepley /*@ 49719a20e4cSMatthew G. Knepley PetscViewerHDF5GetDefaultTimestepping - Get the flag for default timestepping 49819a20e4cSMatthew G. Knepley 49919a20e4cSMatthew G. Knepley Not collective 50019a20e4cSMatthew G. Knepley 50119a20e4cSMatthew G. Knepley Input Parameter: 50219a20e4cSMatthew G. Knepley . viewer - the PetscViewer 50319a20e4cSMatthew G. Knepley 50419a20e4cSMatthew G. Knepley Output Parameter: 50519a20e4cSMatthew G. Knepley . flg - if PETSC_TRUE we will assume that timestepping is on 50619a20e4cSMatthew G. Knepley 50719a20e4cSMatthew G. Knepley Notes: 50819a20e4cSMatthew G. Knepley If the timestepping attribute is not found for an object, then the default timestepping is used 50919a20e4cSMatthew G. Knepley 51019a20e4cSMatthew G. Knepley Level: intermediate 51119a20e4cSMatthew G. Knepley 512db781477SPatrick Sanan .seealso: `PetscViewerHDF5SetDefaultTimestepping()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5GetTimestep()` 51319a20e4cSMatthew G. Knepley @*/ 51419a20e4cSMatthew G. Knepley PetscErrorCode PetscViewerHDF5GetDefaultTimestepping(PetscViewer viewer, PetscBool *flg) 51519a20e4cSMatthew G. Knepley { 51619a20e4cSMatthew G. Knepley PetscFunctionBegin; 51719a20e4cSMatthew G. Knepley PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 1); 518cac4c232SBarry Smith PetscUseMethod(viewer, "PetscViewerHDF5GetDefaultTimestepping_C", (PetscViewer,PetscBool*), (viewer,flg)); 51919a20e4cSMatthew G. Knepley PetscFunctionReturn(0); 52019a20e4cSMatthew G. Knepley } 52119a20e4cSMatthew G. Knepley 5228556b5ebSBarry Smith /*MC 5238556b5ebSBarry Smith PETSCVIEWERHDF5 - A viewer that writes to an HDF5 file 5248556b5ebSBarry Smith 525db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`, 526db781477SPatrick Sanan `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`, 527db781477SPatrick Sanan `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`, 528db781477SPatrick Sanan `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()` 5298556b5ebSBarry Smith 5301b266c99SBarry Smith Level: beginner 5318556b5ebSBarry Smith M*/ 532d1232d7fSToby Isaac 5338cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PetscViewerCreate_HDF5(PetscViewer v) 5345c6c1daeSBarry Smith { 5355c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5; 5365c6c1daeSBarry Smith 5375c6c1daeSBarry Smith PetscFunctionBegin; 53899335c34SVaclav Hapla #if !defined(H5_HAVE_PARALLEL) 53999335c34SVaclav Hapla { 54099335c34SVaclav Hapla PetscMPIInt size; 5419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size)); 5425f80ce2aSJacob 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)"); 54399335c34SVaclav Hapla } 54499335c34SVaclav Hapla #endif 54599335c34SVaclav Hapla 5469566063dSJacob Faibussowitsch PetscCall(PetscNewLog(v,&hdf5)); 5475c6c1daeSBarry Smith 5485c6c1daeSBarry Smith v->data = (void*) hdf5; 5495c6c1daeSBarry Smith v->ops->destroy = PetscViewerDestroy_HDF5; 55082ea9b62SBarry Smith v->ops->setfromoptions = PetscViewerSetFromOptions_HDF5; 551b723ab35SVaclav Hapla v->ops->setup = PetscViewerSetUp_HDF5; 5521b793a25SVaclav Hapla v->ops->view = PetscViewerView_HDF5; 5536226335aSVaclav Hapla v->ops->flush = PetscViewerFlush_HDF5; 5547e4fd573SVaclav Hapla hdf5->btype = FILE_MODE_UNDEFINED; 555908793a3SLisandro Dalcin hdf5->filename = NULL; 5565c6c1daeSBarry Smith hdf5->timestep = -1; 5570298fd71SBarry Smith hdf5->groups = NULL; 5585c6c1daeSBarry Smith 5599c5072fbSVaclav Hapla PetscStackCallHDF5Return(hdf5->dxpl_id,H5Pcreate,(H5P_DATASET_XFER)); 5609c5072fbSVaclav Hapla 5619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_HDF5)); 5629566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_HDF5)); 5639566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_HDF5)); 5649566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_HDF5)); 5659566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetBaseDimension2_C",PetscViewerHDF5SetBaseDimension2_HDF5)); 5669566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetSPOutput_C",PetscViewerHDF5SetSPOutput_HDF5)); 5679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetCollective_C",PetscViewerHDF5SetCollective_HDF5)); 5689566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetCollective_C",PetscViewerHDF5GetCollective_HDF5)); 5699566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5GetDefaultTimestepping_C",PetscViewerHDF5GetDefaultTimestepping_HDF5)); 5709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)v,"PetscViewerHDF5SetDefaultTimestepping_C",PetscViewerHDF5SetDefaultTimestepping_HDF5)); 5715c6c1daeSBarry Smith PetscFunctionReturn(0); 5725c6c1daeSBarry Smith } 5735c6c1daeSBarry Smith 5745c6c1daeSBarry Smith /*@C 5755c6c1daeSBarry Smith PetscViewerHDF5Open - Opens a file for HDF5 input/output. 5765c6c1daeSBarry Smith 577d083f849SBarry Smith Collective 5785c6c1daeSBarry Smith 5795c6c1daeSBarry Smith Input Parameters: 5805c6c1daeSBarry Smith + comm - MPI communicator 5815c6c1daeSBarry Smith . name - name of file 5825c6c1daeSBarry Smith - type - type of file 5835c6c1daeSBarry Smith 5845c6c1daeSBarry Smith Output Parameter: 5855c6c1daeSBarry Smith . hdf5v - PetscViewer for HDF5 input/output to use with the specified file 5865c6c1daeSBarry Smith 58782ea9b62SBarry Smith Options Database: 588a2b725a8SWilliam 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 589a2b725a8SWilliam Gropp - -viewer_hdf5_sp_output - forces (if true) the viewer to write data in single precision independent on the precision of PetscReal 59082ea9b62SBarry Smith 5915c6c1daeSBarry Smith Level: beginner 5925c6c1daeSBarry Smith 5937e4fd573SVaclav Hapla Notes: 5947e4fd573SVaclav Hapla Reading is always available, regardless of the mode. Available modes are 5957e4fd573SVaclav Hapla + FILE_MODE_READ - open existing HDF5 file for read only access, fail if file does not exist [H5Fopen() with H5F_ACC_RDONLY] 5967e4fd573SVaclav Hapla . FILE_MODE_WRITE - if file exists, fully overwrite it, else create new HDF5 file [H5FcreateH5Fcreate() with H5F_ACC_TRUNC] 5977e4fd573SVaclav 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] 5987e4fd573SVaclav Hapla - FILE_MODE_UPDATE - same as FILE_MODE_APPEND 5997e4fd573SVaclav Hapla 6007e4fd573SVaclav 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. 6017e4fd573SVaclav Hapla 6025c6c1daeSBarry Smith This PetscViewer should be destroyed with PetscViewerDestroy(). 6035c6c1daeSBarry Smith 604db781477SPatrick Sanan .seealso: `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`, `PetscViewerHDF5SetBaseDimension2()`, 605db781477SPatrick Sanan `PetscViewerHDF5SetSPOutput()`, `PetscViewerHDF5GetBaseDimension2()`, `VecView()`, `MatView()`, `VecLoad()`, 606db781477SPatrick Sanan `MatLoad()`, `PetscFileMode`, `PetscViewer`, `PetscViewerSetType()`, `PetscViewerFileSetMode()`, `PetscViewerFileSetName()` 6075c6c1daeSBarry Smith @*/ 6085c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5Open(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *hdf5v) 6095c6c1daeSBarry Smith { 6105c6c1daeSBarry Smith PetscFunctionBegin; 6119566063dSJacob Faibussowitsch PetscCall(PetscViewerCreate(comm, hdf5v)); 6129566063dSJacob Faibussowitsch PetscCall(PetscViewerSetType(*hdf5v, PETSCVIEWERHDF5)); 6139566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetMode(*hdf5v, type)); 6149566063dSJacob Faibussowitsch PetscCall(PetscViewerFileSetName(*hdf5v, name)); 6159566063dSJacob Faibussowitsch PetscCall(PetscViewerSetFromOptions(*hdf5v)); 6165c6c1daeSBarry Smith PetscFunctionReturn(0); 6175c6c1daeSBarry Smith } 6185c6c1daeSBarry Smith 6195c6c1daeSBarry Smith /*@C 6205c6c1daeSBarry Smith PetscViewerHDF5GetFileId - Retrieve the file id, this file ID then can be used in direct HDF5 calls 6215c6c1daeSBarry Smith 6225c6c1daeSBarry Smith Not collective 6235c6c1daeSBarry Smith 6245c6c1daeSBarry Smith Input Parameter: 6255c6c1daeSBarry Smith . viewer - the PetscViewer 6265c6c1daeSBarry Smith 6275c6c1daeSBarry Smith Output Parameter: 6285c6c1daeSBarry Smith . file_id - The file id 6295c6c1daeSBarry Smith 6305c6c1daeSBarry Smith Level: intermediate 6315c6c1daeSBarry Smith 632db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()` 6335c6c1daeSBarry Smith @*/ 6345c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetFileId(PetscViewer viewer, hid_t *file_id) 6355c6c1daeSBarry Smith { 6365c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6375c6c1daeSBarry Smith 6385c6c1daeSBarry Smith PetscFunctionBegin; 6395c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 6405c6c1daeSBarry Smith if (file_id) *file_id = hdf5->file_id; 6415c6c1daeSBarry Smith PetscFunctionReturn(0); 6425c6c1daeSBarry Smith } 6435c6c1daeSBarry Smith 6445c6c1daeSBarry Smith /*@C 6455c6c1daeSBarry Smith PetscViewerHDF5PushGroup - Set the current HDF5 group for output 6465c6c1daeSBarry Smith 6475c6c1daeSBarry Smith Not collective 6485c6c1daeSBarry Smith 6495c6c1daeSBarry Smith Input Parameters: 6505c6c1daeSBarry Smith + viewer - the PetscViewer 6515c6c1daeSBarry Smith - name - The group name 6525c6c1daeSBarry Smith 6535c6c1daeSBarry Smith Level: intermediate 6545c6c1daeSBarry Smith 6552e361344SVaclav Hapla Notes: 6562e361344SVaclav Hapla This is designed to mnemonically resemble the Unix cd command. 6572e361344SVaclav 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. 6582e361344SVaclav Hapla . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 6592e361344SVaclav Hapla - "." means the current group is pushed again. 6602e361344SVaclav Hapla 6612e361344SVaclav Hapla Example: 6622e361344SVaclav Hapla Suppose the current group is "/a". 6632e361344SVaclav Hapla + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 6642e361344SVaclav Hapla . If name is ".", then the new group will be "/a". 6652e361344SVaclav Hapla . If name is "b", then the new group will be "/a/b". 6662e361344SVaclav Hapla - If name is "/b", then the new group will be "/b". 6672e361344SVaclav Hapla 6682e361344SVaclav Hapla Developer Notes: 6692e361344SVaclav Hapla The root group "/" is internally stored as NULL. 670820fc2d1SVaclav Hapla 671c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 6725c6c1daeSBarry Smith @*/ 673be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5PushGroup(PetscViewer viewer, const char name[]) 6745c6c1daeSBarry Smith { 6755c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 6767d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 6772e361344SVaclav Hapla size_t i,len; 6782e361344SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]; 6792e361344SVaclav Hapla const char *gname; 6805c6c1daeSBarry Smith 6815c6c1daeSBarry Smith PetscFunctionBegin; 6825c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 683820fc2d1SVaclav Hapla if (name) PetscValidCharPointer(name,2); 6849566063dSJacob Faibussowitsch PetscCall(PetscStrlen(name, &len)); 6852e361344SVaclav Hapla gname = NULL; 6862e361344SVaclav Hapla if (len) { 6872e361344SVaclav Hapla if (len == 1 && name[0] == '.') { 6882e361344SVaclav Hapla /* use current name */ 6892e361344SVaclav Hapla gname = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : NULL; 6902e361344SVaclav Hapla } else if (name[0] == '/') { 6912e361344SVaclav Hapla /* absolute */ 6922e361344SVaclav Hapla for (i=1; i<len; i++) { 6932e361344SVaclav Hapla if (name[i] != '/') { 6942e361344SVaclav Hapla gname = name; 6952e361344SVaclav Hapla break; 6962e361344SVaclav Hapla } 6972e361344SVaclav Hapla } 6982e361344SVaclav Hapla } else { 6992e361344SVaclav Hapla /* relative */ 7002e361344SVaclav Hapla const char *parent = (hdf5->groups && hdf5->groups->name) ? hdf5->groups->name : ""; 7019566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf, sizeof(buf), "%s/%s", parent, name)); 7022e361344SVaclav Hapla gname = buf; 7032e361344SVaclav Hapla } 7042e361344SVaclav Hapla } 7059566063dSJacob Faibussowitsch PetscCall(PetscNew(&groupNode)); 7069566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(gname, (char**) &groupNode->name)); 7075c6c1daeSBarry Smith groupNode->next = hdf5->groups; 7085c6c1daeSBarry Smith hdf5->groups = groupNode; 7095c6c1daeSBarry Smith PetscFunctionReturn(0); 7105c6c1daeSBarry Smith } 7115c6c1daeSBarry Smith 712*1ad4d0e4SMatthew G. Knepley /*@C 713*1ad4d0e4SMatthew G. Knepley PetscViewerHDF5PushGroupRelative - Set the current HDF5 group for output, appending the path in the argument 714*1ad4d0e4SMatthew G. Knepley 715*1ad4d0e4SMatthew G. Knepley Not collective 716*1ad4d0e4SMatthew G. Knepley 717*1ad4d0e4SMatthew G. Knepley Input Parameters: 718*1ad4d0e4SMatthew G. Knepley + viewer - the PetscViewer 719*1ad4d0e4SMatthew G. Knepley - name - The group name to append 720*1ad4d0e4SMatthew G. Knepley 721*1ad4d0e4SMatthew G. Knepley Level: intermediate 722*1ad4d0e4SMatthew G. Knepley 723*1ad4d0e4SMatthew G. Knepley Notes: 724*1ad4d0e4SMatthew G. Knepley This is designed to mnemonically resemble the Unix cd command. 725*1ad4d0e4SMatthew G. Knepley + 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. 726*1ad4d0e4SMatthew G. Knepley . NULL, empty string, or any sequence of all slashes (e.g. "///") is interpreted as the root group "/". 727*1ad4d0e4SMatthew G. Knepley - "." means the current group is pushed again. 728*1ad4d0e4SMatthew G. Knepley 729*1ad4d0e4SMatthew G. Knepley Example: 730*1ad4d0e4SMatthew G. Knepley Suppose the current group is "/a". 731*1ad4d0e4SMatthew G. Knepley + If name is NULL, empty string, or a sequence of all slashes (e.g. "///"), then the new group will be "/". 732*1ad4d0e4SMatthew G. Knepley . If name is ".", then the new group will be "/a". 733*1ad4d0e4SMatthew G. Knepley . If name is "b", then the new group will be "/a/b". 734*1ad4d0e4SMatthew G. Knepley - If name is "/b", then the new group will be "/b". 735*1ad4d0e4SMatthew G. Knepley 736*1ad4d0e4SMatthew G. Knepley Developer Notes: 737*1ad4d0e4SMatthew G. Knepley The root group "/" is internally stored as NULL. 738*1ad4d0e4SMatthew G. Knepley 739*1ad4d0e4SMatthew G. Knepley .seealso: `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()`, `PetscViewerHDF5Open()` 740*1ad4d0e4SMatthew G. Knepley @*/ 741*1ad4d0e4SMatthew G. Knepley PetscErrorCode PetscViewerHDF5PushGroupRelative(PetscViewer viewer, const char name[]) 742*1ad4d0e4SMatthew G. Knepley { 743*1ad4d0e4SMatthew G. Knepley PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 744*1ad4d0e4SMatthew G. Knepley PetscViewerHDF5GroupList *groupNode; 745*1ad4d0e4SMatthew G. Knepley char groupname[PETSC_MAX_PATH_LEN]; 746*1ad4d0e4SMatthew G. Knepley 747*1ad4d0e4SMatthew G. Knepley PetscFunctionBegin; 748*1ad4d0e4SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 749*1ad4d0e4SMatthew G. Knepley if (name) PetscValidCharPointer(name,2); 750*1ad4d0e4SMatthew G. Knepley if (name && name[0]) { 751*1ad4d0e4SMatthew G. Knepley size_t i,len; 752*1ad4d0e4SMatthew G. Knepley PetscCall(PetscStrlen(name, &len)); 753*1ad4d0e4SMatthew G. Knepley for (i=0; i<len; i++) if (name[i] != '/') break; 754*1ad4d0e4SMatthew G. Knepley if (i == len) name = NULL; 755*1ad4d0e4SMatthew G. Knepley } else name = NULL; 756*1ad4d0e4SMatthew G. Knepley PetscCall(PetscNew(&groupNode)); 757*1ad4d0e4SMatthew G. Knepley PetscCall(PetscStrncpy(groupname, hdf5->groups->name, sizeof(groupname))); 758*1ad4d0e4SMatthew G. Knepley if (name[0] != '/') PetscCall(PetscStrlcat(groupname, "/", sizeof(groupname))); 759*1ad4d0e4SMatthew G. Knepley PetscCall(PetscStrlcat(groupname, name, sizeof(groupname))); 760*1ad4d0e4SMatthew G. Knepley PetscCall(PetscStrallocpy(groupname, (char**) &groupNode->name)); 761*1ad4d0e4SMatthew G. Knepley groupNode->next = hdf5->groups; 762*1ad4d0e4SMatthew G. Knepley hdf5->groups = groupNode; 763*1ad4d0e4SMatthew G. Knepley PetscFunctionReturn(0); 764*1ad4d0e4SMatthew G. Knepley } 765*1ad4d0e4SMatthew G. Knepley 7663ef9c667SSatish Balay /*@ 7675c6c1daeSBarry Smith PetscViewerHDF5PopGroup - Return the current HDF5 group for output to the previous value 7685c6c1daeSBarry Smith 7695c6c1daeSBarry Smith Not collective 7705c6c1daeSBarry Smith 7715c6c1daeSBarry Smith Input Parameter: 7725c6c1daeSBarry Smith . viewer - the PetscViewer 7735c6c1daeSBarry Smith 7745c6c1daeSBarry Smith Level: intermediate 7755c6c1daeSBarry Smith 776c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 7775c6c1daeSBarry Smith @*/ 7785c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5PopGroup(PetscViewer viewer) 7795c6c1daeSBarry Smith { 7805c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 7817d964842SVaclav Hapla PetscViewerHDF5GroupList *groupNode; 7825c6c1daeSBarry Smith 7835c6c1daeSBarry Smith PetscFunctionBegin; 7845c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 7855f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->groups,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "HDF5 group stack is empty, cannot pop"); 7865c6c1daeSBarry Smith groupNode = hdf5->groups; 7875c6c1daeSBarry Smith hdf5->groups = hdf5->groups->next; 7889566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode->name)); 7899566063dSJacob Faibussowitsch PetscCall(PetscFree(groupNode)); 7905c6c1daeSBarry Smith PetscFunctionReturn(0); 7915c6c1daeSBarry Smith } 7925c6c1daeSBarry Smith 7935c6c1daeSBarry Smith /*@C 794874270d9SVaclav Hapla PetscViewerHDF5GetGroup - Get the current HDF5 group name (full path), set with PetscViewerHDF5PushGroup()/PetscViewerHDF5PopGroup(). 795874270d9SVaclav Hapla If none has been assigned, returns NULL. 7965c6c1daeSBarry Smith 7975c6c1daeSBarry Smith Not collective 7985c6c1daeSBarry Smith 7995c6c1daeSBarry Smith Input Parameter: 8005c6c1daeSBarry Smith . viewer - the PetscViewer 8015c6c1daeSBarry Smith 8025c6c1daeSBarry Smith Output Parameter: 8035c6c1daeSBarry Smith . name - The group name 8045c6c1daeSBarry Smith 8055c6c1daeSBarry Smith Level: intermediate 8065c6c1daeSBarry Smith 807c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5OpenGroup()` 8085c6c1daeSBarry Smith @*/ 809be7aa430SVaclav Hapla PetscErrorCode PetscViewerHDF5GetGroup(PetscViewer viewer, const char *name[]) 8105c6c1daeSBarry Smith { 8115c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5 *) viewer->data; 8125c6c1daeSBarry Smith 8135c6c1daeSBarry Smith PetscFunctionBegin; 8145c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 815c959eef4SJed Brown PetscValidPointer(name,2); 816a297a907SKarl Rupp if (hdf5->groups) *name = hdf5->groups->name; 8170298fd71SBarry Smith else *name = NULL; 8185c6c1daeSBarry Smith PetscFunctionReturn(0); 8195c6c1daeSBarry Smith } 8205c6c1daeSBarry Smith 8218c557b5aSVaclav Hapla /*@ 822874270d9SVaclav Hapla PetscViewerHDF5OpenGroup - Open the HDF5 group with the name (full path) returned by PetscViewerHDF5GetGroup(), 823874270d9SVaclav Hapla and return this group's ID and file ID. 824874270d9SVaclav Hapla If PetscViewerHDF5GetGroup() yields NULL, then group ID is file ID. 825874270d9SVaclav Hapla 826874270d9SVaclav Hapla Not collective 827874270d9SVaclav Hapla 828874270d9SVaclav Hapla Input Parameter: 829874270d9SVaclav Hapla . viewer - the PetscViewer 830874270d9SVaclav Hapla 831d8d19677SJose E. Roman Output Parameters: 832874270d9SVaclav Hapla + fileId - The HDF5 file ID 833874270d9SVaclav Hapla - groupId - The HDF5 group ID 834874270d9SVaclav Hapla 835e74616beSVaclav Hapla Notes: 836e74616beSVaclav Hapla If the viewer is writable, the group is created if it doesn't exist yet. 837e74616beSVaclav Hapla 838874270d9SVaclav Hapla Level: intermediate 839874270d9SVaclav Hapla 840c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 841874270d9SVaclav Hapla @*/ 84254dbf706SVaclav Hapla PetscErrorCode PetscViewerHDF5OpenGroup(PetscViewer viewer, hid_t *fileId, hid_t *groupId) 84354dbf706SVaclav Hapla { 84490e03537SVaclav Hapla hid_t file_id; 84590e03537SVaclav Hapla H5O_type_t type; 84676d59af2SVaclav Hapla const char *groupName = NULL, *fileName = NULL; 84776d59af2SVaclav Hapla PetscBool writable, has; 84854dbf706SVaclav Hapla 84954dbf706SVaclav Hapla PetscFunctionBegin; 8509566063dSJacob Faibussowitsch PetscCall(PetscViewerWritable(viewer, &writable)); 8519566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &file_id)); 8529566063dSJacob Faibussowitsch PetscCall(PetscViewerFileGetName(viewer, &fileName)); 8539566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetGroup(viewer, &groupName)); 8549566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, groupName, writable, &has, &type)); 85576d59af2SVaclav Hapla if (!has) { 8565f80ce2aSJacob 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); 857f7d195e4SLawrence Mitchell SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_LIB, "HDF5 failed to create group %s although file %s is open for writing", groupName, fileName); 85876d59af2SVaclav Hapla } 8595f80ce2aSJacob 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); 86090e03537SVaclav Hapla PetscStackCallHDF5Return(*groupId,H5Gopen2,(file_id, groupName ? groupName : "/", H5P_DEFAULT)); 86154dbf706SVaclav Hapla *fileId = file_id; 86254dbf706SVaclav Hapla PetscFunctionReturn(0); 86354dbf706SVaclav Hapla } 86454dbf706SVaclav Hapla 8653ef9c667SSatish Balay /*@ 866d7dd068bSVaclav Hapla PetscViewerHDF5PushTimestepping - Activate timestepping mode for subsequent HDF5 reading and writing. 8675c6c1daeSBarry Smith 8685c6c1daeSBarry Smith Not collective 8695c6c1daeSBarry Smith 8705c6c1daeSBarry Smith Input Parameter: 8715c6c1daeSBarry Smith . viewer - the PetscViewer 8725c6c1daeSBarry Smith 8735c6c1daeSBarry Smith Level: intermediate 8745c6c1daeSBarry Smith 875d7dd068bSVaclav Hapla Notes: 876d7dd068bSVaclav Hapla On first PetscViewerHDF5PushTimestepping(), the initial time step is set to 0. 877d7dd068bSVaclav Hapla Next timesteps can then be set using PetscViewerHDF5IncrementTimestep() or PetscViewerHDF5SetTimestep(). 878d7dd068bSVaclav 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()]. 879d7dd068bSVaclav Hapla Use PetscViewerHDF5PopTimestepping() to deactivate timestepping mode; calling it by the end of the program is NOT mandatory. 880d7dd068bSVaclav Hapla Current timestep is remembered between PetscViewerHDF5PopTimestepping() and the next PetscViewerHDF5PushTimestepping(). 881d7dd068bSVaclav Hapla 882d7dd068bSVaclav Hapla If a dataset was stored with timestepping, it can be loaded only in the timestepping mode again. 883d7dd068bSVaclav Hapla Loading a timestepped dataset with timestepping disabled, or vice-versa results in an error. 884d7dd068bSVaclav Hapla 885d7dd068bSVaclav Hapla Developer notes: 886d7dd068bSVaclav Hapla Timestepped HDF5 dataset has an extra dimension and attribute "timestepping" set to true. 887d7dd068bSVaclav Hapla 888db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 889d7dd068bSVaclav Hapla @*/ 890d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PushTimestepping(PetscViewer viewer) 891d7dd068bSVaclav Hapla { 892d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 893d7dd068bSVaclav Hapla 894d7dd068bSVaclav Hapla PetscFunctionBegin; 895d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 8965f80ce2aSJacob Faibussowitsch PetscCheck(!hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping is already pushed"); 897d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_TRUE; 898d7dd068bSVaclav Hapla if (hdf5->timestep < 0) hdf5->timestep = 0; 899d7dd068bSVaclav Hapla PetscFunctionReturn(0); 900d7dd068bSVaclav Hapla } 901d7dd068bSVaclav Hapla 902d7dd068bSVaclav Hapla /*@ 903d7dd068bSVaclav Hapla PetscViewerHDF5PopTimestepping - Deactivate timestepping mode for subsequent HDF5 reading and writing. 904d7dd068bSVaclav Hapla 905d7dd068bSVaclav Hapla Not collective 906d7dd068bSVaclav Hapla 907d7dd068bSVaclav Hapla Input Parameter: 908d7dd068bSVaclav Hapla . viewer - the PetscViewer 909d7dd068bSVaclav Hapla 910d7dd068bSVaclav Hapla Level: intermediate 911d7dd068bSVaclav Hapla 912d7dd068bSVaclav Hapla Notes: 913d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 914d7dd068bSVaclav Hapla 915db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IsTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 916d7dd068bSVaclav Hapla @*/ 917d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5PopTimestepping(PetscViewer viewer) 918d7dd068bSVaclav Hapla { 919d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 920d7dd068bSVaclav Hapla 921d7dd068bSVaclav Hapla PetscFunctionBegin; 922d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9235f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 924d7dd068bSVaclav Hapla hdf5->timestepping = PETSC_FALSE; 925d7dd068bSVaclav Hapla PetscFunctionReturn(0); 926d7dd068bSVaclav Hapla } 927d7dd068bSVaclav Hapla 928d7dd068bSVaclav Hapla /*@ 929d7dd068bSVaclav Hapla PetscViewerHDF5IsTimestepping - Ask the viewer whether it is in timestepping mode currently. 930d7dd068bSVaclav Hapla 931d7dd068bSVaclav Hapla Not collective 932d7dd068bSVaclav Hapla 933d7dd068bSVaclav Hapla Input Parameter: 934d7dd068bSVaclav Hapla . viewer - the PetscViewer 935d7dd068bSVaclav Hapla 936d7dd068bSVaclav Hapla Output Parameter: 937d7dd068bSVaclav Hapla . flg - is timestepping active? 938d7dd068bSVaclav Hapla 939d7dd068bSVaclav Hapla Level: intermediate 940d7dd068bSVaclav Hapla 941d7dd068bSVaclav Hapla Notes: 942d7dd068bSVaclav Hapla See PetscViewerHDF5PushTimestepping() for details. 943d7dd068bSVaclav Hapla 944db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5PopTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 945d7dd068bSVaclav Hapla @*/ 946d7dd068bSVaclav Hapla PetscErrorCode PetscViewerHDF5IsTimestepping(PetscViewer viewer, PetscBool *flg) 947d7dd068bSVaclav Hapla { 948d7dd068bSVaclav Hapla PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 949d7dd068bSVaclav Hapla 950d7dd068bSVaclav Hapla PetscFunctionBegin; 951d7dd068bSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 952d7dd068bSVaclav Hapla *flg = hdf5->timestepping; 953d7dd068bSVaclav Hapla PetscFunctionReturn(0); 954d7dd068bSVaclav Hapla } 955d7dd068bSVaclav Hapla 956d7dd068bSVaclav Hapla /*@ 957d7dd068bSVaclav Hapla PetscViewerHDF5IncrementTimestep - Increments current timestep for the HDF5 output. Fields are stacked in time. 958d7dd068bSVaclav Hapla 959d7dd068bSVaclav Hapla Not collective 960d7dd068bSVaclav Hapla 961d7dd068bSVaclav Hapla Input Parameter: 962d7dd068bSVaclav Hapla . viewer - the PetscViewer 963d7dd068bSVaclav Hapla 964d7dd068bSVaclav Hapla Level: intermediate 965d7dd068bSVaclav Hapla 966d7dd068bSVaclav Hapla Notes: 967d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 968d7dd068bSVaclav Hapla 969db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5SetTimestep()`, `PetscViewerHDF5GetTimestep()` 9705c6c1daeSBarry Smith @*/ 9715c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5IncrementTimestep(PetscViewer viewer) 9725c6c1daeSBarry Smith { 9735c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 9745c6c1daeSBarry Smith 9755c6c1daeSBarry Smith PetscFunctionBegin; 9765c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 9775f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 9785c6c1daeSBarry Smith ++hdf5->timestep; 9795c6c1daeSBarry Smith PetscFunctionReturn(0); 9805c6c1daeSBarry Smith } 9815c6c1daeSBarry Smith 9823ef9c667SSatish Balay /*@ 983d7dd068bSVaclav Hapla PetscViewerHDF5SetTimestep - Set the current timestep for the HDF5 output. Fields are stacked in time. 9845c6c1daeSBarry Smith 985d7dd068bSVaclav Hapla Logically collective 9865c6c1daeSBarry Smith 9875c6c1daeSBarry Smith Input Parameters: 9885c6c1daeSBarry Smith + viewer - the PetscViewer 989d7dd068bSVaclav Hapla - timestep - The timestep 9905c6c1daeSBarry Smith 9915c6c1daeSBarry Smith Level: intermediate 9925c6c1daeSBarry Smith 993d7dd068bSVaclav Hapla Notes: 994d7dd068bSVaclav Hapla This can be called only if the viewer is in timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 995d7dd068bSVaclav Hapla 996db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5GetTimestep()` 9975c6c1daeSBarry Smith @*/ 9985c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5SetTimestep(PetscViewer viewer, PetscInt timestep) 9995c6c1daeSBarry Smith { 10005c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 10015c6c1daeSBarry Smith 10025c6c1daeSBarry Smith PetscFunctionBegin; 10035c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1004d7dd068bSVaclav Hapla PetscValidLogicalCollectiveInt(viewer, timestep, 2); 10055f80ce2aSJacob Faibussowitsch PetscCheck(timestep >= 0,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestep %" PetscInt_FMT " is negative", timestep); 10065f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 10075c6c1daeSBarry Smith hdf5->timestep = timestep; 10085c6c1daeSBarry Smith PetscFunctionReturn(0); 10095c6c1daeSBarry Smith } 10105c6c1daeSBarry Smith 10113ef9c667SSatish Balay /*@ 10125c6c1daeSBarry Smith PetscViewerHDF5GetTimestep - Get the current timestep for the HDF5 output. Fields are stacked in time. 10135c6c1daeSBarry Smith 10145c6c1daeSBarry Smith Not collective 10155c6c1daeSBarry Smith 10165c6c1daeSBarry Smith Input Parameter: 10175c6c1daeSBarry Smith . viewer - the PetscViewer 10185c6c1daeSBarry Smith 10195c6c1daeSBarry Smith Output Parameter: 1020d7dd068bSVaclav Hapla . timestep - The timestep 10215c6c1daeSBarry Smith 10225c6c1daeSBarry Smith Level: intermediate 10235c6c1daeSBarry Smith 1024d7dd068bSVaclav Hapla Notes: 1025d7dd068bSVaclav Hapla This can be called only if the viewer is in the timestepping mode. See PetscViewerHDF5PushTimestepping() for details. 1026d7dd068bSVaclav Hapla 1027db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5PushTimestepping()`, `PetscViewerHDF5IncrementTimestep()`, `PetscViewerHDF5SetTimestep()` 10285c6c1daeSBarry Smith @*/ 10295c6c1daeSBarry Smith PetscErrorCode PetscViewerHDF5GetTimestep(PetscViewer viewer, PetscInt *timestep) 10305c6c1daeSBarry Smith { 10315c6c1daeSBarry Smith PetscViewer_HDF5 *hdf5 = (PetscViewer_HDF5*) viewer->data; 10325c6c1daeSBarry Smith 10335c6c1daeSBarry Smith PetscFunctionBegin; 10345c6c1daeSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1035d7dd068bSVaclav Hapla PetscValidIntPointer(timestep,2); 10365f80ce2aSJacob Faibussowitsch PetscCheck(hdf5->timestepping,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "Timestepping has not been pushed yet. Call PetscViewerHDF5PushTimestepping() first"); 10375c6c1daeSBarry Smith *timestep = hdf5->timestep; 10385c6c1daeSBarry Smith PetscFunctionReturn(0); 10395c6c1daeSBarry Smith } 10405c6c1daeSBarry Smith 104136bce6e8SMatthew G. Knepley /*@C 104236bce6e8SMatthew G. Knepley PetscDataTypeToHDF5DataType - Converts the PETSc name of a datatype to its HDF5 name. 104336bce6e8SMatthew G. Knepley 104436bce6e8SMatthew G. Knepley Not collective 104536bce6e8SMatthew G. Knepley 104636bce6e8SMatthew G. Knepley Input Parameter: 104736bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 104836bce6e8SMatthew G. Knepley 104936bce6e8SMatthew G. Knepley Output Parameter: 105036bce6e8SMatthew G. Knepley . mtype - the MPI datatype (for example MPI_DOUBLE, ...) 105136bce6e8SMatthew G. Knepley 105236bce6e8SMatthew G. Knepley Level: advanced 105336bce6e8SMatthew G. Knepley 1054db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 105536bce6e8SMatthew G. Knepley @*/ 105636bce6e8SMatthew G. Knepley PetscErrorCode PetscDataTypeToHDF5DataType(PetscDataType ptype, hid_t *htype) 105736bce6e8SMatthew G. Knepley { 105836bce6e8SMatthew G. Knepley PetscFunctionBegin; 105936bce6e8SMatthew G. Knepley if (ptype == PETSC_INT) 106036bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 106136bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_LLONG; 106236bce6e8SMatthew G. Knepley #else 106336bce6e8SMatthew G. Knepley *htype = H5T_NATIVE_INT; 106436bce6e8SMatthew G. Knepley #endif 106536bce6e8SMatthew G. Knepley else if (ptype == PETSC_DOUBLE) *htype = H5T_NATIVE_DOUBLE; 106636bce6e8SMatthew G. Knepley else if (ptype == PETSC_LONG) *htype = H5T_NATIVE_LONG; 106736bce6e8SMatthew G. Knepley else if (ptype == PETSC_SHORT) *htype = H5T_NATIVE_SHORT; 1068c9450970SBarry Smith else if (ptype == PETSC_ENUM) *htype = H5T_NATIVE_INT; 1069de3ba262SVaclav Hapla else if (ptype == PETSC_BOOL) *htype = H5T_NATIVE_INT; 107036bce6e8SMatthew G. Knepley else if (ptype == PETSC_FLOAT) *htype = H5T_NATIVE_FLOAT; 107136bce6e8SMatthew G. Knepley else if (ptype == PETSC_CHAR) *htype = H5T_NATIVE_CHAR; 107236bce6e8SMatthew G. Knepley else if (ptype == PETSC_BIT_LOGICAL) *htype = H5T_NATIVE_UCHAR; 10737e97c476SMatthew G. Knepley else if (ptype == PETSC_STRING) *htype = H5Tcopy(H5T_C_S1); 107436bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported PETSc datatype"); 107536bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 107636bce6e8SMatthew G. Knepley } 107736bce6e8SMatthew G. Knepley 107836bce6e8SMatthew G. Knepley /*@C 107936bce6e8SMatthew G. Knepley PetscHDF5DataTypeToPetscDataType - Finds the PETSc name of a datatype from its HDF5 name 108036bce6e8SMatthew G. Knepley 108136bce6e8SMatthew G. Knepley Not collective 108236bce6e8SMatthew G. Knepley 108336bce6e8SMatthew G. Knepley Input Parameter: 108436bce6e8SMatthew G. Knepley . htype - the HDF5 datatype (for example H5T_NATIVE_DOUBLE, ...) 108536bce6e8SMatthew G. Knepley 108636bce6e8SMatthew G. Knepley Output Parameter: 108736bce6e8SMatthew G. Knepley . ptype - the PETSc datatype name (for example PETSC_DOUBLE) 108836bce6e8SMatthew G. Knepley 108936bce6e8SMatthew G. Knepley Level: advanced 109036bce6e8SMatthew G. Knepley 1091db781477SPatrick Sanan .seealso: `PetscDataType`, `PetscHDF5DataTypeToPetscDataType()` 109236bce6e8SMatthew G. Knepley @*/ 109336bce6e8SMatthew G. Knepley PetscErrorCode PetscHDF5DataTypeToPetscDataType(hid_t htype, PetscDataType *ptype) 109436bce6e8SMatthew G. Knepley { 109536bce6e8SMatthew G. Knepley PetscFunctionBegin; 109636bce6e8SMatthew G. Knepley #if defined(PETSC_USE_64BIT_INDICES) 109736bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_LONG; 109836bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LLONG) *ptype = PETSC_INT; 109936bce6e8SMatthew G. Knepley #else 110036bce6e8SMatthew G. Knepley if (htype == H5T_NATIVE_INT) *ptype = PETSC_INT; 110136bce6e8SMatthew G. Knepley #endif 110236bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_DOUBLE) *ptype = PETSC_DOUBLE; 110336bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_LONG) *ptype = PETSC_LONG; 110436bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_SHORT) *ptype = PETSC_SHORT; 110536bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_FLOAT) *ptype = PETSC_FLOAT; 110636bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_CHAR) *ptype = PETSC_CHAR; 110736bce6e8SMatthew G. Knepley else if (htype == H5T_NATIVE_UCHAR) *ptype = PETSC_CHAR; 11087e97c476SMatthew G. Knepley else if (htype == H5T_C_S1) *ptype = PETSC_STRING; 110936bce6e8SMatthew G. Knepley else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported HDF5 datatype"); 111036bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 111136bce6e8SMatthew G. Knepley } 111236bce6e8SMatthew G. Knepley 1113df863907SAlex Fikl /*@C 1114b8ee85a1SVaclav Hapla PetscViewerHDF5WriteAttribute - Write an attribute 111536bce6e8SMatthew G. Knepley 1116343a20b2SVaclav Hapla Collective 1117343a20b2SVaclav Hapla 111836bce6e8SMatthew G. Knepley Input Parameters: 111936bce6e8SMatthew G. Knepley + viewer - The HDF5 viewer 11204302210dSVaclav Hapla . parent - The parent dataset/group name 112136bce6e8SMatthew G. Knepley . name - The attribute name 112236bce6e8SMatthew G. Knepley . datatype - The attribute type 112336bce6e8SMatthew G. Knepley - value - The attribute value 112436bce6e8SMatthew G. Knepley 112536bce6e8SMatthew G. Knepley Level: advanced 112636bce6e8SMatthew G. Knepley 11274302210dSVaclav Hapla Notes: 1128343a20b2SVaclav 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. 11294302210dSVaclav Hapla 1130c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 113136bce6e8SMatthew G. Knepley @*/ 11324302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5WriteAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *value) 113336bce6e8SMatthew G. Knepley { 11344302210dSVaclav Hapla char *parentAbsPath; 113560568592SMatthew G. Knepley hid_t h5, dataspace, obj, attribute, dtype; 1136080f144cSVaclav Hapla PetscBool has; 113736bce6e8SMatthew G. Knepley 113836bce6e8SMatthew G. Knepley PetscFunctionBegin; 11395cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 11404302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1141c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 11424302210dSVaclav Hapla PetscValidLogicalCollectiveEnum(viewer,datatype,4); 1143b7e81f0fSVaclav Hapla PetscValidPointer(value, 5); 11449566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 11459566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_TRUE, NULL, NULL)); 11469566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 11479566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 11487e97c476SMatthew G. Knepley if (datatype == PETSC_STRING) { 11497e97c476SMatthew G. Knepley size_t len; 11509566063dSJacob Faibussowitsch PetscCall(PetscStrlen((const char *) value, &len)); 1151729ab6d9SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 11527e97c476SMatthew G. Knepley } 11539566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 1154729ab6d9SBarry Smith PetscStackCallHDF5Return(dataspace,H5Screate,(H5S_SCALAR)); 11554302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 1156080f144cSVaclav Hapla if (has) { 1157080f144cSVaclav Hapla PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1158080f144cSVaclav Hapla } else { 115960568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Acreate2,(obj, name, dtype, dataspace, H5P_DEFAULT, H5P_DEFAULT)); 1160080f144cSVaclav Hapla } 1161729ab6d9SBarry Smith PetscStackCallHDF5(H5Awrite,(attribute, dtype, value)); 1162729ab6d9SBarry Smith if (datatype == PETSC_STRING) PetscStackCallHDF5(H5Tclose,(dtype)); 1163729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 116460568592SMatthew G. Knepley PetscStackCallHDF5(H5Oclose,(obj)); 1165729ab6d9SBarry Smith PetscStackCallHDF5(H5Sclose,(dataspace)); 11669566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 116736bce6e8SMatthew G. Knepley PetscFunctionReturn(0); 116836bce6e8SMatthew G. Knepley } 116936bce6e8SMatthew G. Knepley 1170df863907SAlex Fikl /*@C 1171577e0e04SVaclav Hapla PetscViewerHDF5WriteObjectAttribute - Write an attribute to the dataset matching the given PetscObject by name 1172577e0e04SVaclav Hapla 1173343a20b2SVaclav Hapla Collective 1174343a20b2SVaclav Hapla 1175577e0e04SVaclav Hapla Input Parameters: 1176577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1177577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1178577e0e04SVaclav Hapla . name - The attribute name 1179577e0e04SVaclav Hapla . datatype - The attribute type 1180577e0e04SVaclav Hapla - value - The attribute value 1181577e0e04SVaclav Hapla 1182577e0e04SVaclav Hapla Notes: 1183343a20b2SVaclav 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). 1184577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1185577e0e04SVaclav Hapla 1186577e0e04SVaclav Hapla Level: advanced 1187577e0e04SVaclav Hapla 1188c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1189577e0e04SVaclav Hapla @*/ 1190577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5WriteObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, const void *value) 1191577e0e04SVaclav Hapla { 1192577e0e04SVaclav Hapla PetscFunctionBegin; 1193577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1194577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1195577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1196b7e81f0fSVaclav Hapla PetscValidPointer(value,5); 11979566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 11989566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5WriteAttribute(viewer, obj->name, name, datatype, value)); 1199577e0e04SVaclav Hapla PetscFunctionReturn(0); 1200577e0e04SVaclav Hapla } 1201577e0e04SVaclav Hapla 1202577e0e04SVaclav Hapla /*@C 1203b8ee85a1SVaclav Hapla PetscViewerHDF5ReadAttribute - Read an attribute 1204ad0c61feSMatthew G. Knepley 1205343a20b2SVaclav Hapla Collective 1206343a20b2SVaclav Hapla 1207ad0c61feSMatthew G. Knepley Input Parameters: 1208ad0c61feSMatthew G. Knepley + viewer - The HDF5 viewer 12094302210dSVaclav Hapla . parent - The parent dataset/group name 1210ad0c61feSMatthew G. Knepley . name - The attribute name 1211a2d6be1bSVaclav Hapla . datatype - The attribute type 1212a2d6be1bSVaclav Hapla - defaultValue - The pointer to the default value 1213ad0c61feSMatthew G. Knepley 1214ad0c61feSMatthew G. Knepley Output Parameter: 1215a2d6be1bSVaclav Hapla . value - The pointer to the read HDF5 attribute value 1216ad0c61feSMatthew G. Knepley 1217a2d6be1bSVaclav Hapla Notes: 1218a2d6be1bSVaclav Hapla If defaultValue is NULL and the attribute is not found, an error occurs. 1219a2d6be1bSVaclav Hapla If defaultValue is not NULL and the attribute is not found, defaultValue is copied to value. 1220a2d6be1bSVaclav Hapla The pointers defaultValue and value can be the same; for instance 1221a2d6be1bSVaclav Hapla $ flg = PETSC_FALSE; 12229566063dSJacob Faibussowitsch $ PetscCall(PetscViewerHDF5ReadAttribute(viewer,name,"attr",PETSC_BOOL,&flg,&flg)); 1223a2d6be1bSVaclav Hapla is valid, but make sure the default value is initialized. 1224a2d6be1bSVaclav Hapla 1225a2d6be1bSVaclav Hapla If the datatype is PETSC_STRING, the output string is newly allocated so one must PetscFree() it when no longer needed. 1226708d4cb3SBarry Smith 1227343a20b2SVaclav 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. 1228ad0c61feSMatthew G. Knepley 1229343a20b2SVaclav Hapla Level: advanced 12304302210dSVaclav Hapla 1231c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1232ad0c61feSMatthew G. Knepley @*/ 1233d70ec8faSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadAttribute(PetscViewer viewer, const char parent[], const char name[], PetscDataType datatype, const void *defaultValue, void *value) 1234ad0c61feSMatthew G. Knepley { 12354302210dSVaclav Hapla char *parentAbsPath; 1236a2d6be1bSVaclav Hapla hid_t h5, obj, attribute, dtype; 1237969748fdSVaclav Hapla PetscBool has; 1238ad0c61feSMatthew G. Knepley 1239ad0c61feSMatthew G. Knepley PetscFunctionBegin; 12405cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 12414302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent, 2); 1242c57a1a9aSVaclav Hapla PetscValidCharPointer(name, 3); 1243a2d6be1bSVaclav Hapla if (defaultValue) PetscValidPointer(defaultValue, 5); 1244a2d6be1bSVaclav Hapla PetscValidPointer(value, 6); 12459566063dSJacob Faibussowitsch PetscCall(PetscDataTypeToHDF5DataType(datatype, &dtype)); 12469566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 12479566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, &has, NULL)); 12489566063dSJacob Faibussowitsch if (has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, &has)); 1249a2d6be1bSVaclav Hapla if (!has) { 1250a2d6be1bSVaclav Hapla if (defaultValue) { 1251a2d6be1bSVaclav Hapla if (defaultValue != value) { 1252a2d6be1bSVaclav Hapla if (datatype == PETSC_STRING) { 12539566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(*(char**)defaultValue, (char**)value)); 1254a2d6be1bSVaclav Hapla } else { 1255a2d6be1bSVaclav Hapla size_t len; 1256ece88022SPierre Jolivet PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(dtype)); 12579566063dSJacob Faibussowitsch PetscCall(PetscMemcpy(value, defaultValue, len)); 1258a2d6be1bSVaclav Hapla } 1259a2d6be1bSVaclav Hapla } 12609566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1261a2d6be1bSVaclav Hapla PetscFunctionReturn(0); 126298921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_FILE_UNEXPECTED, "Attribute %s/%s does not exist and default value not provided", parentAbsPath, name); 1263a2d6be1bSVaclav Hapla } 12649566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 12654302210dSVaclav Hapla PetscStackCallHDF5Return(obj,H5Oopen,(h5, parentAbsPath, H5P_DEFAULT)); 126660568592SMatthew G. Knepley PetscStackCallHDF5Return(attribute,H5Aopen_name,(obj, name)); 1267f0b58479SMatthew G. Knepley if (datatype == PETSC_STRING) { 1268f0b58479SMatthew G. Knepley size_t len; 1269a2d6be1bSVaclav Hapla hid_t atype; 1270e90c5075SMarek Pecha PetscStackCallHDF5Return(atype,H5Aget_type,(attribute)); 1271ece88022SPierre Jolivet PetscStackCallHDF5ReturnNoCheck(len,H5Tget_size,(atype)); 12729566063dSJacob Faibussowitsch PetscCall(PetscMalloc((len+1) * sizeof(char), value)); 1273708d4cb3SBarry Smith PetscStackCallHDF5(H5Tset_size,(dtype, len+1)); 1274708d4cb3SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, *(char**)value)); 1275708d4cb3SBarry Smith } else { 127670efba33SBarry Smith PetscStackCallHDF5(H5Aread,(attribute, dtype, value)); 1277708d4cb3SBarry Smith } 1278729ab6d9SBarry Smith PetscStackCallHDF5(H5Aclose,(attribute)); 1279e90c5075SMarek Pecha /* H5Oclose can be used to close groups, datasets, or committed datatypes */ 1280e90c5075SMarek Pecha PetscStackCallHDF5(H5Oclose,(obj)); 12819566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 1282ad0c61feSMatthew G. Knepley PetscFunctionReturn(0); 1283ad0c61feSMatthew G. Knepley } 1284ad0c61feSMatthew G. Knepley 1285577e0e04SVaclav Hapla /*@C 1286577e0e04SVaclav Hapla PetscViewerHDF5ReadObjectAttribute - Read an attribute from the dataset matching the given PetscObject by name 1287577e0e04SVaclav Hapla 1288343a20b2SVaclav Hapla Collective 1289343a20b2SVaclav Hapla 1290577e0e04SVaclav Hapla Input Parameters: 1291577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1292577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1293577e0e04SVaclav Hapla . name - The attribute name 1294577e0e04SVaclav Hapla - datatype - The attribute type 1295577e0e04SVaclav Hapla 1296577e0e04SVaclav Hapla Output Parameter: 1297577e0e04SVaclav Hapla . value - The attribute value 1298577e0e04SVaclav Hapla 1299577e0e04SVaclav Hapla Notes: 1300577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1301577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1302577e0e04SVaclav Hapla 1303577e0e04SVaclav Hapla Level: advanced 1304577e0e04SVaclav Hapla 1305c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5ReadAttribute()` `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1306577e0e04SVaclav Hapla @*/ 1307a2d6be1bSVaclav Hapla PetscErrorCode PetscViewerHDF5ReadObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscDataType datatype, void *defaultValue, void *value) 1308577e0e04SVaclav Hapla { 1309577e0e04SVaclav Hapla PetscFunctionBegin; 1310577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1311577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1312577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 1313064a246eSJacob Faibussowitsch PetscValidPointer(value, 6); 13149566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 13159566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5ReadAttribute(viewer, obj->name, name, datatype, defaultValue, value)); 1316577e0e04SVaclav Hapla PetscFunctionReturn(0); 1317577e0e04SVaclav Hapla } 1318577e0e04SVaclav Hapla 13199fbee547SJacob Faibussowitsch static inline PetscErrorCode PetscViewerHDF5Traverse_Inner_Internal(hid_t h5, const char name[], PetscBool createGroup, PetscBool *exists_) 1320a75c3fd4SVaclav Hapla { 1321a75c3fd4SVaclav Hapla htri_t exists; 1322a75c3fd4SVaclav Hapla hid_t group; 1323a75c3fd4SVaclav Hapla 1324a75c3fd4SVaclav Hapla PetscFunctionBegin; 1325a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(exists,H5Lexists,(h5, name, H5P_DEFAULT)); 1326a75c3fd4SVaclav Hapla if (exists) PetscStackCallHDF5Return(exists,H5Oexists_by_name,(h5, name, H5P_DEFAULT)); 1327a75c3fd4SVaclav Hapla if (!exists && createGroup) { 1328a75c3fd4SVaclav Hapla PetscStackCallHDF5Return(group,H5Gcreate2,(h5, name, H5P_DEFAULT, H5P_DEFAULT, H5P_DEFAULT)); 1329a75c3fd4SVaclav Hapla PetscStackCallHDF5(H5Gclose,(group)); 1330a75c3fd4SVaclav Hapla exists = PETSC_TRUE; 1331a75c3fd4SVaclav Hapla } 1332a75c3fd4SVaclav Hapla *exists_ = (PetscBool) exists; 1333a75c3fd4SVaclav Hapla PetscFunctionReturn(0); 1334a75c3fd4SVaclav Hapla } 1335a75c3fd4SVaclav Hapla 1336bb286ee1SVaclav Hapla static PetscErrorCode PetscViewerHDF5Traverse_Internal(PetscViewer viewer, const char name[], PetscBool createGroup, PetscBool *has, H5O_type_t *otype) 13375cdeb108SMatthew G. Knepley { 133890e03537SVaclav Hapla const char rootGroupName[] = "/"; 13395cdeb108SMatthew G. Knepley hid_t h5; 1340e5bf9ebcSVaclav Hapla PetscBool exists=PETSC_FALSE; 134115b861d2SVaclav Hapla PetscInt i; 134215b861d2SVaclav Hapla int n; 134385688ae2SVaclav Hapla char **hierarchy; 134485688ae2SVaclav Hapla char buf[PETSC_MAX_PATH_LEN]=""; 13455cdeb108SMatthew G. Knepley 13465cdeb108SMatthew G. Knepley PetscFunctionBegin; 13475cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 134890e03537SVaclav Hapla if (name) PetscValidCharPointer(name, 2); 134990e03537SVaclav Hapla else name = rootGroupName; 1350ccd66a83SVaclav Hapla if (has) { 1351064a246eSJacob Faibussowitsch PetscValidBoolPointer(has, 4); 13525cdeb108SMatthew G. Knepley *has = PETSC_FALSE; 1353ccd66a83SVaclav Hapla } 1354ccd66a83SVaclav Hapla if (otype) { 1355064a246eSJacob Faibussowitsch PetscValidIntPointer(otype, 5); 135656cc0592SVaclav Hapla *otype = H5O_TYPE_UNKNOWN; 1357ccd66a83SVaclav Hapla } 13589566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 135985688ae2SVaclav Hapla 136085688ae2SVaclav Hapla /* 136185688ae2SVaclav Hapla Unfortunately, H5Oexists_by_name() fails if any object in hierarchy is missing. 136285688ae2SVaclav Hapla Hence, each of them needs to be tested separately: 136385688ae2SVaclav Hapla 1) whether it's a valid link 136485688ae2SVaclav Hapla 2) whether this link resolves to an object 136585688ae2SVaclav Hapla See H5Oexists_by_name() documentation. 136685688ae2SVaclav Hapla */ 13679566063dSJacob Faibussowitsch PetscCall(PetscStrToArray(name,'/',&n,&hierarchy)); 136885688ae2SVaclav Hapla if (!n) { 136985688ae2SVaclav Hapla /* Assume group "/" always exists in accordance with HDF5 >= 1.10.0. See H5Lexists() documentation. */ 1370ccd66a83SVaclav Hapla if (has) *has = PETSC_TRUE; 1371ccd66a83SVaclav Hapla if (otype) *otype = H5O_TYPE_GROUP; 13729566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n,hierarchy)); 137385688ae2SVaclav Hapla PetscFunctionReturn(0); 137485688ae2SVaclav Hapla } 137585688ae2SVaclav Hapla for (i=0; i<n; i++) { 13769566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf,"/")); 13779566063dSJacob Faibussowitsch PetscCall(PetscStrcat(buf,hierarchy[i])); 13789566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Inner_Internal(h5, buf, createGroup, &exists)); 1379a75c3fd4SVaclav Hapla if (!exists) break; 138085688ae2SVaclav Hapla } 13819566063dSJacob Faibussowitsch PetscCall(PetscStrToArrayDestroy(n,hierarchy)); 138285688ae2SVaclav Hapla 138385688ae2SVaclav Hapla /* If the object exists, get its type */ 1384ccd66a83SVaclav Hapla if (exists && otype) { 13855cdeb108SMatthew G. Knepley H5O_info_t info; 13865cdeb108SMatthew G. Knepley 138776276748SVaclav Hapla /* We could use H5Iget_type() here but that would require opening the object. This way we only need its name. */ 138804633f7fSVaclav Hapla PetscStackCallHDF5(H5Oget_info_by_name,(h5, name, &info, H5P_DEFAULT)); 138956cc0592SVaclav Hapla *otype = info.type; 13905cdeb108SMatthew G. Knepley } 1391ccd66a83SVaclav Hapla if (has) *has = exists; 13925cdeb108SMatthew G. Knepley PetscFunctionReturn(0); 13935cdeb108SMatthew G. Knepley } 13945cdeb108SMatthew G. Knepley 13954302210dSVaclav Hapla /*@C 13960a9f38efSVaclav Hapla PetscViewerHDF5HasGroup - Check whether the current (pushed) group exists in the HDF5 file 13970a9f38efSVaclav Hapla 1398343a20b2SVaclav Hapla Collective 1399343a20b2SVaclav Hapla 14000a9f38efSVaclav Hapla Input Parameters: 14014302210dSVaclav Hapla + viewer - The HDF5 viewer 14024302210dSVaclav Hapla - path - The group path 14030a9f38efSVaclav Hapla 14040a9f38efSVaclav Hapla Output Parameter: 14050a9f38efSVaclav Hapla . has - Flag for group existence 14060a9f38efSVaclav Hapla 14070a9f38efSVaclav Hapla Level: advanced 14080a9f38efSVaclav Hapla 14094302210dSVaclav Hapla Notes: 1410785c443eSVaclav 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. 1411785c443eSVaclav Hapla So NULL or empty path means the current pushed group. 1412343a20b2SVaclav Hapla If path exists but is not a group, PETSC_FALSE is returned. 14134302210dSVaclav Hapla 1414db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()`, `PetscViewerHDF5OpenGroup()` 14150a9f38efSVaclav Hapla @*/ 14164302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasGroup(PetscViewer viewer, const char path[], PetscBool *has) 14170a9f38efSVaclav Hapla { 14180a9f38efSVaclav Hapla H5O_type_t type; 14194302210dSVaclav Hapla char *abspath; 14200a9f38efSVaclav Hapla 14210a9f38efSVaclav Hapla PetscFunctionBegin; 14220a9f38efSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 14234302210dSVaclav Hapla if (path) PetscValidCharPointer(path,2); 14244302210dSVaclav Hapla PetscValidBoolPointer(has,3); 14259566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath)); 14269566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 14274302210dSVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_GROUP); 14289566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 14290a9f38efSVaclav Hapla PetscFunctionReturn(0); 14300a9f38efSVaclav Hapla } 14310a9f38efSVaclav Hapla 143289e0ef10SVaclav Hapla /*@C 143389e0ef10SVaclav Hapla PetscViewerHDF5HasDataset - Check whether a given dataset exists in the HDF5 file 143489e0ef10SVaclav Hapla 1435343a20b2SVaclav Hapla Collective 1436343a20b2SVaclav Hapla 143789e0ef10SVaclav Hapla Input Parameters: 143889e0ef10SVaclav Hapla + viewer - The HDF5 viewer 143989e0ef10SVaclav Hapla - path - The dataset path 144089e0ef10SVaclav Hapla 144189e0ef10SVaclav Hapla Output Parameter: 144289e0ef10SVaclav Hapla . has - Flag whether dataset exists 144389e0ef10SVaclav Hapla 144489e0ef10SVaclav Hapla Level: advanced 144589e0ef10SVaclav Hapla 144689e0ef10SVaclav Hapla Notes: 1447343a20b2SVaclav 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. 144889e0ef10SVaclav Hapla If path is NULL or empty, has is set to PETSC_FALSE. 1449343a20b2SVaclav Hapla If path exists but is not a dataset, has is set to PETSC_FALSE as well. 145089e0ef10SVaclav Hapla 1451db781477SPatrick Sanan .seealso: `PetscViewerHDF5HasObject()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5HasGroup()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 145289e0ef10SVaclav Hapla @*/ 145389e0ef10SVaclav Hapla PetscErrorCode PetscViewerHDF5HasDataset(PetscViewer viewer, const char path[], PetscBool *has) 145489e0ef10SVaclav Hapla { 145589e0ef10SVaclav Hapla H5O_type_t type; 145689e0ef10SVaclav Hapla char *abspath; 145789e0ef10SVaclav Hapla 145889e0ef10SVaclav Hapla PetscFunctionBegin; 145989e0ef10SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 146089e0ef10SVaclav Hapla if (path) PetscValidCharPointer(path,2); 146189e0ef10SVaclav Hapla PetscValidBoolPointer(has,3); 14629566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, path, &abspath)); 14639566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, abspath, PETSC_FALSE, NULL, &type)); 146489e0ef10SVaclav Hapla *has = (PetscBool)(type == H5O_TYPE_DATASET); 14659566063dSJacob Faibussowitsch PetscCall(PetscFree(abspath)); 146689e0ef10SVaclav Hapla PetscFunctionReturn(0); 146789e0ef10SVaclav Hapla } 146889e0ef10SVaclav Hapla 14690a9f38efSVaclav Hapla /*@ 1470e3f143f7SVaclav Hapla PetscViewerHDF5HasObject - Check whether a dataset with the same name as given object exists in the HDF5 file under current group 1471ecce1506SVaclav Hapla 1472343a20b2SVaclav Hapla Collective 1473343a20b2SVaclav Hapla 1474ecce1506SVaclav Hapla Input Parameters: 1475ecce1506SVaclav Hapla + viewer - The HDF5 viewer 1476ecce1506SVaclav Hapla - obj - The named object 1477ecce1506SVaclav Hapla 1478ecce1506SVaclav Hapla Output Parameter: 147989e0ef10SVaclav Hapla . has - Flag for dataset existence 1480ecce1506SVaclav Hapla 1481e3f143f7SVaclav Hapla Notes: 148289e0ef10SVaclav Hapla If the object is unnamed, an error occurs. 1483343a20b2SVaclav Hapla If the path current_group/object_name exists but is not a dataset, has is set to PETSC_FALSE as well. 1484e3f143f7SVaclav Hapla 1485ecce1506SVaclav Hapla Level: advanced 1486ecce1506SVaclav Hapla 1487db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasDataset()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1488ecce1506SVaclav Hapla @*/ 1489ecce1506SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObject(PetscViewer viewer, PetscObject obj, PetscBool *has) 1490ecce1506SVaclav Hapla { 149189e0ef10SVaclav Hapla size_t len; 1492ecce1506SVaclav Hapla 1493ecce1506SVaclav Hapla PetscFunctionBegin; 1494c57a1a9aSVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1495c57a1a9aSVaclav Hapla PetscValidHeader(obj,2); 14964302210dSVaclav Hapla PetscValidBoolPointer(has,3); 14979566063dSJacob Faibussowitsch PetscCall(PetscStrlen(obj->name, &len)); 14985f80ce2aSJacob Faibussowitsch PetscCheck(len,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Object must be named"); 14999566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasDataset(viewer, obj->name, has)); 1500ecce1506SVaclav Hapla PetscFunctionReturn(0); 1501ecce1506SVaclav Hapla } 1502ecce1506SVaclav Hapla 1503df863907SAlex Fikl /*@C 1504ece83b6cSVaclav Hapla PetscViewerHDF5HasAttribute - Check whether an attribute exists 1505e7bdbf8eSMatthew G. Knepley 1506343a20b2SVaclav Hapla Collective 1507343a20b2SVaclav Hapla 1508e7bdbf8eSMatthew G. Knepley Input Parameters: 1509e7bdbf8eSMatthew G. Knepley + viewer - The HDF5 viewer 15104302210dSVaclav Hapla . parent - The parent dataset/group name 1511e7bdbf8eSMatthew G. Knepley - name - The attribute name 1512e7bdbf8eSMatthew G. Knepley 1513e7bdbf8eSMatthew G. Knepley Output Parameter: 1514e7bdbf8eSMatthew G. Knepley . has - Flag for attribute existence 1515e7bdbf8eSMatthew G. Knepley 1516e7bdbf8eSMatthew G. Knepley Level: advanced 1517e7bdbf8eSMatthew G. Knepley 15184302210dSVaclav Hapla Notes: 1519343a20b2SVaclav 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. 15204302210dSVaclav Hapla 1521c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasObjectAttribute()`, `PetscViewerHDF5WriteAttribute()`, `PetscViewerHDF5ReadAttribute()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1522e7bdbf8eSMatthew G. Knepley @*/ 15234302210dSVaclav Hapla PetscErrorCode PetscViewerHDF5HasAttribute(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 1524e7bdbf8eSMatthew G. Knepley { 15254302210dSVaclav Hapla char *parentAbsPath; 1526e7bdbf8eSMatthew G. Knepley 1527e7bdbf8eSMatthew G. Knepley PetscFunctionBegin; 15285cdeb108SMatthew G. Knepley PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 15294302210dSVaclav Hapla if (parent) PetscValidCharPointer(parent,2); 1530c57a1a9aSVaclav Hapla PetscValidCharPointer(name,3); 15314302210dSVaclav Hapla PetscValidBoolPointer(has,4); 15329566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetAbsolutePath_Internal(viewer, parent, &parentAbsPath)); 15339566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5Traverse_Internal(viewer, parentAbsPath, PETSC_FALSE, has, NULL)); 15349566063dSJacob Faibussowitsch if (*has) PetscCall(PetscViewerHDF5HasAttribute_Internal(viewer, parentAbsPath, name, has)); 15359566063dSJacob Faibussowitsch PetscCall(PetscFree(parentAbsPath)); 153606db490cSVaclav Hapla PetscFunctionReturn(0); 153706db490cSVaclav Hapla } 153806db490cSVaclav Hapla 1539577e0e04SVaclav Hapla /*@C 1540577e0e04SVaclav Hapla PetscViewerHDF5HasObjectAttribute - Check whether an attribute is attached to the dataset matching the given PetscObject by name 1541577e0e04SVaclav Hapla 1542343a20b2SVaclav Hapla Collective 1543343a20b2SVaclav Hapla 1544577e0e04SVaclav Hapla Input Parameters: 1545577e0e04SVaclav Hapla + viewer - The HDF5 viewer 1546577e0e04SVaclav Hapla . obj - The object whose name is used to lookup the parent dataset, relative to the current group. 1547577e0e04SVaclav Hapla - name - The attribute name 1548577e0e04SVaclav Hapla 1549577e0e04SVaclav Hapla Output Parameter: 1550577e0e04SVaclav Hapla . has - Flag for attribute existence 1551577e0e04SVaclav Hapla 1552577e0e04SVaclav Hapla Notes: 1553577e0e04SVaclav Hapla This fails if current_group/object_name doesn't resolve to a dataset (the path doesn't exist or is not a dataset). 1554577e0e04SVaclav Hapla You might want to check first if it does using PetscViewerHDF5HasObject(). 1555577e0e04SVaclav Hapla 1556577e0e04SVaclav Hapla Level: advanced 1557577e0e04SVaclav Hapla 1558c2e3fba1SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerHDF5HasAttribute()`, `PetscViewerHDF5WriteObjectAttribute()`, `PetscViewerHDF5ReadObjectAttribute()`, `PetscViewerHDF5HasObject()`, `PetscViewerHDF5PushGroup()`, `PetscViewerHDF5PopGroup()`, `PetscViewerHDF5GetGroup()` 1559577e0e04SVaclav Hapla @*/ 1560577e0e04SVaclav Hapla PetscErrorCode PetscViewerHDF5HasObjectAttribute(PetscViewer viewer, PetscObject obj, const char name[], PetscBool *has) 1561577e0e04SVaclav Hapla { 1562577e0e04SVaclav Hapla PetscFunctionBegin; 1563577e0e04SVaclav Hapla PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,1); 1564577e0e04SVaclav Hapla PetscValidHeader(obj,2); 1565577e0e04SVaclav Hapla PetscValidCharPointer(name,3); 15664302210dSVaclav Hapla PetscValidBoolPointer(has,4); 15679566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5CheckNamedObject_Internal(viewer, obj)); 15689566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5HasAttribute(viewer, obj->name, name, has)); 1569577e0e04SVaclav Hapla PetscFunctionReturn(0); 1570577e0e04SVaclav Hapla } 1571577e0e04SVaclav Hapla 157206db490cSVaclav Hapla static PetscErrorCode PetscViewerHDF5HasAttribute_Internal(PetscViewer viewer, const char parent[], const char name[], PetscBool *has) 157306db490cSVaclav Hapla { 157406db490cSVaclav Hapla hid_t h5; 157506db490cSVaclav Hapla htri_t hhas; 157606db490cSVaclav Hapla 157706db490cSVaclav Hapla PetscFunctionBegin; 15789566063dSJacob Faibussowitsch PetscCall(PetscViewerHDF5GetFileId(viewer, &h5)); 15792f936e54SVaclav Hapla PetscStackCallHDF5Return(hhas,H5Aexists_by_name,(h5, parent, name, H5P_DEFAULT)); 15805cdeb108SMatthew G. Knepley *has = hhas ? PETSC_TRUE : PETSC_FALSE; 1581e7bdbf8eSMatthew G. Knepley PetscFunctionReturn(0); 1582e7bdbf8eSMatthew G. Knepley } 1583e7bdbf8eSMatthew G. Knepley 1584a75e6a4aSMatthew G. Knepley /* 1585a75e6a4aSMatthew G. Knepley The variable Petsc_Viewer_HDF5_keyval is used to indicate an MPI attribute that 1586a75e6a4aSMatthew G. Knepley is attached to a communicator, in this case the attribute is a PetscViewer. 1587a75e6a4aSMatthew G. Knepley */ 1588d4c7638eSBarry Smith PetscMPIInt Petsc_Viewer_HDF5_keyval = MPI_KEYVAL_INVALID; 1589a75e6a4aSMatthew G. Knepley 1590a75e6a4aSMatthew G. Knepley /*@C 1591a75e6a4aSMatthew G. Knepley PETSC_VIEWER_HDF5_ - Creates an HDF5 PetscViewer shared by all processors in a communicator. 1592a75e6a4aSMatthew G. Knepley 1593d083f849SBarry Smith Collective 1594a75e6a4aSMatthew G. Knepley 1595a75e6a4aSMatthew G. Knepley Input Parameter: 1596a75e6a4aSMatthew G. Knepley . comm - the MPI communicator to share the HDF5 PetscViewer 1597a75e6a4aSMatthew G. Knepley 1598a75e6a4aSMatthew G. Knepley Level: intermediate 1599a75e6a4aSMatthew G. Knepley 1600a75e6a4aSMatthew G. Knepley Options Database Keys: 160110699b91SBarry Smith . -viewer_hdf5_filename <name> - name of the HDF5 file 1602a75e6a4aSMatthew G. Knepley 1603a75e6a4aSMatthew G. Knepley Environmental variables: 160410699b91SBarry Smith . PETSC_VIEWER_HDF5_FILENAME - name of the HDF5 file 1605a75e6a4aSMatthew G. Knepley 1606a75e6a4aSMatthew G. Knepley Notes: 1607a75e6a4aSMatthew G. Knepley Unlike almost all other PETSc routines, PETSC_VIEWER_HDF5_ does not return 1608a75e6a4aSMatthew G. Knepley an error code. The HDF5 PetscViewer is usually used in the form 1609a75e6a4aSMatthew G. Knepley $ XXXView(XXX object, PETSC_VIEWER_HDF5_(comm)); 1610a75e6a4aSMatthew G. Knepley 1611db781477SPatrick Sanan .seealso: `PetscViewerHDF5Open()`, `PetscViewerCreate()`, `PetscViewerDestroy()` 1612a75e6a4aSMatthew G. Knepley @*/ 1613a75e6a4aSMatthew G. Knepley PetscViewer PETSC_VIEWER_HDF5_(MPI_Comm comm) 1614a75e6a4aSMatthew G. Knepley { 1615a75e6a4aSMatthew G. Knepley PetscErrorCode ierr; 1616a75e6a4aSMatthew G. Knepley PetscBool flg; 1617a75e6a4aSMatthew G. Knepley PetscViewer viewer; 1618a75e6a4aSMatthew G. Knepley char fname[PETSC_MAX_PATH_LEN]; 1619a75e6a4aSMatthew G. Knepley MPI_Comm ncomm; 1620a75e6a4aSMatthew G. Knepley 1621a75e6a4aSMatthew G. Knepley PetscFunctionBegin; 1622908793a3SLisandro 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);} 1623a75e6a4aSMatthew G. Knepley if (Petsc_Viewer_HDF5_keyval == MPI_KEYVAL_INVALID) { 1624908793a3SLisandro Dalcin ierr = MPI_Comm_create_keyval(MPI_COMM_NULL_COPY_FN,MPI_COMM_NULL_DELETE_FN,&Petsc_Viewer_HDF5_keyval,NULL); 1625908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1626a75e6a4aSMatthew G. Knepley } 162747435625SJed Brown ierr = MPI_Comm_get_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void**)&viewer,(int*)&flg); 1628908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1629a75e6a4aSMatthew G. Knepley if (!flg) { /* PetscViewer not yet created */ 1630a75e6a4aSMatthew G. Knepley ierr = PetscOptionsGetenv(ncomm,"PETSC_VIEWER_HDF5_FILENAME",fname,PETSC_MAX_PATH_LEN,&flg); 16312cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1632a75e6a4aSMatthew G. Knepley if (!flg) { 1633a75e6a4aSMatthew G. Knepley ierr = PetscStrcpy(fname,"output.h5"); 16342cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1635a75e6a4aSMatthew G. Knepley } 1636a75e6a4aSMatthew G. Knepley ierr = PetscViewerHDF5Open(ncomm,fname,FILE_MODE_WRITE,&viewer); 16372cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1638a75e6a4aSMatthew G. Knepley ierr = PetscObjectRegisterDestroy((PetscObject)viewer); 16392cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 164047435625SJed Brown ierr = MPI_Comm_set_attr(ncomm,Petsc_Viewer_HDF5_keyval,(void*)viewer); 1641908793a3SLisandro Dalcin if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_INITIAL," ");PetscFunctionReturn(NULL);} 1642a75e6a4aSMatthew G. Knepley } 1643a75e6a4aSMatthew G. Knepley ierr = PetscCommDestroy(&ncomm); 16442cb5e1ccSBarry Smith if (ierr) {PetscError(PETSC_COMM_SELF,__LINE__,"PETSC_VIEWER_HDF5_",__FILE__,PETSC_ERR_PLIB,PETSC_ERROR_REPEAT," ");PetscFunctionReturn(NULL);} 1645a75e6a4aSMatthew G. Knepley PetscFunctionReturn(viewer); 1646a75e6a4aSMatthew G. Knepley } 1647