3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/komega/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_KOMEGA_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_KOMEGA_STAGGERED_FLUXVARIABLES_HH
26
27#include <numeric>
33
34namespace Dumux {
35
36// forward declaration
37template<class TypeTag, class BaseFluxVariables, DiscretizationMethod discMethod>
38class KOmegaFluxVariablesImpl;
39
44template<class TypeTag, class BaseFluxVariables>
45class KOmegaFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethod::staggered>
46: public BaseFluxVariables
47{
48 using ParentType = BaseFluxVariables;
49
51
52 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
53 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
54 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
55
56 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
57 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
58
59 using GridFaceVariables = typename GridVariables::GridFaceVariables;
60 using ElementFaceVariables = typename GridFaceVariables::LocalView;
61 using FaceVariables = typename GridFaceVariables::FaceVariables;
62
67 using Element = typename GridView::template Codim<0>::Entity;
68 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
70 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
71 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
73
74 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
75 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
76
77public:
78
82 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
83 const Element &element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const ElementFaceVariables& elemFaceVars,
87 const SubControlVolumeFace &scvf,
88 const FluxVariablesCache& fluxVarsCache)
89 {
90 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
91 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
92
93 // calculate advective flux
94 auto upwindTermK = [](const auto& volVars)
95 {
96 return volVars.turbulentKineticEnergy();
97 };
98 auto upwindTermOmega = [](const auto& volVars)
99 {
100 return volVars.dissipation();
101 };
102
103 flux[turbulentKineticEnergyEqIdx]
104 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
105 flux[dissipationEqIdx]
106 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermOmega);
107
108 // calculate diffusive flux
109 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
110 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
111 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
112 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
113
114 // effective diffusion coefficients
115 Scalar insideCoeff_k = insideVolVars.kinematicViscosity()
116 + ( insideVolVars.sigmaK() * insideVolVars.turbulentKineticEnergy() / insideVolVars.dissipation() );
117 Scalar outsideCoeff_k = outsideVolVars.kinematicViscosity()
118 + ( outsideVolVars.sigmaK() * outsideVolVars.turbulentKineticEnergy() / outsideVolVars.dissipation() );
119 Scalar insideCoeff_w = insideVolVars.kinematicViscosity()
120 + ( insideVolVars.sigmaOmega() * insideVolVars.turbulentKineticEnergy() / insideVolVars.dissipation() );
121 Scalar outsideCoeff_w = outsideVolVars.kinematicViscosity()
122 + ( outsideVolVars.sigmaOmega() * outsideVolVars.turbulentKineticEnergy() / outsideVolVars.dissipation() );
123
124 // scale by extrusion factor
125 insideCoeff_k *= insideVolVars.extrusionFactor();
126 outsideCoeff_k *= outsideVolVars.extrusionFactor();
127 insideCoeff_w *= insideVolVars.extrusionFactor();
128 outsideCoeff_w *= outsideVolVars.extrusionFactor();
129
130 // average and distance
131 Scalar coeff_k = arithmeticMean(insideCoeff_k, outsideCoeff_k,
132 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
133 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
134 Scalar coeff_w = arithmeticMean(insideCoeff_w, outsideCoeff_w,
135 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
136 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
137 Scalar distance = 0.0;
138
139 // adapt boundary handling
140 if (scvf.boundary())
141 {
142 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
143 coeff_k = insideCoeff_k;
144 coeff_w = insideCoeff_w;
145 }
146 else
147 {
148 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
149 }
150
151 const auto bcTypes = problem.boundaryTypes(element, scvf);
152 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
153 || bcTypes.isSymmetry())))
154 {
155 flux[turbulentKineticEnergyEqIdx]
156 += coeff_k / distance
157 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
158 * scvf.area();
159 }
160 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
161 || bcTypes.isSymmetry())))
162 {
163 flux[dissipationEqIdx]
164 += coeff_w / distance
165 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
166 * scvf.area();
167 }
168 return flux;
169 }
170
174 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
175 const Element& element,
176 const SubControlVolumeFace& scvf,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const ElementFaceVariables& elemFaceVars,
180 const GridFluxVariablesCache& gridFluxVarsCache)
181 {
182 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
183
184 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
185 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
186 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
187 * scvf.area() * scvf.directionSign() * insideVolVars.extrusionFactor();
188 }
189};
190
191} // end namespace
192
193#endif
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
Definition: freeflow/rans/twoeq/komega/fluxvariables.hh:34
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/komega/staggered/fluxvariables.hh:174
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/komega/staggered/fluxvariables.hh:82
Declares all properties used in Dumux.