24#ifndef DUMUX_STAGGERED_ONEEQ_LOCAL_RESIDUAL_HH
25#define DUMUX_STAGGERED_ONEEQ_LOCAL_RESIDUAL_HH
27#include <dune/common/hybridutilities.hh>
41template<
class TypeTag,
class BaseLocalRes
idual,
class DiscretizationMethod>
42class OneEqResidualImpl;
44template<
class TypeTag,
class BaseLocalRes
idual>
46:
public BaseLocalResidual
48 using ParentType = BaseLocalResidual;
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 ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
60 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
61 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
66 using Element =
typename GridView::template Codim<0>::Entity;
68 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
73 static constexpr int viscosityTildeEqIdx = Indices::viscosityTildeEqIdx - ModelTraits::dim();
76 using ParentType::ParentType;
80 const SubControlVolume& scv,
81 const VolumeVariables& volVars)
const
83 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
84 storage[viscosityTildeEqIdx] = volVars.viscosityTilde() * volVars.density();
89 const Element &element,
90 const FVElementGeometry& fvGeometry,
91 const ElementVolumeVariables& elemVolVars,
92 const ElementFaceVariables& elemFaceVars,
93 const SubControlVolume &scv)
const
95 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
96 elemVolVars, elemFaceVars, scv);
98 const auto& volVars = elemVolVars[scv];
100 source[viscosityTildeEqIdx] += volVars.cb1() * (1.0 - volVars.ft2())
101 * volVars.stressTensorScalarProductTilde()
102 * volVars.viscosityTilde() * volVars.density();
104 source[viscosityTildeEqIdx] -= (volVars.cw1() * volVars.fW()
105 - volVars.cb1() * volVars.ft2() / problem.karmanConstant() / problem.karmanConstant())
106 * volVars.viscosityTilde() * volVars.viscosityTilde()
107 / volVars.wallDistance() / volVars.wallDistance() * volVars.density();;
109 for (
unsigned int axisIdx = 0; axisIdx < ModelTraits::dim(); ++axisIdx)
111 source[viscosityTildeEqIdx] += volVars.cb2() / volVars.sigma()
112 * volVars.storedViscosityTildeGradient()[axisIdx]
113 * volVars.storedViscosityTildeGradient()[axisIdx]
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 one-equation turbulence models using the staggered discr...
Definition: freeflow/rans/oneeq/localresidual.hh:36
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/oneeq/staggered/localresidual.hh:79
CellCenterPrimaryVariables computeSourceForCellCenter(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolume &scv) const
Definition: freeflow/rans/oneeq/staggered/localresidual.hh:88
Declares all properties used in Dumux.