24#ifndef DUMUX_LOWREKEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_LOWREKEPSILON_STAGGERED_FLUXVARIABLES_HH
37template<
class TypeTag,
class BaseFluxVariables, DiscretizationMethod discMethod>
38class LowReKEpsilonFluxVariablesImpl;
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.dissipationTilde();
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.kinematicViscosity()
116 + insideVolVars.kinematicEddyViscosity() / insideVolVars.sigmaK();
117 Scalar outsideCoeff_k = outsideVolVars.kinematicViscosity()
118 + outsideVolVars.kinematicEddyViscosity() / outsideVolVars.sigmaK();
119 Scalar insideCoeff_e = insideVolVars.kinematicViscosity()
120 + insideVolVars.kinematicEddyViscosity() / insideVolVars.sigmaEpsilon();
121 Scalar outsideCoeff_e = outsideVolVars.kinematicViscosity()
122 + outsideVolVars.kinematicEddyViscosity() / outsideVolVars.sigmaEpsilon();
125 insideCoeff_k *= insideVolVars.extrusionFactor();
126 outsideCoeff_k *= outsideVolVars.extrusionFactor();
127 insideCoeff_e *= insideVolVars.extrusionFactor();
128 outsideCoeff_e *= outsideVolVars.extrusionFactor();
132 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
133 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
135 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
136 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
137 Scalar distance = 0.0;
142 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
143 coeff_k = insideCoeff_k;
144 coeff_e = insideCoeff_e;
148 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
151 const auto bcTypes = problem.boundaryTypes(element, scvf);
152 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
153 || bcTypes.isSymmetry())))
155 flux[turbulentKineticEnergyEqIdx]
156 += coeff_k / distance
157 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
160 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
161 || bcTypes.isSymmetry())))
163 flux[dissipationEqIdx]
164 += coeff_e / distance
165 * (insideVolVars.dissipationTilde() - outsideVolVars.dissipationTilde())
175 const Element& element,
176 const SubControlVolumeFace& scvf,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const ElementFaceVariables& elemFaceVars,
180 const GridFluxVariablesCache& gridFluxVarsCache)
182 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
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();
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/lowrekepsilon/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/lowrekepsilon/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/lowrekepsilon/staggered/fluxvariables.hh:82
Declares all properties used in Dumux.