24#ifndef DUMUX_STAGGERED_LOWREKEPSILON_LOCAL_RESIDUAL_HH
25#define DUMUX_STAGGERED_LOWREKEPSILON_LOCAL_RESIDUAL_HH
27#include <dune/common/hybridutilities.hh>
40template<
class TypeTag,
class BaseLocalRes
idual,
class DiscretizationMethod>
41class LowReKEpsilonResidualImpl;
43template<
class TypeTag,
class BaseLocalRes
idual>
45:
public BaseLocalResidual
47 using ParentType = BaseLocalResidual;
51 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
52 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
53 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
55 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
56 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
59 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
60 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
65 using Element =
typename GridView::template Codim<0>::Entity;
67 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
72 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
73 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
76 using ParentType::ParentType;
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars)
const
83 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
85 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy() * volVars.density();
86 storage[dissipationEqIdx] = volVars.dissipationTilde() * volVars.density();
92 const Element &element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const ElementFaceVariables& elemFaceVars,
96 const SubControlVolume &scv)
const
98 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
99 elemVolVars, elemFaceVars, scv);
101 const auto& volVars = elemVolVars[scv];
104 source[turbulentKineticEnergyEqIdx] += 2.0 * volVars.dynamicEddyViscosity()
105 * volVars.stressTensorScalarProduct();
106 source[dissipationEqIdx] += volVars.cOneEpsilon() * volVars.fOne()
107 * volVars.dissipationTilde() / volVars.turbulentKineticEnergy()
108 * 2.0 * volVars.dynamicEddyViscosity()
109 * volVars.stressTensorScalarProduct();
112 source[turbulentKineticEnergyEqIdx] -= volVars.dissipationTilde() * volVars.density();
113 source[dissipationEqIdx] -= volVars.cTwoEpsilon() * volVars.fTwo() * volVars.density()
114 * volVars.dissipationTilde() * volVars.dissipationTilde()
115 / volVars.turbulentKineticEnergy();
118 source[turbulentKineticEnergyEqIdx] -= volVars.dValue() * volVars.density();
119 source[dissipationEqIdx] += volVars.eValue() * volVars.density();
The available discretization methods in Dumux.
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Element-wise calculation of the residual for low-Reynolds k-epsilon models using the staggered discre...
Definition: freeflow/rans/twoeq/lowrekepsilon/localresidual.hh:36
CellCenterPrimaryVariables computeSourceForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolume &scv) const
Definition: freeflow/rans/twoeq/lowrekepsilon/staggered/localresidual.hh:91
CellCenterPrimaryVariables computeStorageForCellCenter(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate fluxes entering or leaving the cell center control volume.
Definition: freeflow/rans/twoeq/lowrekepsilon/staggered/localresidual.hh:79
Declares all properties used in Dumux.