24#ifndef DUMUX_LOWREKEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_LOWREKEPSILON_STAGGERED_FLUXVARIABLES_HH
43template<
class TypeTag,
class BaseFluxVariables,
class DiscretizationMethod>
44class LowReKEpsilonFluxVariablesImpl;
46template<
class TypeTag,
class BaseFluxVariables>
48:
public BaseFluxVariables
50 using ParentType = BaseFluxVariables;
54 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
55 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
56 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
58 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
59 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
61 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
62 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
63 using FaceVariables =
typename GridFaceVariables::FaceVariables;
68 using FVElementGeometry =
typename GridGeometry::LocalView;
69 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
71 using GridView =
typename GridGeometry::GridView;
73 using Element =
typename GridView::template Codim<0>::Entity;
78 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
79 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
87 const Element &element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const ElementFaceVariables& elemFaceVars,
91 const SubControlVolumeFace &scvf,
92 const FluxVariablesCache& fluxVarsCache)
94 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
95 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
98 auto upwindTermK = [](
const auto& volVars)
100 return volVars.turbulentKineticEnergy() * volVars.density();
102 auto upwindTermEpsilon = [](
const auto& volVars)
104 return volVars.dissipationTilde() * volVars.density();
107 flux[turbulentKineticEnergyEqIdx]
108 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermK);
109 flux[dissipationEqIdx]
110 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermEpsilon);
113 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
114 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
115 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
116 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
119 Scalar insideCoeff_k = insideVolVars.viscosity() + insideVolVars.kinematicEddyViscosity()
120 * insideVolVars.density() / insideVolVars.sigmaK();
121 Scalar outsideCoeff_k = outsideVolVars.viscosity() + outsideVolVars.kinematicEddyViscosity()
122 * outsideVolVars.density() / outsideVolVars.sigmaK();
123 Scalar insideCoeff_e = insideVolVars.viscosity() + insideVolVars.kinematicEddyViscosity()
124 * insideVolVars.density() / insideVolVars.sigmaEpsilon();
125 Scalar outsideCoeff_e = outsideVolVars.viscosity() + outsideVolVars.kinematicEddyViscosity()
126 * outsideVolVars.density() / outsideVolVars.sigmaEpsilon();
129 insideCoeff_k *= insideVolVars.extrusionFactor();
130 outsideCoeff_k *= outsideVolVars.extrusionFactor();
131 insideCoeff_e *= insideVolVars.extrusionFactor();
132 outsideCoeff_e *= outsideVolVars.extrusionFactor();
134 Scalar coeff_k = 0.0;
135 Scalar coeff_e = 0.0;
139 coeff_k = insideCoeff_k;
140 coeff_e = insideCoeff_e;
141 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
147 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
148 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
150 (outsideScv.dofPosition() - scvf.ipGlobal()).two_norm(),
151 (insideScv.dofPosition() - scvf.ipGlobal()).two_norm());
152 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
155 const auto bcTypes = problem.boundaryTypes(element, scvf);
156 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
157 || bcTypes.isSymmetry())))
159 flux[turbulentKineticEnergyEqIdx]
161 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
162 * Extrusion::area(fvGeometry, scvf);
164 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
165 || bcTypes.isSymmetry())))
167 flux[dissipationEqIdx]
169 * (insideVolVars.dissipationTilde() - outsideVolVars.dissipationTilde())
170 * Extrusion::area(fvGeometry, scvf);
179 const Element& element,
180 const SubControlVolumeFace& scvf,
181 const FVElementGeometry& fvGeometry,
182 const ElementVolumeVariables& elemVolVars,
183 const ElementFaceVariables& elemFaceVars,
184 const GridFluxVariablesCache& gridFluxVarsCache)
186 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
188 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
189 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
190 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
191 * Extrusion::area(fvGeometry, 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:294
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
The flux variables class for the low-Reynolds k-epsilon model using the staggered grid discretization...
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:178
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:86
Declares all properties used in Dumux.