From 9663689f8e5502975d17c2fd9ba189b073ea60c5 Mon Sep 17 00:00:00 2001 From: Timo Koch <timo.koch@iws.uni-stuttgart.de> Date: Wed, 10 Jun 2020 17:17:24 +0200 Subject: [PATCH] Add classes for face-related frontal/lateral scvfs and scvs --- dumux/assembly/staggeredlocalresidual.hh | 30 ++++--- .../freeflow/fvgridgeometrytraits.hh | 60 ++++++++++++++ .../staggered/fvgridgeometry.hh | 58 +++++++------- .../navierstokes/staggered/fluxvariables.hh | 78 +++++++++++++------ .../test_staggered_free_flow_geometry.cc | 26 +------ .../staggered/test_staggeredfvgeometry.cc | 9 +++ 6 files changed, 177 insertions(+), 84 deletions(-) diff --git a/dumux/assembly/staggeredlocalresidual.hh b/dumux/assembly/staggeredlocalresidual.hh index 15ef2a4a9a..a9933129c7 100644 --- a/dumux/assembly/staggeredlocalresidual.hh +++ b/dumux/assembly/staggeredlocalresidual.hh @@ -47,9 +47,11 @@ class StaggeredLocalResidual using ElementBoundaryTypes = GetPropType<TypeTag, Properties::ElementBoundaryTypes>; using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView; - using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView; - using SubControlVolume = typename FVElementGeometry::SubControlVolume; - using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; + using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; + using FVElementGeometry = typename GridGeometry::LocalView; + using SubControlVolume = typename GridGeometry::SubControlVolume; + using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; + using FaceSubControlVolume = typename GridGeometry::Traits::FaceSubControlVolume; using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>; using CellCenterResidual = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>; @@ -241,8 +243,13 @@ public: const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); const auto extrusionFactor = elemVolVars[scv].extrusionFactor(); + // contruct staggered scv (half of the element) + auto scvCenter = scvf.center() - scv.center(); + scvCenter *= 0.5; + FaceSubControlVolume faceScv(scvCenter, 0.5*scv.volume()); + // multiply by 0.5 because we only consider half of a staggered control volume here - source *= 0.5*scv.volume()*extrusionFactor; + source *= faceScv.volume()*extrusionFactor; residual -= source; } @@ -272,18 +279,19 @@ public: const ElementFaceVariables& curElemFaceVars, const SubControlVolumeFace& scvf) const { - FaceResidualValue storage(0.0); const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); - auto prevFaceStorage = asImp_().computeStorageForFace(problem, scvf, prevElemVolVars[scv], prevElemFaceVars); - auto curFaceStorage = asImp_().computeStorageForFace(problem, scvf, curElemVolVars[scv], curElemFaceVars); - storage = std::move(curFaceStorage); - storage -= std::move(prevFaceStorage); + auto storage = asImp_().computeStorageForFace(problem, scvf, curElemVolVars[scv], curElemFaceVars); + storage -= asImp_().computeStorageForFace(problem, scvf, prevElemVolVars[scv], prevElemFaceVars); const auto extrusionFactor = curElemVolVars[scv].extrusionFactor(); - // multiply by 0.5 because we only consider half of a staggered control volume here - storage *= 0.5*scv.volume()*extrusionFactor; + // contruct staggered scv (half of the element) + auto scvCenter = scvf.center() - scv.center(); + scvCenter *= 0.5; + FaceSubControlVolume faceScv(scvCenter, 0.5*scv.volume()); + + storage *= faceScv.volume()*extrusionFactor; storage /= timeLoop_->timeStepSize(); residual += storage; diff --git a/dumux/discretization/staggered/freeflow/fvgridgeometrytraits.hh b/dumux/discretization/staggered/freeflow/fvgridgeometrytraits.hh index 051ca1c413..eb40284cd8 100644 --- a/dumux/discretization/staggered/freeflow/fvgridgeometrytraits.hh +++ b/dumux/discretization/staggered/freeflow/fvgridgeometrytraits.hh @@ -35,6 +35,57 @@ namespace Dumux { +template<class GridView> +class FreeflowStaggeredSCV +{ + using ctype = typename GridView::ctype; + using GlobalPosition = Dune::FieldVector<ctype, GridView::dimensionworld>; +public: + struct Traits + { + using Scalar = ctype; + using Geometry = typename GridView::template Codim<0>::Geometry; + }; + + FreeflowStaggeredSCV(const GlobalPosition& center, const ctype volume) + : center_(center), volume_(volume) {} + + const GlobalPosition& center() const + { return center_; } + + ctype volume() const + { return volume_; } +private: + GlobalPosition center_; + ctype volume_; +}; + +template<class GridView> +class FreeflowStaggeredSCVF +{ + using ctype = typename GridView::ctype; + using GlobalPosition = Dune::FieldVector<ctype, GridView::dimensionworld>; +public: + struct Traits + { + using Scalar = ctype; + using Geometry = typename GridView::template Codim<1>::Geometry; + }; + + FreeflowStaggeredSCVF(const GlobalPosition& center, const ctype area) + : center_(center), area_(area) {} + + const GlobalPosition& center() const + { return center_; } + + ctype area() const + { return area_; } +private: + GlobalPosition center_; + ctype area_; +}; + + /*! * \ingroup StaggeredDiscretization * \brief Default traits for the finite volume grid geometry. @@ -60,6 +111,15 @@ struct StaggeredFreeFlowDefaultFVGridGeometryTraits template<class GridGeometry, bool cachingEnabled> using LocalView = StaggeredFVElementGeometry<GridGeometry, cachingEnabled>; + + struct PublicTraits + { + using CellSubControlVolume = SubControlVolume; + using CellSubControlVolumeFace = SubControlVolumeFace; + using FaceSubControlVolume = FreeflowStaggeredSCV<GridView>; + using FaceLateralSubControlVolumeFace = FreeflowStaggeredSCVF<GridView>; + using FaceFrontalSubControlVolumeFace = FreeflowStaggeredSCVF<GridView>; + }; }; } //end namespace Dumux diff --git a/dumux/discretization/staggered/fvgridgeometry.hh b/dumux/discretization/staggered/fvgridgeometry.hh index feca485740..e6119c7c43 100644 --- a/dumux/discretization/staggered/fvgridgeometry.hh +++ b/dumux/discretization/staggered/fvgridgeometry.hh @@ -174,37 +174,40 @@ class StaggeredFVGridGeometry; * This builds up the sub control volumes and sub control volume faces * for each element. Specialization in case the FVElementGeometries are stored. */ -template<class GV, class Traits> -class StaggeredFVGridGeometry<GV, true, Traits> -: public BaseGridGeometry<GV, Traits> +template<class GV, class T> +class StaggeredFVGridGeometry<GV, true, T> +: public BaseGridGeometry<GV, T> { - using ThisType = StaggeredFVGridGeometry<GV, true, Traits>; - using ParentType = BaseGridGeometry<GV, Traits>; + using ThisType = StaggeredFVGridGeometry<GV, true, T>; + using ParentType = BaseGridGeometry<GV, T>; using GridIndexType = typename IndexTraits<GV>::GridIndex; using LocalIndexType = typename IndexTraits<GV>::LocalIndex; using Element = typename GV::template Codim<0>::Entity; - using IntersectionMapper = typename Traits::IntersectionMapper; - using GeometryHelper = typename Traits::GeometryHelper; - using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>; + using IntersectionMapper = typename T::IntersectionMapper; + using GeometryHelper = typename T::GeometryHelper; + using ConnectivityMap = typename T::template ConnectivityMap<ThisType>; public: + //! export the traits + using Traits = typename T::PublicTraits; + //! export discretization method static constexpr DiscretizationMethod discMethod = DiscretizationMethod::staggered; - static constexpr int upwindSchemeOrder = Traits::upwindSchemeOrder; + static constexpr int upwindSchemeOrder = T::upwindSchemeOrder; static constexpr bool useHigherOrder = upwindSchemeOrder > 1; static constexpr bool cachingEnabled = true; //! export the type of the fv element geometry (the local view type) - using LocalView = typename Traits::template LocalView<ThisType, true>; + using LocalView = typename T::template LocalView<ThisType, true>; //! export the type of sub control volume - using SubControlVolume = typename Traits::SubControlVolume; + using SubControlVolume = typename T::SubControlVolume; //! export the type of sub control volume - using SubControlVolumeFace = typename Traits::SubControlVolumeFace; + using SubControlVolumeFace = typename T::SubControlVolumeFace; //! export the grid view type using GridView = GV; //! export the dof type indices - using DofTypeIndices = typename Traits::DofTypeIndices; + using DofTypeIndices = typename T::DofTypeIndices; //! return a integral constant for cell center dofs static constexpr auto cellCenterIdx() @@ -433,38 +436,41 @@ private: * This builds up the sub control volumes and sub control volume faces * for each element. Specialization in case the FVElementGeometries are stored. */ -template<class GV, class Traits> -class StaggeredFVGridGeometry<GV, false, Traits> -: public BaseGridGeometry<GV, Traits> +template<class GV, class T> +class StaggeredFVGridGeometry<GV, false, T> +: public BaseGridGeometry<GV, T> { - using ThisType = StaggeredFVGridGeometry<GV, false, Traits>; - using ParentType = BaseGridGeometry<GV, Traits>; + using ThisType = StaggeredFVGridGeometry<GV, false, T>; + using ParentType = BaseGridGeometry<GV, T>; using GridIndexType = typename IndexTraits<GV>::GridIndex; using LocalIndexType = typename IndexTraits<GV>::LocalIndex; using Element = typename GV::template Codim<0>::Entity; - using IntersectionMapper = typename Traits::IntersectionMapper; - using ConnectivityMap = typename Traits::template ConnectivityMap<ThisType>; + using IntersectionMapper = typename T::IntersectionMapper; + using ConnectivityMap = typename T::template ConnectivityMap<ThisType>; public: + //! export the traits + using Traits = typename T::PublicTraits; + //! export discretization method static constexpr DiscretizationMethod discMethod = DiscretizationMethod::staggered; - static constexpr int upwindSchemeOrder = Traits::upwindSchemeOrder; + static constexpr int upwindSchemeOrder = T::upwindSchemeOrder; static constexpr bool useHigherOrder = upwindSchemeOrder > 1; static constexpr bool cachingEnabled = false; - using GeometryHelper = typename Traits::GeometryHelper; + using GeometryHelper = typename T::GeometryHelper; //! export the type of the fv element geometry (the local view type) - using LocalView = typename Traits::template LocalView<ThisType, false>; + using LocalView = typename T::template LocalView<ThisType, false>; //! export the type of sub control volume - using SubControlVolume = typename Traits::SubControlVolume; + using SubControlVolume = typename T::SubControlVolume; //! export the type of sub control volume - using SubControlVolumeFace = typename Traits::SubControlVolumeFace; + using SubControlVolumeFace = typename T::SubControlVolumeFace; //! export the grid view type using GridView = GV; //! export the dof type indices - using DofTypeIndices = typename Traits::DofTypeIndices; + using DofTypeIndices = typename T::DofTypeIndices; //! return a integral constant for cell center dofs static constexpr auto cellCenterIdx() diff --git a/dumux/freeflow/navierstokes/staggered/fluxvariables.hh b/dumux/freeflow/navierstokes/staggered/fluxvariables.hh index 6c8ef859a6..4363a3269a 100644 --- a/dumux/freeflow/navierstokes/staggered/fluxvariables.hh +++ b/dumux/freeflow/navierstokes/staggered/fluxvariables.hh @@ -78,6 +78,8 @@ class NavierStokesFluxVariablesImpl<TypeTag, DiscretizationMethod::staggered> using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; using Indices = typename ModelTraits::Indices; using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; + using FaceFrontalSubControlVolumeFace = typename GridGeometry::Traits::FaceFrontalSubControlVolumeFace; + using FaceLateralSubControlVolumeFace = typename GridGeometry::Traits::FaceLateralSubControlVolumeFace; using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>; using FacePrimaryVariables = GetPropType<TypeTag, Properties::FacePrimaryVariables>; using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>; @@ -246,8 +248,10 @@ public: frontalFlux += pressure * -1.0 * scvf.directionSign(); // Account for the staggered face's area. For rectangular elements, this equals the area of the scvf - // our velocity dof of interest lives on. - return frontalFlux * scvf.area() * insideVolVars.extrusionFactor(); + // our velocity dof of interest lives on but with adjusted centroid + const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); + FaceFrontalSubControlVolumeFace frontalFace(scv.center(), scvf.area()); + return frontalFlux * frontalFace.area() * insideVolVars.extrusionFactor(); } /*! @@ -290,14 +294,14 @@ public: { const auto eIdx = scvf.insideScvIdx(); // Get the face normal to the face the dof lives on. The staggered sub face conincides with half of this lateral face. - const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx); + const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx); // Create a boundaryTypes object (will be empty if not at a boundary). std::optional<BoundaryTypes> lateralFaceBoundaryTypes; // Check if there is face/element parallel to our face of interest where the dof lives on. If there is no parallel neighbor, // we are on a boundary where we have to check for boundary conditions. - if (lateralScvf.boundary()) + if (lateralFace.boundary()) { // Retrieve the boundary types that correspond to the center of the lateral scvf. As a convention, we always query // the type of BCs at the center of the element's "actual" lateral scvf (not the face of the staggered control volume). @@ -310,7 +314,7 @@ public: // | || | (center of lateral scvf) // ---------------- V position at which the value of the boundary conditions will be evaluated // (center of the staggered lateral face) - lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralScvf)); + lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralFace)); } // If the current scvf is on a bounary and if there is a Neumann or Beavers-Joseph-(Saffmann) BC for the stress in tangential direction, @@ -319,14 +323,15 @@ public: if (currentScvfBoundaryTypes) { // Handle Neumann BCs. - if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralScvf.directionIndex()))) + if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralFace.directionIndex()))) { // Get the location of the lateral staggered face's center. + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx); lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, - scvf.makeBoundaryFace(lateralStaggeredFaceCenter))[Indices::velocity(lateralScvf.directionIndex())] - * extrusionFactor_(elemVolVars, lateralScvf) * lateralScvf.area() * 0.5 - * lateralScvf.directionSign(); + scvf.makeBoundaryFace(lateralStaggeredFaceCenter))[Indices::velocity(lateralFace.directionIndex())] + * extrusionFactor_(elemVolVars, lateralFace) * lateralScvf.area() + * lateralFace.directionSign(); continue; } @@ -334,21 +339,22 @@ public: // It is not clear how to evaluate the BJ condition here. // For symmetry reasons, our own scvf should then have the same Neumann flux as the lateral face. // TODO: We should clarify if this is the correct approach. - if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralScvf.directionIndex())) && lateralFaceBoundaryTypes && + if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())) && lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex()))) { + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx); lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, - lateralScvf.makeBoundaryFace(lateralStaggeredFaceCenter))[Indices::velocity(scvf.directionIndex())] - * extrusionFactor_(elemVolVars, lateralScvf) * lateralScvf.area() * 0.5 - * lateralScvf.directionSign(); + lateralFace.makeBoundaryFace(lateralStaggeredFaceCenter))[Indices::velocity(scvf.directionIndex())] + * extrusionFactor_(elemVolVars, lateralFace) * lateralScvf.area() + * lateralFace.directionSign(); continue; } } // Check if the lateral face (perpendicular to our current scvf) lies on a boundary. If yes, boundary conditions might need to be treated // and further calculations can be skipped. - if (lateralScvf.boundary()) + if (lateralFace.boundary()) { // Check if we have a symmetry boundary condition. If yes, the tangential part of the momentum flux can be neglected // and we may skip any further calculations for the given sub face. @@ -360,10 +366,11 @@ public: { // Construct a temporary scvf which corresponds to the staggered sub face, featuring the location // the staggered faces's center. + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx); - const auto lateralBoundaryFace = lateralScvf.makeBoundaryFace(lateralStaggeredFaceCenter); + const auto lateralBoundaryFace = lateralFace.makeBoundaryFace(lateralStaggeredFaceCenter); lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())] - * elemVolVars[lateralScvf.insideScvIdx()].extrusionFactor() * lateralScvf.area() * 0.5; + * elemVolVars[lateralFace.insideScvIdx()].extrusionFactor() * lateralScvf.area(); continue; } @@ -518,8 +525,9 @@ private: const bool selfIsUpstream = lateralFace.directionSign() == sign(transportingVelocity); StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods()); + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); return upwindHelper.computeUpwindLateralMomentum(selfIsUpstream, lateralFace, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes) - * transportingVelocity * lateralFace.directionSign() * lateralFace.area() * 0.5 * extrusionFactor_(elemVolVars, lateralFace); + * transportingVelocity * lateralFace.directionSign() * lateralScvf.area() * extrusionFactor_(elemVolVars, lateralFace); } /*! @@ -596,7 +604,8 @@ private: } // Account for the area of the staggered lateral face (0.5 of the coinciding scfv). - return lateralDiffusiveFlux * lateralFace.area() * 0.5 * extrusionFactor_(elemVolVars, lateralFace); + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); + return lateralDiffusiveFlux * lateralScvf.area() * extrusionFactor_(elemVolVars, lateralFace); } /*! @@ -618,6 +627,29 @@ private: return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter; }; + /*! + * \brief Get the location of the lateral staggered sub control volume face center. + * Only needed for boundary conditions if the current scvf or the lateral one is on a bounary. + * + * \verbatim + * --------###o#### || frontal face of staggered half-control-volume + * | || | current scvf # lateral staggered face of interest (may lie on a boundary) + * | || | x dof position + * | || x~~~~> vel.Self -- element boundaries, current scvf may lie on a boundary + * | || | o center of the lateral staggered scvf + * | || | + * ---------------- + * \endverbatim + */ + GlobalPosition lateralStaggeredSCVFCenter_(const SubControlVolumeFace& lateralFace, + const SubControlVolumeFace& currentFace, + const int localSubFaceIdx) const + { + auto pos = currentFace.pairData(localSubFaceIdx).lateralStaggeredFaceCenter - lateralFace.center(); + pos *= 0.5; + return pos; + }; + //! helper function to get the averaged extrusion factor for a face static Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf) { @@ -643,14 +675,14 @@ private: const std::size_t localSubFaceIdx) const { const auto eIdx = scvf.insideScvIdx(); - const auto& lateralScvf = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx); + const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx); - if (problem.useWallFunction(element, lateralScvf, Indices::velocity(scvf.directionIndex()))) + if (problem.useWallFunction(element, lateralFace, Indices::velocity(scvf.directionIndex()))) { - const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx); - const auto lateralBoundaryFace = lateralScvf.makeBoundaryFace(lateralStaggeredFaceCenter); + FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area()); + const auto lateralBoundaryFace = lateralFace.makeBoundaryFace(lateralStaggeredFaceCenter_(scvf, localSubFaceIdx)); lateralFlux += problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())] - * extrusionFactor_(elemVolVars, lateralScvf) * lateralScvf.area() * 0.5; + * extrusionFactor_(elemVolVars, lateralFace) * lateralScvf.area(); return true; } else diff --git a/test/discretization/staggered/test_staggered_free_flow_geometry.cc b/test/discretization/staggered/test_staggered_free_flow_geometry.cc index b1cdde06fc..5df69efaa5 100644 --- a/test/discretization/staggered/test_staggered_free_flow_geometry.cc +++ b/test/discretization/staggered/test_staggered_free_flow_geometry.cc @@ -36,6 +36,7 @@ #include <dumux/common/defaultmappertraits.hh> #include <dumux/discretization/cellcentered/subcontrolvolume.hh> #include <dumux/discretization/staggered/fvelementgeometry.hh> +#include <dumux/discretization/staggered/freeflow/fvgridgeometrytraits.hh> #include <dumux/discretization/staggered/fvgridgeometry.hh> #include <dumux/discretization/staggered/freeflow/subcontrolvolumeface.hh> #include <dumux/discretization/staggered/freeflow/staggeredgeometryhelper.hh> @@ -52,29 +53,6 @@ public: }; } // end namespace Detail -//! the fv grid geometry traits for this test -template<class GridView, int upwOrder> -struct TestFVGGTraits : public DefaultMapperTraits<GridView> -{ - using SubControlVolume = CCSubControlVolume<GridView>; - using SubControlVolumeFace = FreeFlowStaggeredSubControlVolumeFace<GridView, upwOrder>; - using IntersectionMapper = ConformingGridIntersectionMapper<GridView>; - using GeometryHelper = FreeFlowStaggeredGeometryHelper<GridView, upwOrder>; - static constexpr int upwindSchemeOrder = upwOrder; - - struct DofTypeIndices - { - using CellCenterIdx = Dune::index_constant<0>; - using FaceIdx = Dune::index_constant<1>; - }; - - template<class GridGeometry> - using ConnectivityMap = StaggeredFreeFlowConnectivityMap<GridGeometry>; - - template<class GridGeometry, bool cachingEnabled> - using LocalView = StaggeredFVElementGeometry<GridGeometry, cachingEnabled>; -}; - } // end namespace Dumux #endif @@ -96,7 +74,7 @@ int main (int argc, char *argv[]) try static constexpr int upwindSchemeOrder = 2; using GridGeometry = StaggeredFVGridGeometry<typename Grid::LeafGridView, /*enable caching=*/ true, - TestFVGGTraits<typename Grid::LeafGridView, upwindSchemeOrder> >; + StaggeredFreeFlowDefaultFVGridGeometryTraits<typename Grid::LeafGridView, upwindSchemeOrder> >; using FVElementGeometry = typename GridGeometry::LocalView; using SubControlVolume = typename FVElementGeometry::SubControlVolume; using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; diff --git a/test/discretization/staggered/test_staggeredfvgeometry.cc b/test/discretization/staggered/test_staggeredfvgeometry.cc index e8a067704c..cfbb08cb42 100644 --- a/test/discretization/staggered/test_staggeredfvgeometry.cc +++ b/test/discretization/staggered/test_staggeredfvgeometry.cc @@ -78,6 +78,15 @@ struct TestFVGGTraits : public DefaultMapperTraits<GridView> template<class GridGeometry, bool enableCache> using LocalView = StaggeredFVElementGeometry<GridGeometry, enableCache>; + + struct PublicTraits + { + using CellSubControlVolume = SubControlVolume; + using CellSubControlVolumeFace = SubControlVolumeFace; + using FaceSubControlVolume = SubControlVolume; + using FaceLateralSubControlVolumeFace = SubControlVolumeFace; + using FaceFrontalSubControlVolumeFace = SubControlVolumeFace; + }; }; } // end namespace Dumux -- GitLab