version 3.10-dev
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// 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_KOMEGA_PROBLEM_HH
13#define DUMUX_KOMEGA_PROBLEM_HH
14
22
23#include "model.hh"
24
25namespace Dumux {
26
33template<class TypeTag>
34class RANSProblemImpl<TypeTag, TurbulenceModel::komega> : public RANSProblemBase<TypeTag>
35{
37 using Implementation = GetPropType<TypeTag, Properties::Problem>;
39
41 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
42
45 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
47
48 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
49 using DimVector = typename Element::Geometry::GlobalCoordinate;
50
51public:
52 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
53 : ParentType(gridGeometry, paramGroup)
54 { }
55
60 {
61 ParentType::updateStaticWallProperties();
62
63 // update size and initial values of the global vectors
64 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
65 storedDissipation_.resize(this->gridGeometry().elementMapper().size(), 0.0);
66 storedDissipationGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
67 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
68 storedTurbulentKineticEnergyGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
69 }
70
76 template<class SolutionVector>
77 void updateDynamicWallProperties(const SolutionVector& curSol)
78 {
79 ParentType::updateDynamicWallProperties(curSol);
80
81 auto fvGeometry = localView(this->gridGeometry());
82 for (const auto& element : elements(this->gridGeometry().gridView()))
83 {
84 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
85 fvGeometry.bindElement(element);
86
87 for (auto&& scv : scvs(fvGeometry))
88 {
89 const int dofIdx = scv.dofIndex();
90 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
91 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
92 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
93 // NOTE: first update the turbulence quantities
94 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
95 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
96 // NOTE: then update the volVars
97 VolumeVariables volVars;
98 volVars.update(elemSol, asImp_(), element, scv);
99 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*this);
100 }
101 }
102
103 // calculate cell-centered gradients
104 for (const auto& element : elements(this->gridGeometry().gridView()))
105 {
106 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
107
108 for (unsigned int axisIdx = 0; axisIdx < DimVector::dimension; ++axisIdx)
109 {
110 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
111 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
112
113 // Cell centered TKE Gradient
114 storedTurbulentKineticEnergyGradient_[elementIdx][axisIdx]
115 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
116 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
117 // Cell centered Omega Gradient
118 storedDissipationGradient_[elementIdx][axisIdx]
119 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
120 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
121 }
122 }
123 }
124
126 const Scalar betaOmega() const
127 { return 0.0708; }
128
130 {
131 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
132 return useStoredEddyViscosity;
133 }
134
135 Scalar storedDynamicEddyViscosity(const int elementIdx) const
136 { return storedDynamicEddyViscosity_[elementIdx]; }
137
138 Scalar storedTurbulentKineticEnergy(const int elementIdx) const
139 { return storedTurbulentKineticEnergy_[elementIdx]; }
140
141 Scalar storedDissipation(const int elementIdx) const
142 { return storedDissipation_[elementIdx]; }
143
144 DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
145 { return storedTurbulentKineticEnergyGradient_[elementIdx]; }
146
147 DimVector storedDissipationGradient(const int elementIdx) const
148 { return storedDissipationGradient_[elementIdx]; }
149
150private:
151 std::vector<Scalar> storedDynamicEddyViscosity_;
152 std::vector<Scalar> storedTurbulentKineticEnergy_;
153 std::vector<Scalar> storedDissipation_;
154 std::vector<DimVector> storedDissipationGradient_;
155 std::vector<DimVector> storedTurbulentKineticEnergyGradient_;
156
158 Implementation &asImp_()
159 { return *static_cast<Implementation *>(this); }
160
162 const Implementation &asImp_() const
163 { return *static_cast<const Implementation *>(this); }
164};
165
166} // end namespace Dumux
167
168#endif
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:47
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:138
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/komega/problem.hh:129
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:135
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/komega/problem.hh:59
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:144
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:147
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/komega/problem.hh:77
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
Definition: freeflow/rans/twoeq/komega/problem.hh:52
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:141
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/twoeq/komega/problem.hh:126
forward declare
Definition: freeflow/rans/problem.hh:31
Defines all properties used in Dumux.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
TurbulenceModel
The available free flow turbulence models in Dumux.
Definition: turbulencemodel.hh:26
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
Definition: adapt.hh:17
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.
Base class for all staggered fv problems.
The available free flow turbulence models in Dumux.