24#ifndef DUMUX_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
25#define DUMUX_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
27#include <dune/common/hybridutilities.hh>
28#include <dune/common/exceptions.hh>
42template<
class TypeTag,
class BaseLocalRes
idual,
class DiscretizationMethod>
45template<
class TypeTag,
class BaseLocalRes
idual>
46class SSTResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
47:
public BaseLocalResidual
49 using ParentType = BaseLocalResidual;
53 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
54 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
55 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
57 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
58 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
61 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
62 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
67 using Element =
typename GridView::template Codim<0>::Entity;
69 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
71 using CellCenterResidual = CellCenterPrimaryVariables;
75 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
76 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
79 using ParentType::ParentType;
83 const SubControlVolume& scv,
84 const VolumeVariables& volVars)
const
86 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
88 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy()*volVars.density();
89 storage[dissipationEqIdx] = volVars.dissipation()*volVars.density();
95 const Element &element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const ElementFaceVariables& elemFaceVars,
99 const SubControlVolume &scv)
const
101 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
102 elemVolVars, elemFaceVars, scv);
105 const auto& volVars = elemVolVars[scv];
110 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
111 source[turbulentKineticEnergyEqIdx] += productionTerm;
112 source[dissipationEqIdx] += volVars.gammaBSL() / volVars.kinematicEddyViscosity() * productionTerm;
115 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarBSL() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
116 source[dissipationEqIdx] -= volVars.betaBSL() * volVars.density() * volVars.dissipation() * volVars.dissipation();
119 Scalar gradientProduct = 0.0;
120 for (
unsigned int i = 0; i < ModelTraits::dim(); ++i)
121 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
122 * volVars.storedDissipationGradient()[i];
124 source[dissipationEqIdx] += 2.0 * volVars.density() * (1.0-volVars.F1()) * volVars.sigmaOmega2() / volVars.dissipation() * gradientProduct;
129 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
130 source[turbulentKineticEnergyEqIdx] += productionTerm;
131 source[dissipationEqIdx] += volVars.gammaSST() / volVars.kinematicEddyViscosity() * productionTerm;
134 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarSST() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
135 source[dissipationEqIdx] -= volVars.betaSST() * volVars.density() * volVars.dissipation() * volVars.dissipation();
138 Scalar gradientProduct = 0.0;
139 for (
unsigned int i = 0; i < ModelTraits::dim(); ++i)
140 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
141 * volVars.storedDissipationGradient()[i];
142 source[dissipationEqIdx] += 2.0 * volVars.density() * (1.0-volVars.F1()) * volVars.sigmaOmega2() / volVars.dissipation() * gradientProduct;
145 DUNE_THROW(Dune::NotImplemented,
"\nThis SST Model is not implemented.\n");
The available discretization methods in Dumux.
The available free flow turbulence models in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
Element-wise calculation of the residual for SST models using the staggered discretization.
Definition: freeflow/rans/twoeq/sst/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/sst/staggered/localresidual.hh:94
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/sst/staggered/localresidual.hh:82
Declares all properties used in Dumux.