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; }
73template<
class NumEqVector,
typename std::enable_if_t<NumEqVector::size() == 3,
int> = 0>
86 template<
class Problem,
class FVElementGeometry,
class ElementVolumeVariables>
88 const typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity& element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const typename FVElementGeometry::SubControlVolumeFace& scvf)
93 using Scalar =
typename NumEqVector::value_type;
94 using FluidSystem =
typename ElementVolumeVariables::VolumeVariables::FluidSystem;
99 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
100 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
102 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
103 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
105 const auto gradVelocity = [&]()
108 const auto velocityXLeft = insideVolVars.velocity(0);
109 const auto velocityYLeft = insideVolVars.velocity(1);
110 const auto velocityXRight = outsideVolVars.velocity(0);
111 const auto velocityYRight = outsideVolVars.velocity(1);
115 const auto& cellCenterToCellCenter = outsideScv.center() - insideScv.center();
116 const auto distance = cellCenterToCellCenter.two_norm();
117 const auto&
unitNormal = scvf.unitOuterNormal();
119 return std::array<Scalar, 2>{
120 (velocityXRight-velocityXLeft)*direction/
distance,
121 (velocityYRight-velocityYLeft)*direction/
distance
126 const auto waterDepthLeft = insideVolVars.waterDepth();
127 const auto waterDepthRight = outsideVolVars.waterDepth();
128 const auto averageDepth = 2.0*(waterDepthLeft*waterDepthRight)/(waterDepthLeft + waterDepthRight);
131 const Scalar turbViscosity = [&]()
134 static const auto turbBGViscosity = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.TurbulentViscosity", 1.0e-6);
137 static const auto useMixingLengthTurbulenceModel = getParamFromGroup<bool>(problem.paramGroup(),
"ShallowWater.UseMixingLengthTurbulenceModel",
false);
140 if (!useMixingLengthTurbulenceModel)
141 return turbBGViscosity;
143 using SpatialParams =
typename Problem::SpatialParams;
144 using E =
typename FVElementGeometry::GridGeometry::GridView::template Codim<0>::Entity;
145 using SCV =
typename FVElementGeometry::SubControlVolume;
146 if constexpr (!Detail::implementsFrictionLaw<SpatialParams, E, SCV>())
147 DUNE_THROW(Dune::IOError,
"Mixing length turbulence model enabled but spatial parameters do not implement the frictionLaw interface!");
153 static const auto turbConstV = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.VerticalCoefficientOfMixingLengthModel", 1.0);
154 static const auto turbConstH = getParamFromGroup<Scalar>(problem.paramGroup(),
"ShallowWater.HorizontalCoefficientOfMixingLengthModel", 0.1);
162 constexpr Scalar kappa = 0.41;
164 const Scalar ustar = [&]()
168 const auto& bottomShearStressInside = problem.spatialParams().frictionLaw(element, insideScv).bottomShearStress(insideVolVars);
169 const auto& bottomShearStressOutside = problem.spatialParams().frictionLaw(element, outsideScv).bottomShearStress(outsideVolVars);
170 const auto bottomShearStressInsideMag = bottomShearStressInside.two_norm();
171 const auto bottomShearStressOutsideMag = bottomShearStressOutside.two_norm();
175 const auto averageBottomShearStress = 2.0*(bottomShearStressInsideMag*bottomShearStressOutsideMag)
176 / max(1.0e-8,bottomShearStressInsideMag + bottomShearStressOutsideMag);
180 static_assert(!FluidSystem::isCompressible(0),
181 "The shallow water model assumes incompressible fluids"
185 return sqrt(averageBottomShearStress/insideVolVars.density());
188 const auto turbViscosityV = turbConstV * (kappa/6.0) * ustar * averageDepth;
215 const auto [gradU, gradV] = gradVelocity;
216 const auto mixingLengthSquared = turbConstH * turbConstH * averageDepth * averageDepth;
217 const auto turbViscosityH = mixingLengthSquared * sqrt(2.0*gradU*gradU + 2.0*gradV*gradV);
220 return turbBGViscosity + sqrt(turbViscosityV*turbViscosityV + turbViscosityH*turbViscosityH);
227 const auto freeSurfaceInside = insideVolVars.waterDepth() + insideVolVars.bedSurface();
228 const auto freeSurfaceOutside = outsideVolVars.waterDepth() + outsideVolVars.bedSurface();
229 const auto interfaceWaterDepth = max(min(freeSurfaceInside , freeSurfaceOutside) - max(insideVolVars.bedSurface(),outsideVolVars.bedSurface()),0.0);
230 const auto& [gradU, gradV] = gradVelocity;
231 const auto uViscousFlux = turbViscosity * interfaceWaterDepth * gradU;
232 const auto vViscousFlux = turbViscosity * interfaceWaterDepth * gradV;
235 localFlux[1] = -uViscousFlux * scvf.area();
236 localFlux[2] = -vViscousFlux * scvf.area();
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Function to limit the fluxes.
Classes related to flux variables caching.
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:294
Vector unitNormal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:70
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:87
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
The empty filler class corresponding to EmptyCache.
Definition: fluxvariablescaching.hh:32
Definition: fluxvariablescaching.hh:67
Compute the shallow water viscous momentum flux due to (turbulent) viscosity.
Definition: shallowwaterviscousflux.hh:75