3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
shallowwaterviscousflux.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_FLUX_SHALLOW_WATER_VISCOUS_FLUX_HH
25#define DUMUX_FLUX_SHALLOW_WATER_VISCOUS_FLUX_HH
26
27#include <cmath>
28#include <algorithm>
29#include <utility>
30#include <type_traits>
31#include <array>
32
33#include <dune/common/std/type_traits.hh>
34#include <dune/common/exceptions.hh>
35
39
40namespace Dumux {
41
42#ifndef DOXYGEN
43namespace Detail {
44// helper struct detecting if the user-defined spatial params class has a frictionLaw function
45template <typename T, typename ...Ts>
46using FrictionLawDetector = decltype(std::declval<T>().frictionLaw(std::declval<Ts>()...));
47
48template<class T, typename ...Args>
49static constexpr bool implementsFrictionLaw()
50{ return Dune::Std::is_detected<FrictionLawDetector, T, Args...>::value; }
51} // end namespace Detail
52#endif
53
60template<class PrimaryVariables, class NumEqVector,
61 typename std::enable_if_t<NumEqVector::size() == 3, int> = 0>
63{
64
65public:
66
83 template<class Problem, class FVElementGeometry, class ElementVolumeVariables>
84 static NumEqVector flux(const Problem& problem,
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)
89 {
90 using Scalar = typename PrimaryVariables::value_type;
91
92 NumEqVector localFlux(0.0);
93
94 // Get the inside and outside volume variables
95 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
96 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
97
98 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
99 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
100
101 const auto gradVelocity = [&]()
102 {
103 // The left (inside) and right (outside) states
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);
108
109 // Compute the velocity gradients normal to the face
110 // Factor that takes the direction of the unit vector into account
111 const auto& cellCenterToCellCenter = outsideScv.center() - insideScv.center();
112 const auto distance = cellCenterToCellCenter.two_norm();
113 const auto& unitNormal = scvf.unitOuterNormal();
114 const auto direction = (unitNormal*cellCenterToCellCenter)/distance;
115 return std::array<Scalar, 2>{
116 (velocityXRight-velocityXLeft)*direction/distance,
117 (velocityYRight-velocityYLeft)*direction/distance
118 };
119 }();
120
121 // Use a harmonic average of the depth at the interface.
122 const auto waterDepthLeft = insideVolVars.waterDepth();
123 const auto waterDepthRight = outsideVolVars.waterDepth();
124 const auto averageDepth = 2.0*(waterDepthLeft*waterDepthRight)/(waterDepthLeft + waterDepthRight);
125
126 // compute the turbulent viscosity contribution
127 const Scalar turbViscosity = [&]()
128 {
129 // The (constant) background turbulent viscosity
130 static const auto turbBGViscosity = getParamFromGroup<Scalar>(problem.paramGroup(), "ShallowWater.TurbulentViscosity", 1.0e-6);
131
132 // Check whether the mixing-length turbulence model is used
133 static const auto useMixingLengthTurbulenceModel = getParamFromGroup<bool>(problem.paramGroup(), "ShallowWater.UseMixingLengthTurbulenceModel", false);
134
135 // constant eddy viscosity equal to the prescribed background eddy viscosity
136 if (!useMixingLengthTurbulenceModel)
137 return turbBGViscosity;
138
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!");
144 else
145 {
146 // turbulence model based on mixing length
147 // Compute the turbulent viscosity using a combined horizonal/vertical mixing length approach
148 // Turbulence coefficients: vertical (Elder like) and horizontal (Smagorinsky like)
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);
151
158 constexpr Scalar kappa = 0.41;
159 // Compute the average shear velocity on the face
160 const Scalar ustar = [&]()
161 {
162 // Get the bottom shear stress in the two adjacent cells
163 // Note that element is not needed in spatialParams().frictionLaw (should be removed). For now we simply pass the same element twice
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();
168
169 // Use a harmonic average of the viscosity and the depth at the interface.
170 using std::max;
171 const auto averageBottomShearStress = 2.0*(bottomShearStressInsideMag*bottomShearStressOutsideMag)
172 / max(1.0e-8,bottomShearStressInsideMag + bottomShearStressOutsideMag);
173
174 // Shear velocity possibly needed in mixing-length turbulence model in the computation of the diffusion flux
175 using std::sqrt;
176 return sqrt(averageBottomShearStress);
177 }();
178
179 const auto turbViscosityV = turbConstV * (kappa/6.0) * ustar * averageDepth;
180
205 using std::sqrt;
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);
209
210 // Total turbulent viscosity
211 return turbBGViscosity + sqrt(turbViscosityV*turbViscosityV + turbViscosityH*turbViscosityH);
212 }
213 }();
214
215 // Compute the viscous momentum fluxes
216 const auto [gradU, gradV] = gradVelocity;
217 const auto uViscousFlux = turbViscosity * averageDepth * gradU;
218 const auto vViscousFlux = turbViscosity * averageDepth * gradV;
219
220 // compute the mobility of the flux with the fluxlimiter
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);
223
224 const auto limitingDepth = (waterDepthLeft + waterDepthRight) * 0.5;
225 const auto mobility = ShallowWater::fluxLimiterLET(limitingDepth,
226 limitingDepth,
227 upperWaterDepthFluxLimiting,
228 lowerWaterDepthFluxLimiting);
229
230 localFlux[0] = 0.0;
231 localFlux[1] = -uViscousFlux * mobility * scvf.area();
232 localFlux[2] = -vViscousFlux * mobility * scvf.area();
233
234 return localFlux;
235 }
236};
237
238} // end namespace Dumux
239
240#endif
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
Definition: adapt.hh:29
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