3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/komega/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_KOMEGA_PROBLEM_HH
25#define DUMUX_KOMEGA_PROBLEM_HH
26
34
35#include "model.hh"
36
37namespace Dumux {
38
45template<class TypeTag>
46class RANSProblemImpl<TypeTag, TurbulenceModel::komega> : public RANSProblemBase<TypeTag>
47{
48 using ParentType = RANSProblemBase<TypeTag>;
49 using Implementation = GetPropType<TypeTag, Properties::Problem>;
51
53 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
54
58 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
60
61 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
62 using DimVector = typename Element::Geometry::GlobalCoordinate;
63
64public:
65 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
66 : ParentType(gridGeometry, paramGroup)
67 { }
68
73 {
74 ParentType::updateStaticWallProperties();
75
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 void updateDynamicWallProperties(const SolutionVector& curSol)
90 {
91 ParentType::updateDynamicWallProperties(curSol);
92
93 for (const auto& element : elements(this->gridGeometry().gridView()))
94 {
95 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
96
97 auto fvGeometry = localView(this->gridGeometry());
98 fvGeometry.bindElement(element);
99 for (auto&& scv : scvs(fvGeometry))
100 {
101 const int dofIdx = scv.dofIndex();
102 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
103 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
104 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
105 // NOTE: first update the turbulence quantities
106 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
107 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
108 // NOTE: then update the volVars
109 VolumeVariables volVars;
110 volVars.update(elemSol, asImp_(), element, scv);
111 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*this);
112 }
113 }
114
115 // calculate cell-centered gradients
116 for (const auto& element : elements(this->gridGeometry().gridView()))
117 {
118 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
119
120 for (unsigned int dimIdx = 0; dimIdx < DimVector::dimension; ++dimIdx)
121 {
122 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, dimIdx, 0);
123 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, dimIdx, 1);
124
125 // Cell centered TKE Gradient
126 storedTurbulentKineticEnergyGradient_[elementIdx][dimIdx]
127 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
128 / (ParentType::cellCenter(neighborIdx1)[dimIdx] - ParentType::cellCenter(neighborIdx0)[dimIdx]);
129 // Cell centered Omega Gradient
130 storedDissipationGradient_[elementIdx][dimIdx]
131 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
132 / (ParentType::cellCenter(neighborIdx1)[dimIdx] - ParentType::cellCenter(neighborIdx0)[dimIdx]);
133 }
134 }
135 }
136
138 const Scalar betaOmega() const
139 { return 0.0708; }
140
142 {
143 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
144 return useStoredEddyViscosity;
145 }
146
147 Scalar storedDynamicEddyViscosity(const int elementIdx) const
148 { return storedDynamicEddyViscosity_[elementIdx]; }
149
150 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
151 { return storedTurbulentKineticEnergy_[elementIdx]; }
152
153 Scalar storedDissipation(const int elementIdx) const
154 { return storedDissipation_[elementIdx]; }
155
156 DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
157 { return storedTurbulentKineticEnergyGradient_[elementIdx]; }
158
159 DimVector storedDissipationGradient(const int elementIdx) const
160 { return storedDissipationGradient_[elementIdx]; }
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
170 Implementation &asImp_()
171 { return *static_cast<Implementation *>(this); }
172
174 const Implementation &asImp_() const
175 { return *static_cast<const Implementation *>(this); }
176};
177
178} // end namespace Dumux
179
180#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:37
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:150
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/komega/problem.hh:89
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/komega/problem.hh:141
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:147
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/komega/problem.hh:72
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:156
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:159
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Definition: freeflow/rans/twoeq/komega/problem.hh:65
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:153
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/twoeq/komega/problem.hh:138
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.