3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/sst/problem.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_PROBLEM_HH
25#define DUMUX_SST_PROBLEM_HH
26
34
35#include "model.hh"
36
37namespace Dumux {
38
45template<class TypeTag>
46class RANSProblemImpl<TypeTag, TurbulenceModel::sst> : public RANSProblemBase<TypeTag>
47{
49 using Implementation = GetPropType<TypeTag, Properties::Problem>;
51
53 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
54
57 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
59
60 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
61 using DimVector = typename Element::Geometry::GlobalCoordinate;
62
63public:
64 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
65 : ParentType(gridGeometry, paramGroup)
66 {
67 sstModelVersion_ = sstModelFromString(getParamFromGroup<std::string>(paramGroup, "RANS.SSTModelVersion", "SST"));
68 }
69
74 {
75 ParentType::updateStaticWallProperties();
76 // update size and initial values of the global vectors
77 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
78 storedDissipation_.resize(this->gridGeometry().elementMapper().size(), 0.0);
79 storedDissipationGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
80 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
81 storedTurbulentKineticEnergyGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
82 }
83
89 template<class SolutionVector>
90 void updateDynamicWallProperties(const SolutionVector& curSol)
91 {
92 ParentType::updateDynamicWallProperties(curSol);
93
94 for (const auto& element : elements(this->gridGeometry().gridView()))
95 {
96 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
97
98 auto fvGeometry = localView(this->gridGeometry());
99 fvGeometry.bindElement(element);
100 for (auto&& scv : scvs(fvGeometry))
101 {
102 const int dofIdx = scv.dofIndex();
103 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
104 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
105 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
106 // NOTE: first update the turbulence quantities
107 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
108 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
109 // NOTE: then update the volVars
110 VolumeVariables volVars;
111 volVars.update(elemSol, asImp_(), element, scv);
112 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*this);
113 }
114 }
115
116 // calculate cell-centered gradients
117 for (const auto& element : elements(this->gridGeometry().gridView()))
118 {
119 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
120
121 for (unsigned int axisIdx = 0; axisIdx < DimVector::dimension; ++axisIdx)
122 {
123 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
124 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
125
126 // Cell centered TKE Gradient
127 storedTurbulentKineticEnergyGradient_[elementIdx][axisIdx]
128 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
129 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
130 // Cell centered Omega Gradient
131 storedDissipationGradient_[elementIdx][axisIdx]
132 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
133 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
134 }
135 }
136 }
137
139 {
140 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
141 return useStoredEddyViscosity;
142 }
143
144 Scalar storedDynamicEddyViscosity(const int elementIdx) const
145 { return storedDynamicEddyViscosity_[elementIdx]; }
146
147 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
148 { return storedTurbulentKineticEnergy_[elementIdx]; }
149
150 Scalar storedDissipation(const int elementIdx) const
151 { return storedDissipation_[elementIdx]; }
152
153 DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
154 { return storedTurbulentKineticEnergyGradient_[elementIdx]; }
155
156 DimVector storedDissipationGradient(const int elementIdx) const
157 { return storedDissipationGradient_[elementIdx]; }
158
160 { return sstModelVersion_; }
161
162private:
163 std::vector<Scalar> storedDynamicEddyViscosity_;
164 std::vector<Scalar> storedTurbulentKineticEnergy_;
165 std::vector<Scalar> storedDissipation_;
166 std::vector<DimVector> storedDissipationGradient_;
167 std::vector<DimVector> storedTurbulentKineticEnergyGradient_;
168
169 SSTModel sstModelVersion_;
170
172 Implementation &asImp_()
173 { return *static_cast<Implementation *>(this); }
174
176 const Implementation &asImp_() const
177 { return *static_cast<const Implementation *>(this); }
178};
179
180} // end namespace Dumux
181
182#endif
Base class for all staggered fv problems.
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
The available free flow turbulence models in Dumux.
TurbulenceModel
The available free flow turbulence models in Dumux.
Definition: turbulencemodel.hh:38
SSTModel
The available variations of the SST Turbulence Model.
Definition: turbulencemodel.hh:75
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
SSTModel sstModelFromString(const std::string &sstModel)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: turbulencemodel.hh:94
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:55
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:156
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/sst/problem.hh:138
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/sst/problem.hh:73
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:144
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:147
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:153
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Definition: freeflow/rans/twoeq/sst/problem.hh:64
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:150
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/sst/problem.hh:90
SSTModel sstModelVersion() const
Definition: freeflow/rans/twoeq/sst/problem.hh:159
Declares all properties used in Dumux.
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.