24#ifndef DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_KEPSILON_STAGGERED_FLUXVARIABLES_HH
43template<
class TypeTag,
class BaseFluxVariables,
class DiscretizationMethod>
44class KEpsilonFluxVariablesImpl;
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;
70 using GlobalPosition =
typename SubControlVolumeFace::GlobalPosition;
72 using GridView =
typename GridGeometry::GridView;
74 using Element =
typename GridView::template Codim<0>::Entity;
79 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
80 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
88 const Element &element,
89 const FVElementGeometry& fvGeometry,
90 const ElementVolumeVariables& elemVolVars,
91 const ElementFaceVariables& elemFaceVars,
92 const SubControlVolumeFace &scvf,
93 const FluxVariablesCache& fluxVarsCache)
95 CellCenterPrimaryVariables flux = ParentType::computeMassFlux(problem, element, fvGeometry,
96 elemVolVars, elemFaceVars, scvf, fluxVarsCache);
99 auto upwindTermK = [](
const auto& volVars)
101 return volVars.turbulentKineticEnergy() * volVars.density();
103 auto upwindTermEpsilon = [](
const auto& volVars)
105 return volVars.dissipation() * volVars.density();
108 flux[turbulentKineticEnergyEqIdx]
109 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermK);
110 flux[dissipationEqIdx ]
111 = ParentType::advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTermEpsilon);
114 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
115 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
116 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
117 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
120 Scalar insideCoeff_k = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaK()) + insideVolVars.viscosity();
121 Scalar outsideCoeff_k = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaK()) + outsideVolVars.viscosity();
122 Scalar insideCoeff_e = (insideVolVars.dynamicEddyViscosity() / insideVolVars.sigmaEpsilon()) + insideVolVars.viscosity();
123 Scalar outsideCoeff_e = (outsideVolVars.dynamicEddyViscosity() / outsideVolVars.sigmaEpsilon()) + outsideVolVars.viscosity();
126 insideCoeff_k *= insideVolVars.extrusionFactor();
127 outsideCoeff_k *= outsideVolVars.extrusionFactor();
128 insideCoeff_e *= insideVolVars.extrusionFactor();
129 outsideCoeff_e *= outsideVolVars.extrusionFactor();
131 Scalar coeff_k = 0.0;
132 Scalar coeff_e = 0.0;
136 coeff_k = insideCoeff_k;
137 coeff_e = insideCoeff_e;
138 distance = (insideScv.dofPosition() - scvf.ipGlobal()).two_norm();
146 distance = (outsideScv.dofPosition() - insideScv.dofPosition()).two_norm();
149 const auto bcTypes = problem.boundaryTypes(element, scvf);
151 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::turbulentKineticEnergyEqIdx)
152 || bcTypes.isSymmetry()
153 || bcTypes.hasWall())))
155 if (!(insideVolVars.isMatchingPoint() && outsideVolVars.isMatchingPoint())
156 || !(insideVolVars.isMatchingPoint() && outsideVolVars.inNearWallRegion())
157 || !(insideVolVars.inNearWallRegion() && outsideVolVars.isMatchingPoint()))
159 flux[turbulentKineticEnergyEqIdx]
161 * (insideVolVars.turbulentKineticEnergy() - outsideVolVars.turbulentKineticEnergy())
162 * Extrusion::area(fvGeometry, scvf);
166 if (!(scvf.boundary() && (bcTypes.isOutflow(Indices::dissipationEqIdx)
167 || bcTypes.isSymmetry())))
169 flux[dissipationEqIdx]
171 * (insideVolVars.dissipation() - outsideVolVars.dissipation())
172 * Extrusion::area(fvGeometry, scvf);
181 const Element& element,
182 const SubControlVolumeFace& scvf,
183 const FVElementGeometry& fvGeometry,
184 const ElementVolumeVariables& elemVolVars,
185 const ElementFaceVariables& elemFaceVars,
186 const GridFluxVariablesCache& gridFluxVarsCache)
188 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
190 return ParentType::computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
191 + ParentType::computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache)
192 + 2.0 / ModelTraits::dim() * insideVolVars.density() * insideVolVars.turbulentKineticEnergy()
193 * 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 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:87
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:180
Declares all properties used in Dumux.