3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/sst/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_SST_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_SST_STAGGERED_FLUXVARIABLES_HH
26
27#include <numeric>
28#include <dune/common/exceptions.hh>
36
37namespace Dumux {
38
44// forward declaration
45template<class TypeTag, class BaseFluxVariables, class DiscretizationMethod>
46class SSTFluxVariablesImpl;
47
48template<class TypeTag, class BaseFluxVariables>
49class SSTFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethods::Staggered>
50: public BaseFluxVariables
51{
52 using ParentType = BaseFluxVariables;
53
55
56 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
57 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
58 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
59
60 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
61 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
62
63 using GridFaceVariables = typename GridVariables::GridFaceVariables;
64 using ElementFaceVariables = typename GridFaceVariables::LocalView;
65 using FaceVariables = typename GridFaceVariables::FaceVariables;
66
70 using FVElementGeometry = typename GridGeometry::LocalView;
71 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
72 using Extrusion = Extrusion_t<GridGeometry>;
73 using GridView = typename GridGeometry::GridView;
75 using Element = typename GridView::template Codim<0>::Entity;
77 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
79
80 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
81 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
82
83public:
84
88 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
89 const Element &element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const ElementFaceVariables& elemFaceVars,
93 const SubControlVolumeFace &scvf,
94 const FluxVariablesCache& fluxVarsCache)
95 {
96 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
97 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
98
99 // calculate advective flux
100 auto upwindTermK = [](const auto& volVars)
101 { return volVars.turbulentKineticEnergy() * volVars.density(); };
102 auto upwindTermOmega = [](const auto& volVars)
103 { return volVars.dissipation() * volVars.density(); };
104
105 flux[turbulentKineticEnergyEqIdx]
106 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
107 flux[dissipationEqIdx]
108 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermOmega);
109
110 // calculate diffusive flux
111 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
112 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
113 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
114 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
115 Scalar insideCoeff_k = 0.0, insideCoeff_w = 0.0, outsideCoeff_k = 0.0, outsideCoeff_w = 0.0;
116
117 if(problem.sstModelVersion() == SSTModel::BSL)
118 {
119 insideCoeff_k = insideVolVars.viscosity()
120 + ( insideVolVars.sigmaKBSL() * insideVolVars.dynamicEddyViscosity() );
121 outsideCoeff_k = outsideVolVars.viscosity()
122 + ( outsideVolVars.sigmaKBSL() * outsideVolVars.dynamicEddyViscosity() );
123 insideCoeff_w = insideVolVars.viscosity()
124 + ( insideVolVars.sigmaOmegaBSL() * insideVolVars.dynamicEddyViscosity() );
125 outsideCoeff_w = outsideVolVars.viscosity()
126 + ( outsideVolVars.sigmaOmegaBSL() * outsideVolVars.dynamicEddyViscosity() );
127 }
128 else if(problem.sstModelVersion() == SSTModel::SST)
129 {
130 insideCoeff_k = insideVolVars.viscosity()
131 + ( insideVolVars.sigmaKSST() * insideVolVars.dynamicEddyViscosity() );
132 outsideCoeff_k = outsideVolVars.viscosity()
133 + ( outsideVolVars.sigmaKSST() * outsideVolVars.dynamicEddyViscosity() );
134 insideCoeff_w = insideVolVars.viscosity()
135 + ( insideVolVars.sigmaOmegaSST() * insideVolVars.dynamicEddyViscosity() );
136 outsideCoeff_w = outsideVolVars.viscosity()
137 + ( outsideVolVars.sigmaOmegaSST() * outsideVolVars.dynamicEddyViscosity() );
138 }
139 else
140 DUNE_THROW(Dune::NotImplemented, "\nThis SST Model is not implemented.\n");
141
142
143 // scale by extrusion factor
144 insideCoeff_k *= insideVolVars.extrusionFactor();
145 outsideCoeff_k *= outsideVolVars.extrusionFactor();
146 insideCoeff_w *= insideVolVars.extrusionFactor();
147 outsideCoeff_w *= outsideVolVars.extrusionFactor();
148
149 Scalar distance = 0.0;
150 Scalar coeff_k = 0.0;
151 Scalar coeff_w = 0.0;
152 if (scvf.boundary())
153 {
154 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
155 coeff_k = insideCoeff_k;
156 coeff_w = insideCoeff_w;
157 }
158 else
159 {
160 // average and distance
161 coeff_k = arithmeticMean(insideCoeff_k, outsideCoeff_k,
162 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
163 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
164 coeff_w = arithmeticMean(insideCoeff_w, outsideCoeff_w,
165 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
166 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
167 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
168 }
169
170 const auto bcTypes = problem.boundaryTypes(element, scvf);
171 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
172 || bcTypes.isSymmetry())))
173 {
174 flux[turbulentKineticEnergyEqIdx]
175 += coeff_k / distance
176 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
177 * Extrusion::area(scvf);
178 }
179 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
180 || bcTypes.isSymmetry())))
181 {
182 flux[dissipationEqIdx]
183 += coeff_w / distance
184 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
185 * Extrusion::area(scvf);
186 }
187 return flux;
188 }
189
193 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
194 const Element& element,
195 const SubControlVolumeFace& scvf,
196 const FVElementGeometry& fvGeometry,
197 const ElementVolumeVariables& elemVolVars,
198 const ElementFaceVariables& elemFaceVars,
199 const GridFluxVariablesCache& gridFluxVarsCache)
200 {
201 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
202
203 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
204 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
205 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
206 * Extrusion::area(scvf) * scvf.directionSign() * insideVolVars.extrusionFactor();
207 }
208};
209
210} // end namespace
211
212#endif
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
The available free flow turbulence models in Dumux.
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:292
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:50
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
Definition: propertysystem.hh:150
The flux variables class for the SST model using the staggered grid discretization.
Definition: freeflow/rans/twoeq/sst/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/sst/staggered/fluxvariables.hh:88
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/sst/staggered/fluxvariables.hh:193
Declares all properties used in Dumux.