3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
25#define DUMUX_SST_STAGGERED_KOMEGA_LOCAL_RESIDUAL_HH
26
27#include <dune/common/hybridutilities.hh>
28#include <dune/common/exceptions.hh>
33
34namespace Dumux {
35
41// forward declaration
42template<class TypeTag, class BaseLocalResidual, class DiscretizationMethod>
43class SSTResidualImpl;
44
45template<class TypeTag, class BaseLocalResidual>
46class SSTResidualImpl<TypeTag, BaseLocalResidual, DiscretizationMethods::Staggered>
47: public BaseLocalResidual
48{
49 using ParentType = BaseLocalResidual;
50
52
53 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
54 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
55 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
56
57 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
58 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
60
61 using GridFaceVariables = typename GridVariables::GridFaceVariables;
62 using ElementFaceVariables = typename GridFaceVariables::LocalView;
63
67 using Element = typename GridView::template Codim<0>::Entity;
68 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
69 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
70 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
71 using CellCenterResidual = CellCenterPrimaryVariables;
74
75 static constexpr int turbulentKineticEnergyEqIdx = Indices::turbulentKineticEnergyEqIdx - ModelTraits::dim();
76 static constexpr int dissipationEqIdx = Indices::dissipationEqIdx - ModelTraits::dim();
77
78public:
79 using ParentType::ParentType;
80
82 CellCenterPrimaryVariables computeStorageForCellCenter(const Problem& problem,
83 const SubControlVolume& scv,
84 const VolumeVariables& volVars) const
85 {
86 CellCenterPrimaryVariables storage = ParentType::computeStorageForCellCenter(problem, scv, volVars);
87
88 storage[turbulentKineticEnergyEqIdx] = volVars.turbulentKineticEnergy()*volVars.density();
89 storage[dissipationEqIdx] = volVars.dissipation()*volVars.density();
90
91 return storage;
92 }
93
94 CellCenterPrimaryVariables computeSourceForCellCenter(const Problem& problem,
95 const Element &element,
96 const FVElementGeometry& fvGeometry,
97 const ElementVolumeVariables& elemVolVars,
98 const ElementFaceVariables& elemFaceVars,
99 const SubControlVolume &scv) const
100 {
101 CellCenterPrimaryVariables source = ParentType::computeSourceForCellCenter(problem, element, fvGeometry,
102 elemVolVars, elemFaceVars, scv);
103
104 using std::min;
105 const auto& volVars = elemVolVars[scv];
106
107 if(problem.sstModelVersion() == SSTModel::BSL)
108 {
109 // production
110 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
111 source[turbulentKineticEnergyEqIdx] += productionTerm;
112 source[dissipationEqIdx] += volVars.gammaBSL() / volVars.kinematicEddyViscosity() * productionTerm;
113
114 // destruction
115 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarBSL() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
116 source[dissipationEqIdx] -= volVars.betaBSL() * volVars.density() * volVars.dissipation() * volVars.dissipation();
117
118 // cross-diffusion term
119 Scalar gradientProduct = 0.0;
120 for (unsigned int i = 0; i < ModelTraits::dim(); ++i)
121 gradientProduct += volVars.storedTurbulentKineticEnergyGradient()[i]
122 * volVars.storedDissipationGradient()[i];
123
124 source[dissipationEqIdx] += 2.0 * volVars.density() * (1.0-volVars.F1()) * volVars.sigmaOmega2() / volVars.dissipation() * gradientProduct;
125 }
126 else if(problem.sstModelVersion() == SSTModel::SST)
127 {
128 // production
129 Scalar productionTerm = 2.0 * volVars.dynamicEddyViscosity() * volVars.stressTensorScalarProduct();
130 source[turbulentKineticEnergyEqIdx] += productionTerm;
131 source[dissipationEqIdx] += volVars.gammaSST() / volVars.kinematicEddyViscosity() * productionTerm;
132
133 // destruction
134 source[turbulentKineticEnergyEqIdx] -= volVars.betaStarSST() * volVars.density() * volVars.turbulentKineticEnergy() * volVars.dissipation();
135 source[dissipationEqIdx] -= volVars.betaSST() * volVars.density() * volVars.dissipation() * volVars.dissipation();
136
137 // cross-diffusion term
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;
143 }
144 else
145 DUNE_THROW(Dune::NotImplemented, "\nThis SST Model is not implemented.\n");
146
147 return source;
148 }
149};
150} // end namespace Dumux
151
152#endif
The available discretization methods in Dumux.
The available free flow turbulence models in Dumux.
Definition: adapt.hh:29
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 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.