3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/kepsilon/staggered/fluxvariables.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_KEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
26
27#include <numeric>
34
35namespace Dumux {
36
42// forward declaration
43template<class TypeTag, class BaseFluxVariables, DiscretizationMethod discMethod>
44class KEpsilonFluxVariablesImpl;
45
46template<class TypeTag, class BaseFluxVariables>
47class KEpsilonFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethod::staggered>
48: public BaseFluxVariables
49{
50 using ParentType = BaseFluxVariables;
51
53
54 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
55 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
56 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
57
58 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
59 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
60
61 using GridFaceVariables = typename GridVariables::GridFaceVariables;
62 using ElementFaceVariables = typename GridFaceVariables::LocalView;
63 using FaceVariables = typename GridFaceVariables::FaceVariables;
64
68 using FVElementGeometry = typename GridGeometry::LocalView;
69 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
70 using Extrusion = Extrusion_t<GridGeometry>;
71 using GridView = typename GridGeometry::GridView;
73 using Element = typename GridView::template Codim<0>::Entity;
75 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
77
78 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
79 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
80
81public:
82
86 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
87 const Element &element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const ElementFaceVariables& elemFaceVars,
91 const SubControlVolumeFace &scvf,
92 const FluxVariablesCache& fluxVarsCache)
93 {
94 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
95 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
96
97 // calculate advective flux
98 auto upwindTermK = [](const auto& volVars)
99 {
100 return volVars.turbulentKineticEnergy() * volVars.density();
101 };
102 auto upwindTermEpsilon = [](const auto& volVars)
103 {
104 return volVars.dissipation() * volVars.density();
105 };
106
107 flux[turbulentKineticEnergyEqIdx]
108 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
109 flux[dissipationEqIdx ]
110 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermEpsilon);
111
112 // calculate diffusive flux
113 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
114 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
115 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
116 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
117
118 // effective diffusion coefficients
119 Scalar insideCoeff_k = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaK()) + insideVolVars.viscosity();
120 Scalar outsideCoeff_k = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaK()) + outsideVolVars.viscosity();
121 Scalar insideCoeff_e = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaEpsilon()) + insideVolVars.viscosity();
122 Scalar outsideCoeff_e = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaEpsilon()) + outsideVolVars.viscosity();
123
124 // scale by extrusion factor
125 insideCoeff_k *= insideVolVars.extrusionFactor();
126 outsideCoeff_k *= outsideVolVars.extrusionFactor();
127 insideCoeff_e *= insideVolVars.extrusionFactor();
128 outsideCoeff_e *= outsideVolVars.extrusionFactor();
129
130 Scalar coeff_k = 0.0;
131 Scalar coeff_e = 0.0;
132 Scalar distance = 0.0;
133 if (scvf.boundary())
134 {
135 coeff_k = insideCoeff_k;
136 coeff_e = insideCoeff_e;
137 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
138 }
139 else
140 {
141 // average and distance
142 // is more stable with simple/unweighted arithmetic mean
143 coeff_k = arithmeticMean(insideCoeff_k, outsideCoeff_k);
144 coeff_e = arithmeticMean(insideCoeff_e, outsideCoeff_e);
145 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
146 }
147
148 const auto bcTypes = problem.boundaryTypes(element, scvf);
149 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
150 || bcTypes.isSymmetry()
151 || problem.isOnWall(scvf))))
152 {
153 if (!(insideVolVars.isMatchingPoint() && outsideVolVars.isMatchingPoint())
154 || !(insideVolVars.isMatchingPoint() && outsideVolVars.inNearWallRegion())
155 || !(insideVolVars.inNearWallRegion() && outsideVolVars.isMatchingPoint()))
156 {
157 flux[turbulentKineticEnergyEqIdx]
158 += coeff_k / distance
159 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
160 * Extrusion::area(scvf);
161 }
162 }
163 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
164 || bcTypes.isSymmetry())))
165 {
166 flux[dissipationEqIdx]
167 += coeff_e / distance
168 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
169 * Extrusion::area(scvf);
170 }
171 return flux;
172 }
173
177 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
178 const Element& element,
179 const SubControlVolumeFace& scvf,
180 const FVElementGeometry& fvGeometry,
181 const ElementVolumeVariables& elemVolVars,
182 const ElementFaceVariables& elemFaceVars,
183 const GridFluxVariablesCache& gridFluxVarsCache)
184 {
185 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
186
187 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
188 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
189 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
190 * Extrusion::area(scvf) * scvf.directionSign() * insideVolVars.extrusionFactor();
191 }
192};
193
194} // end namespace
195
196#endif
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
Base class for the flux variables living on a sub control volume face.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
constexpr Scalar arithmeticMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) arithmetic mean of two scalar values.
Definition: math.hh:49
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
The flux variables class for the k-epsilon model using the staggered grid discretization.
Definition: freeflow/rans/twoeq/kepsilon/fluxvariables.hh:34
CellCenterPrimaryVariables computeMassFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Computes the flux for the cell center residual.
Definition: freeflow/rans/twoeq/kepsilon/staggered/fluxvariables.hh:86
FacePrimaryVariables computeMomentumFlux(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
Returns the momentum flux over all staggered faces.
Definition: freeflow/rans/twoeq/kepsilon/staggered/fluxvariables.hh:177
Declares all properties used in Dumux.