3.6-git
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
73template<class NumEqVector, typename std::enable_if_t<NumEqVector::size() == 3, int> = 0>
75{
76
77public:
78
86 template<class Problem, class FVElementGeometry, class ElementVolumeVariables>
87 static NumEqVector flux(const Problem& problem,
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)
92 {
93 using Scalar = typename NumEqVector::value_type;
94 using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
95
96 NumEqVector localFlux(0.0);
97
98 // Get the inside and outside volume variables
99 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
100 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
101
102 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
103 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
104
105 const auto gradVelocity = [&]()
106 {
107 // The left (inside) and right (outside) states
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);
112
113 // Compute the velocity gradients normal to the face
114 // Factor that takes the direction of the unit vector into account
115 const auto& cellCenterToCellCenter = outsideScv.center() - insideScv.center();
116 const auto distance = cellCenterToCellCenter.two_norm();
117 const auto& unitNormal = scvf.unitOuterNormal();
118 const auto direction = (unitNormal*cellCenterToCellCenter)/distance;
119 return std::array<Scalar, 2>{
120 (velocityXRight-velocityXLeft)*direction/distance,
121 (velocityYRight-velocityYLeft)*direction/distance
122 };
123 }();
124
125 // Use a harmonic average of the depth at the interface.
126 const auto waterDepthLeft = insideVolVars.waterDepth();
127 const auto waterDepthRight = outsideVolVars.waterDepth();
128 const auto averageDepth = 2.0*(waterDepthLeft*waterDepthRight)/(waterDepthLeft + waterDepthRight);
129
130 // compute the turbulent viscosity contribution
131 const Scalar turbViscosity = [&]()
132 {
133 // The (constant) background turbulent viscosity
134 static const auto turbBGViscosity = getParamFromGroup<Scalar>(problem.paramGroup(), "ShallowWater.TurbulentViscosity", 1.0e-6);
135
136 // Check whether the mixing-length turbulence model is used
137 static const auto useMixingLengthTurbulenceModel = getParamFromGroup<bool>(problem.paramGroup(), "ShallowWater.UseMixingLengthTurbulenceModel", false);
138
139 // constant eddy viscosity equal to the prescribed background eddy viscosity
140 if (!useMixingLengthTurbulenceModel)
141 return turbBGViscosity;
142
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!");
148 else
149 {
150 // turbulence model based on mixing length
151 // Compute the turbulent viscosity using a combined horizontal/vertical mixing length approach
152 // Turbulence coefficients: vertical (Elder like) and horizontal (Smagorinsky like)
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);
155
162 constexpr Scalar kappa = 0.41;
163 // Compute the average shear velocity on the face
164 const Scalar ustar = [&]()
165 {
166 // Get the bottom shear stress in the two adjacent cells
167 // Note that element is not needed in spatialParams().frictionLaw (should be removed). For now we simply pass the same element twice
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();
172
173 // Use a harmonic average of the viscosity and the depth at the interface.
174 using std::max;
175 const auto averageBottomShearStress = 2.0*(bottomShearStressInsideMag*bottomShearStressOutsideMag)
176 / max(1.0e-8,bottomShearStressInsideMag + bottomShearStressOutsideMag);
177
178 // Shear velocity possibly needed in mixing-length turbulence model in the computation of the diffusion flux
179 using std::sqrt;
180 static_assert(!FluidSystem::isCompressible(0),
181 "The shallow water model assumes incompressible fluids"
182 );
183 // Since the shallow water equations are incompressible, the density is constant.
184 // Therefore, insideVolVars.density() is equal to outsideVolVars.density().
185 return sqrt(averageBottomShearStress/insideVolVars.density());
186 }();
187
188 const auto turbViscosityV = turbConstV * (kappa/6.0) * ustar * averageDepth;
189
214 using std::sqrt;
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);
218
219 // Total turbulent viscosity
220 return turbBGViscosity + sqrt(turbViscosityV*turbViscosityV + turbViscosityH*turbViscosityH);
221 }
222 }();
223
224 // Compute the viscous momentum fluxes with connecting water depth at the interface
225 using std::min;
226 using std::max;
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;
233
234 localFlux[0] = 0.0;
235 localFlux[1] = -uViscousFlux * scvf.area();
236 localFlux[2] = -vViscousFlux * scvf.area();
237
238 return localFlux;
239 }
240};
241
242} // end namespace Dumux
243
244#endif
Function to limit the fluxes.
Classes related to flux variables caching.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
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