3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/rans/twoeq/lowrekepsilon/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_LOWREKEPSILON_PROBLEM_HH
25#define DUMUX_LOWREKEPSILON_PROBLEM_HH
26
34
35#include "model.hh"
36
37namespace Dumux {
38
45template<class 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
61public:
63 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
64 : ParentType(gridGeometry, paramGroup)
65 { }
66
71 {
72 ParentType::updateStaticWallProperties();
73
74 // update size and initial values of the global vectors
75 storedDissipationTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
76 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
77 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
78 }
79
85 void updateDynamicWallProperties(const SolutionVector& curSol)
86 {
87 ParentType::updateDynamicWallProperties(curSol);
88
89 for (const auto& element : elements(this->gridGeometry().gridView()))
90 {
91 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
92
93 auto fvGeometry = localView(this->gridGeometry());
94 fvGeometry.bindElement(element);
95 for (auto&& scv : scvs(fvGeometry))
96 {
97 const int dofIdx = scv.dofIndex();
98 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
99 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
100 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
101 // NOTE: first update the turbulence quantities
102 storedDissipationTilde_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
103 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
104 // NOTE: then update the volVars
105 VolumeVariables volVars;
106 volVars.update(elemSol, asImp_(), element, scv);
107 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
108 }
109 }
110 }
111
113 {
114 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", true);
115 return useStoredEddyViscosity;
116 }
117
118 Scalar storedDissipationTilde(const int elementIdx) const
119 { return storedDissipationTilde_[elementIdx]; }
120
121 Scalar storedDynamicEddyViscosity(const int elementIdx) const
122 { return storedDynamicEddyViscosity_[elementIdx]; }
123
124 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
125 { return storedTurbulentKineticEnergy_[elementIdx]; }
126
127private:
128 std::vector<Scalar> storedDissipationTilde_;
129 std::vector<Scalar> storedDynamicEddyViscosity_;
130 std::vector<Scalar> storedTurbulentKineticEnergy_;
131
133 Implementation &asImp_()
134 { return *static_cast<Implementation *>(this); }
135
137 const Implementation &asImp_() const
138 { return *static_cast<const Implementation *>(this); }
139};
140
141}
142
143#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 storedDissipationTilde(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:118
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:70
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:63
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:85
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:124
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:121
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:112
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.