version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_KOMEGA_STAGGERED_FLUXVARIABLES_HH
13#define DUMUX_KOMEGA_STAGGERED_FLUXVARIABLES_HH
14
15#include <numeric>
22
23namespace Dumux {
24
30// forward declaration
31template<class TypeTag, class BaseFluxVariables, class DiscretizationMethod>
32class KOmegaFluxVariablesImpl;
33
34template<class TypeTag, class BaseFluxVariables>
35class KOmegaFluxVariablesImpl<TypeTag, BaseFluxVariables, DiscretizationMethods::Staggered>
36: public BaseFluxVariables
37{
38 using ParentType = BaseFluxVariables;
39
41
42 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
43 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
44 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
45
46 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
47 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
48
49 using GridFaceVariables = typename GridVariables::GridFaceVariables;
50 using ElementFaceVariables = typename GridFaceVariables::LocalView;
51 using FaceVariables = typename GridFaceVariables::FaceVariables;
52
56 using FVElementGeometry = typename GridGeometry::LocalView;
57 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
58 using Extrusion = Extrusion_t<GridGeometry>;
59 using GridView = typename GridGeometry::GridView;
61 using Element = typename GridView::template Codim<0>::Entity;
63 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
65
66 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
67 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
68
69public:
70
74 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
75 const Element &element,
76 const FVElementGeometry& fvGeometry,
77 const ElementVolumeVariables& elemVolVars,
78 const ElementFaceVariables& elemFaceVars,
79 const SubControlVolumeFace &scvf,
80 const FluxVariablesCache& fluxVarsCache)
81 {
82 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
83 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
84
85 // calculate advective flux
86 auto upwindTermK = [](const auto& volVars)
87 {
88 return volVars.turbulentKineticEnergy()*volVars.density();
89 };
90 auto upwindTermOmega = [](const auto& volVars)
91 {
92 return volVars.dissipation()*volVars.density();
93 };
94
95 flux[turbulentKineticEnergyEqIdx]
96 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermK);
97 flux[dissipationEqIdx]
98 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermOmega);
99
100 // calculate diffusive flux
101 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
102 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
103 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
104 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
105
106 // effective diffusion coefficients
107 Scalar insideCoeff_k = insideVolVars.viscosity()
108 + ( insideVolVars.sigmaK() * insideVolVars.density()* insideVolVars.turbulentKineticEnergy() / insideVolVars.dissipation() );
109 Scalar outsideCoeff_k = outsideVolVars.viscosity()
110 + ( outsideVolVars.sigmaK() * outsideVolVars.density()* outsideVolVars.turbulentKineticEnergy() / outsideVolVars.dissipation() );
111 Scalar insideCoeff_w = insideVolVars.viscosity()
112 + ( insideVolVars.sigmaOmega() * insideVolVars.density()* insideVolVars.turbulentKineticEnergy() / insideVolVars.dissipation() );
113 Scalar outsideCoeff_w = outsideVolVars.viscosity()
114 + ( outsideVolVars.sigmaOmega() * outsideVolVars.density()* outsideVolVars.turbulentKineticEnergy() / outsideVolVars.dissipation() );
115
116 // scale by extrusion factor
117 insideCoeff_k *= insideVolVars.extrusionFactor();
118 outsideCoeff_k *= outsideVolVars.extrusionFactor();
119 insideCoeff_w *= insideVolVars.extrusionFactor();
120 outsideCoeff_w *= outsideVolVars.extrusionFactor();
121
122 Scalar distance = 0.0;
123 Scalar coeff_k = 0.0;
124 Scalar coeff_w = 0.0;
125 if (scvf.boundary())
126 {
127 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
128 coeff_k = insideCoeff_k;
129 coeff_w = insideCoeff_w;
130 }
131 else
132 {
133 // average and distance
134 coeff_k = arithmeticMean(insideCoeff_k, outsideCoeff_k,
135 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
136 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
137 coeff_w = arithmeticMean(insideCoeff_w, outsideCoeff_w,
138 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
139 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
140 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
141 }
142
143 const auto bcTypes = problem.boundaryTypes(element, scvf);
144 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
145 || bcTypes.isSymmetry())))
146 {
147 flux[turbulentKineticEnergyEqIdx]
148 += coeff_k / distance
149 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
150 * Extrusion::area(fvGeometry, scvf);
151 }
152 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
153 || bcTypes.isSymmetry())))
154 {
155 flux[dissipationEqIdx]
156 += coeff_w / distance
157 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
158 * Extrusion::area(fvGeometry, scvf);
159 }
160 return flux;
161 }
162
166 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
167 const Element& element,
168 const SubControlVolumeFace& scvf,
169 const FVElementGeometry& fvGeometry,
170 const ElementVolumeVariables& elemVolVars,
171 const ElementFaceVariables& elemFaceVars,
172 const GridFluxVariablesCache& gridFluxVarsCache)
173 {
174 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
175
176 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
177 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
178 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
179 * Extrusion::area(fvGeometry, scvf) * scvf.directionSign() * insideVolVars.extrusionFactor();
180 }
181};
182
183} // end namespace
184
185#endif
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:166
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:74
The flux variables class for the k-omega model using the staggered grid discretization.
Definition: freeflow/rans/twoeq/komega/fluxvariables.hh:22
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Base class for the flux variables living on a sub control volume face.
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:38
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:282
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166