24#ifndef DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
37template<
class TypeTag,
class BaseFluxVariables, DiscretizationMethod discMethod>
38class KEpsilonFluxVariablesImpl;
44template<
class TypeTag,
class BaseFluxVariables>
46:
public BaseFluxVariables
48 using ParentType = BaseFluxVariables;
52 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
53 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
54 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
56 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
57 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
59 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
60 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
61 using FaceVariables =
typename GridFaceVariables::FaceVariables;
67 using Element =
typename GridView::template Codim<0>::Entity;
70 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
74 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
75 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
83 const Element &element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const ElementFaceVariables& elemFaceVars,
87 const SubControlVolumeFace &scvf,
88 const FluxVariablesCache& fluxVarsCache)
90 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
91 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
94 auto upwindTermK = [](
const auto& volVars)
96 return volVars.turbulentKineticEnergy();
98 auto upwindTermEpsilon = [](
const auto& volVars)
100 return volVars.dissipation();
103 flux[turbulentKineticEnergyEqIdx]
104 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermK);
105 flux[dissipationEqIdx ]
106 = ParentType::advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTermEpsilon);
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()];
115 Scalar insideCoeff_k = insideVolVars.kinematicEddyViscosity() / insideVolVars.sigmaK();
116 Scalar outsideCoeff_k = outsideVolVars.kinematicEddyViscosity() / outsideVolVars.sigmaK();
117 Scalar insideCoeff_e = insideVolVars.kinematicEddyViscosity() / insideVolVars.sigmaEpsilon();
118 Scalar outsideCoeff_e = outsideVolVars.kinematicEddyViscosity() / outsideVolVars.sigmaEpsilon();
119 static const auto kEpsilonEnableKinematicViscosity_
120 = getParamFromGroup<bool>(problem.paramGroup(),
"KEpsilon.EnableKinematicViscosity",
true);
121 if (kEpsilonEnableKinematicViscosity_)
123 insideCoeff_k += insideVolVars.kinematicViscosity();
124 outsideCoeff_k += outsideVolVars.kinematicViscosity();
125 insideCoeff_e += insideVolVars.kinematicViscosity();
126 outsideCoeff_e += outsideVolVars.kinematicViscosity();
130 insideCoeff_k *= insideVolVars.extrusionFactor();
131 outsideCoeff_k *= outsideVolVars.extrusionFactor();
132 insideCoeff_e *= insideVolVars.extrusionFactor();
133 outsideCoeff_e *= outsideVolVars.extrusionFactor();
139 Scalar distance = 0.0;
144 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
145 coeff_k = insideCoeff_k;
146 coeff_e = insideCoeff_e;
150 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
153 const auto bcTypes = problem.boundaryTypes(element, scvf);
154 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
155 || bcTypes.isSymmetry()
156 || problem.isOnWall(scvf))))
158 if (!(insideVolVars.isMatchingPoint() && outsideVolVars.isMatchingPoint())
159 || !(insideVolVars.isMatchingPoint() && outsideVolVars.inNearWallRegion())
160 || !(insideVolVars.inNearWallRegion() && outsideVolVars.isMatchingPoint()))
162 flux[turbulentKineticEnergyEqIdx]
163 += coeff_k / distance
164 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
168 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
169 || bcTypes.isSymmetry())))
171 flux[dissipationEqIdx]
172 += coeff_e / distance
173 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
183 const Element& element,
184 const SubControlVolumeFace& scvf,
185 const FVElementGeometry& fvGeometry,
186 const ElementVolumeVariables& elemVolVars,
187 const ElementFaceVariables& elemFaceVars,
188 const GridFluxVariablesCache& gridFluxVarsCache)
190 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
192 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
193 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
194 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
195 * scvf.area() * scvf.directionSign() * insideVolVars.extrusionFactor();
Base class for the flux variables living on a sub control volume face.
The available discretization methods in Dumux.
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
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/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:82
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:182
Declares all properties used in Dumux.