version 3.8
freeflow/rans/twoeq/sst/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_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
13#define DUMUX_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
14
15#include <dune/common/hybridutilities.hh>
16#include <dune/common/exceptions.hh>
21
22namespace Dumux {
23
29// forward declaration
30template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
31class SSTResidualImpl;
32
33template<class TypeTag, class BaseLocalResidual>
34class SSTResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
35: public BaseLocalResidual
36{
37 using ParentType = BaseLocalResidual;
38
40
41 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
42 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
43 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
44
45 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
46 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
48
49 using GridFaceVariables = typename GridVariables::GridFaceVariables;
50 using ElementFaceVariables = typename GridFaceVariables::LocalView;
51
55 using Element = typename GridView::template Codim<0>::Entity;
56 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
57 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
58 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
59 using CellCenterResidual = CellCenterPrimaryVariables;
62
63 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
64 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
65
66public:
67 using ParentType::ParentType;
68
70 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
71 const SubControlVolume& scv,
72 const VolumeVariables& volVars) const
73 {
74 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
75
76 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy()*volVars.density();
77 storage[dissipationEqIdx] = volVars.dissipation()*volVars.density();
78
79 return storage;
80 }
81
82 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
83 const Element &element,
84 const FVElementGeometry& fvGeometry,
85 const ElementVolumeVariables& elemVolVars,
86 const ElementFaceVariables& elemFaceVars,
87 const SubControlVolume &scv) const
88 {
89 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
90 elemVolVars, elemFaceVars, scv);
91
92 using std::min;
93 const auto& volVars = elemVolVars[scv];
94
95 if(problem.sstModelVersion() == SSTModel::BSL)
96 {
97 // production
98 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
99 source[turbulentKineticEnergyEqIdx] += productionTerm;
100 source[dissipationEqIdx] += volVars.gammaBSL() / volVars.kinematicEddyViscosity() * productionTerm;
101
102 // destruction
103 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarBSL() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
104 source[dissipationEqIdx] -= volVars.betaBSL() * volVars.density() * volVars.dissipation() * volVars.dissipation();
105
106 // cross-diffusion term
107 Scalar gradientProduct = 0.0;
108 for (unsigned int i = 0; i < ModelTraits::dim(); ++i)
109 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
110 * volVars.storedDissipationGradient()[i];
111
112 source[dissipationEqIdx] += 2.0 * volVars.density() * (1.0-volVars.F1()) * volVars.sigmaOmega2() / volVars.dissipation() * gradientProduct;
113 }
114 else if(problem.sstModelVersion() == SSTModel::SST)
115 {
116 // production
117 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
118 source[turbulentKineticEnergyEqIdx] += productionTerm;
119 source[dissipationEqIdx] += volVars.gammaSST() / volVars.kinematicEddyViscosity() * productionTerm;
120
121 // destruction
122 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarSST() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
123 source[dissipationEqIdx] -= volVars.betaSST() * volVars.density() * volVars.dissipation() * volVars.dissipation();
124
125 // cross-diffusion term
126 Scalar gradientProduct = 0.0;
127 for (unsigned int i = 0; i < ModelTraits::dim(); ++i)
128 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
129 * volVars.storedDissipationGradient()[i];
130 source[dissipationEqIdx] += 2.0 * volVars.density() * (1.0-volVars.F1()) * volVars.sigmaOmega2() / volVars.dissipation() * gradientProduct;
131 }
132 else
133 DUNE_THROW(Dune::NotImplemented, "\nThis SST Model is not implemented.\n");
134
135 return source;
136 }
137};
138} // end namespace Dumux
139
140#endif
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:82
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:70
Element-wise calculation of the residual for SST models using the staggered discretization.
Definition: freeflow/rans/twoeq/sst/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
The available free flow turbulence models in Dumux.