24#ifndef DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
44template<
class TypeTag,
class BaseFluxVariables,
class DiscretizationMethod>
45class KEpsilonFluxVariablesImpl;
47template<
class TypeTag,
class BaseFluxVariables>
49:
public BaseFluxVariables
51 using ParentType = BaseFluxVariables;
55 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
56 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
57 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
59 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
60 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
62 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
63 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
64 using FaceVariables =
typename GridFaceVariables::FaceVariables;
69 using FVElementGeometry =
typename GridGeometry::LocalView;
70 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
71 using GlobalPosition =
typename SubControlVolumeFace::GlobalPosition;
73 using GridView =
typename GridGeometry::GridView;
75 using Element =
typename GridView::template Codim<0>::Entity;
80 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
81 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
89 const Element &element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const ElementFaceVariables& elemFaceVars,
93 const SubControlVolumeFace &scvf,
94 const FluxVariablesCache& fluxVarsCache)
96 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
97 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
100 auto upwindTermK = [](
const auto& volVars)
102 return volVars.turbulentKineticEnergy() * volVars.density();
104 auto upwindTermEpsilon = [](
const auto& volVars)
106 return volVars.dissipation() * volVars.density();
109 flux[turbulentKineticEnergyEqIdx]
110 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
111 flux[dissipationEqIdx ]
112 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermEpsilon);
115 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
116 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
117 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
118 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
121 Scalar insideCoeff_k = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaK()) + insideVolVars.viscosity();
122 Scalar outsideCoeff_k = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaK()) + outsideVolVars.viscosity();
123 Scalar insideCoeff_e = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaEpsilon()) + insideVolVars.viscosity();
124 Scalar outsideCoeff_e = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaEpsilon()) + outsideVolVars.viscosity();
127 insideCoeff_k *= insideVolVars.extrusionFactor();
128 outsideCoeff_k *= outsideVolVars.extrusionFactor();
129 insideCoeff_e *= insideVolVars.extrusionFactor();
130 outsideCoeff_e *= outsideVolVars.extrusionFactor();
132 Scalar coeff_k = 0.0;
133 Scalar coeff_e = 0.0;
137 coeff_k = insideCoeff_k;
138 coeff_e = insideCoeff_e;
139 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
147 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
150 const auto bcTypes = problem.boundaryTypes(element, scvf);
153 if constexpr (Deprecated::hasIsOnWall<Problem, GlobalPosition>())
156 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
157 || bcTypes.isSymmetry()
158 || problem.isOnWall(scvf))))
160 if (!(insideVolVars.isMatchingPoint() && outsideVolVars.isMatchingPoint())
161 || !(insideVolVars.isMatchingPoint() && outsideVolVars.inNearWallRegion())
162 || !(insideVolVars.inNearWallRegion() && outsideVolVars.isMatchingPoint()))
164 flux[turbulentKineticEnergyEqIdx]
166 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
167 * Extrusion::area(scvf);
174 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
175 || bcTypes.isSymmetry()
176 || bcTypes.hasWall())))
178 if (!(insideVolVars.isMatchingPoint() && outsideVolVars.isMatchingPoint())
179 || !(insideVolVars.isMatchingPoint() && outsideVolVars.inNearWallRegion())
180 || !(insideVolVars.inNearWallRegion() && outsideVolVars.isMatchingPoint()))
182 flux[turbulentKineticEnergyEqIdx]
184 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
185 * Extrusion::area(scvf);
190 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
191 || bcTypes.isSymmetry())))
193 flux[dissipationEqIdx]
195 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
196 * Extrusion::area(scvf);
205 const Element& element,
206 const SubControlVolumeFace& scvf,
207 const FVElementGeometry& fvGeometry,
208 const ElementVolumeVariables& elemVolVars,
209 const ElementFaceVariables& elemFaceVars,
210 const GridFluxVariablesCache& gridFluxVarsCache)
212 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
214 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
215 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
216 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
217 * Extrusion::area(scvf) * scvf.directionSign() * insideVolVars.extrusionFactor();
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.
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
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 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: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/kepsilon/staggered/fluxvariables.hh:204
Declares all properties used in Dumux.