24#ifndef DUMUX_FLUX_SHALLOW_WATER_VISCOUS_FLUX_HH
25#define DUMUX_FLUX_SHALLOW_WATER_VISCOUS_FLUX_HH
33#include <dune/common/std/type_traits.hh>
34#include <dune/common/exceptions.hh>
45template <
typename T,
typename ...Ts>
46using FrictionLawDetector =
decltype(std::declval<T>().frictionLaw(std::declval<Ts>()...));
48template<
class T,
typename ...Args>
49static constexpr bool implementsFrictionLaw()
50{
return Dune::Std::is_detected<FrictionLawDetector, T, Args...>::value; }
61 typename std::enable_if_t<NumEqVector::size() == 3,
int> = 0>
83 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables>
85 const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
86 const FVElementGeometry& fvGeometry,
87 const ElementVolumeVariables& elemVolVars,
88 const typename FVElementGeometry::SubControlVolumeFace& scvf)
90 using Scalar =
typename PrimaryVariables::value_type;
95 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
96 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
98 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
99 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
101 const auto gradVelocity = [&]()
104 const auto velocityXLeft = insideVolVars.velocity(0);
105 const auto velocityYLeft = insideVolVars.velocity(1);
106 const auto velocityXRight = outsideVolVars.velocity(0);
107 const auto velocityYRight = outsideVolVars.velocity(1);
111 const auto& cellCenterToCellCenter = outsideScv.center() - insideScv.center();
112 const auto distance = cellCenterToCellCenter.two_norm();
113 const auto&
unitNormal = scvf.unitOuterNormal();
115 return std::array<Scalar, 2>{
116 (velocityXRight-velocityXLeft)*direction/
distance,
117 (velocityYRight-velocityYLeft)*direction/
distance
122 const auto waterDepthLeft = insideVolVars.waterDepth();
123 const auto waterDepthRight = outsideVolVars.waterDepth();
124 const auto averageDepth = 2.0*(waterDepthLeft*waterDepthRight)/(waterDepthLeft + waterDepthRight);
127 const Scalar turbViscosity = [&]()
130 static const auto turbBGViscosity = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.TurbulentViscosity", 1.0e-6);
133 static const auto useMixingLengthTurbulenceModel = getParamFromGroup<bool>(problem.paramGroup(),
"ShallowWater.UseMixingLengthTurbulenceModel",
false);
136 if (!useMixingLengthTurbulenceModel)
137 return turbBGViscosity;
139 using SpatialParams =
typename Problem::SpatialParams;
140 using E =
typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity;
141 using SCV =
typename FVElementGeometry::SubControlVolume;
142 if constexpr (!Detail::implementsFrictionLaw<SpatialParams, E, SCV>())
143 DUNE_THROW(Dune::IOError,
"Mixing length turbulence model enabled but spatial parameters do not implement the frictionLaw interface!");
149 static const auto turbConstV = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.VerticalCoefficientOfMixingLengthModel", 1.0);
150 static const auto turbConstH = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.HorizontalCoefficientOfMixingLengthModel", 0.1);
158 constexpr Scalar kappa = 0.41;
160 const Scalar ustar = [&]()
164 const auto& bottomShearStressInside = problem.spatialParams().frictionLaw(element, insideScv).shearStress(insideVolVars);
165 const auto& bottomShearStressOutside = problem.spatialParams().frictionLaw(element, outsideScv).shearStress(outsideVolVars);
166 const auto bottomShearStressInsideMag = bottomShearStressInside.two_norm();
167 const auto bottomShearStressOutsideMag = bottomShearStressOutside.two_norm();
171 const auto averageBottomShearStress = 2.0*(bottomShearStressInsideMag*bottomShearStressOutsideMag)
172 / max(1.0e-8,bottomShearStressInsideMag + bottomShearStressOutsideMag);
176 return sqrt(averageBottomShearStress);
179 const auto turbViscosityV = turbConstV * (kappa/6.0) * ustar * averageDepth;
206 const auto [gradU, gradV] = gradVelocity;
207 const auto mixingLengthSquared = turbConstH * turbConstH * averageDepth * averageDepth;
208 const auto turbViscosityH = mixingLengthSquared * sqrt(2.0*gradU*gradU + 2.0*gradV*gradV);
211 return turbBGViscosity + sqrt(turbViscosityV*turbViscosityV + turbViscosityH*turbViscosityH);
216 const auto [gradU, gradV] = gradVelocity;
217 const auto uViscousFlux = turbViscosity * averageDepth * gradU;
218 const auto vViscousFlux = turbViscosity * averageDepth * gradV;
221 static const auto upperWaterDepthFluxLimiting = getParamFromGroup<double>(problem.paramGroup(),
"FluxLimiterLET.UpperWaterDepth", 1e-3);
222 static const auto lowerWaterDepthFluxLimiting = getParamFromGroup<double>(problem.paramGroup(),
"FluxLimiterLET.LowerWaterDepth", 1e-5);
224 const auto limitingDepth = (waterDepthLeft + waterDepthRight) * 0.5;
227 upperWaterDepthFluxLimiting,
228 lowerWaterDepthFluxLimiting);
231 localFlux[1] = -uViscousFlux *
mobility * scvf.area();
232 localFlux[2] = -vViscousFlux *
mobility * scvf.area();
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Function to limit the fluxes.
Classes related to flux variables caching.
Vector unitNormal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:68
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:138
static NumEqVector flux(const Problem &problem, const typename FVElementGeometry::GridGeometry::GridView::template Codim< 0 >::Entity &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const typename FVElementGeometry::SubControlVolumeFace &scvf)
Compute the viscous momentum flux contribution from the interface shear stress.
Definition: shallowwaterviscousflux.hh:84
static Scalar fluxLimiterLET(const Scalar valueLeft, const Scalar valueRight, const Scalar upperH, const Scalar lowerH)
Flux limiter function to scale fluxes for small water depths.
Definition: fluxlimiterlet.hh:52
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
std::string mobility(int phaseIdx) noexcept
I/O name of mobility for multiphase systems.
Definition: name.hh:101
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:67
Computes the shallow water viscous momentum flux due to (turbulent) viscosity by adding all surroundi...
Definition: shallowwaterviscousflux.hh:63