version 3.8
freeflow/rans/twoeq/komega/staggered/localresidual.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
13#define DUMUX_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
19
20namespace Dumux {
21
27// forward declaration
28template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
29class KOmegaResidualImpl;
30
31template<class TypeTag, class BaseLocalResidual>
32class KOmegaResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
33: public BaseLocalResidual
34{
35 using ParentType = BaseLocalResidual;
36
38
39 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
40 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
41 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
42
43 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
44 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
46
47 using GridFaceVariables = typename GridVariables::GridFaceVariables;
48 using ElementFaceVariables = typename GridFaceVariables::LocalView;
49
53 using Element = typename GridView::template Codim<0>::Entity;
54 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
55 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
56 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
57 using CellCenterResidual = CellCenterPrimaryVariables;
60
61 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
62 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
63
64public:
65 using ParentType::ParentType;
66
68 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
69 const SubControlVolume& scv,
70 const VolumeVariables& volVars) const
71 {
72 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
73
74 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy()*volVars.density();
75 storage[dissipationEqIdx] = volVars.dissipation()*volVars.density();
76
77 return storage;
78 }
79
80 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
81 const Element &element,
82 const FVElementGeometry& fvGeometry,
83 const ElementVolumeVariables& elemVolVars,
84 const ElementFaceVariables& elemFaceVars,
85 const SubControlVolume &scv) const
86 {
87 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
88 elemVolVars, elemFaceVars, scv);
89
90 using std::min;
91 const auto& volVars = elemVolVars[scv];
92
93 // production
94 static const auto enableKOmegaProductionLimiter
95 = getParamFromGroup<bool>(problem.paramGroup(), "KOmega.EnableProductionLimiter", false);
96 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
97 if (enableKOmegaProductionLimiter)
98 {
99 Scalar productionAlternative = 20.0 * volVars.density() * volVars.betaK() * volVars.turbulentKineticEnergy() * volVars.dissipation();
100 productionTerm = min(productionTerm, productionAlternative);
101 }
102 source[turbulentKineticEnergyEqIdx] += productionTerm;
103 source[dissipationEqIdx] += volVars.alpha() * volVars.dissipation() / volVars.turbulentKineticEnergy() * productionTerm;
104
105 // destruction
106 source[turbulentKineticEnergyEqIdx] -= volVars.betaK() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
107 source[dissipationEqIdx] -= volVars.betaOmega() * volVars.density() * volVars.dissipation() * volVars.dissipation();
108
109 // cross-diffusion term
110 Scalar gradientProduct = 0.0;
111 for (unsigned int i = 0; i < ModelTraits::dim(); ++i)
112 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
113 * volVars.storedDissipationGradient()[i];
114 if (gradientProduct > 0.0)
115 source[dissipationEqIdx] += 0.125 *volVars.density() / volVars.dissipation() * gradientProduct;
116
117 return source;
118 }
119};
120} // end namespace Dumux
121
122#endif
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/komega/staggered/localresidual.hh:68
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/komega/staggered/localresidual.hh:80
Element-wise calculation of the residual for k-omega models using the staggered discretization.
Definition: freeflow/rans/twoeq/komega/localresidual.hh:24
Defines all properties used in Dumux.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
Definition: adapt.hh:17