version 3.8
freeflow/rans/oneeq/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_ONEEQ_PROBLEM_HH
13#define DUMUX_ONEEQ_PROBLEM_HH
14
22
23#include "model.hh"
24
25namespace Dumux {
26
33template<class TypeTag>
34class RANSProblemImpl<TypeTag, TurbulenceModel::oneeq> : public RANSProblemBase<TypeTag>
35{
36
38 using Implementation = GetPropType<TypeTag, Properties::Problem>;
40 using Grid = typename GridView::Grid;
41
43 using DimVector = Dune::FieldVector<Scalar, Grid::dimension>;
44
46 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
47 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
48
51 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
53
54public:
56 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
57 : ParentType(gridGeometry, paramGroup)
58 { }
59
64 {
65 ParentType::updateStaticWallProperties();
66
67 // update size and initial values of the global vectors
68 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
69 storedViscosityTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
70 storedViscosityTildeGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
71 }
72
78 template<class SolutionVector>
79 void updateDynamicWallProperties(const SolutionVector& curSol)
80 {
81 ParentType::updateDynamicWallProperties(curSol);
82
83 auto fvGeometry = localView(this->gridGeometry());
84 for (const auto& element : elements(this->gridGeometry().gridView()))
85 {
86 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
87 fvGeometry.bindElement(element);
88
89 for (auto&& scv : scvs(fvGeometry))
90 {
91 const int dofIdx = scv.dofIndex();
92 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
93 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
94 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
95 // NOTE: first update the turbulence quantities
96 storedViscosityTilde_[elementIdx] = elemSol[0][Indices::viscosityTildeIdx];
97 // NOTE: then update the volVars
98 VolumeVariables volVars;
99 volVars.update(elemSol, asImp_(), element, scv);
100 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
101 }
102 }
103
104 // calculate cell-center-averaged gradient
105 for (const auto& element : elements(this->gridGeometry().gridView()))
106 {
107 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
108 fvGeometry.bindElement(element);
109
110 for (unsigned int axisIdx = 0; axisIdx < Grid::dimension; ++axisIdx)
111 {
112 const unsigned int neighborIndex0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
113 const unsigned int neighborIndex1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
114
115 // calculate cell-centered turbulentEddyViscosity (viscosityTilde) gradient
116 storedViscosityTildeGradient_[elementIdx][axisIdx]
117 = (storedViscosityTilde(neighborIndex1) - storedViscosityTilde(neighborIndex0))
118 / (ParentType::cellCenter(neighborIndex1)[axisIdx] - ParentType::cellCenter(neighborIndex0)[axisIdx]);
119 }
120
121 // Adjust for dirichlet boundary conditions
122 for (auto&& scvf : scvfs(fvGeometry))
123 {
124 const unsigned int normDim = scvf.directionIndex();
125 if (scvf.boundary() && asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::viscosityTildeIdx))
126 {
127 // face Value
128 Scalar dirichletViscosityTilde = asImp_().dirichlet(element, scvf)[Indices::viscosityTildeIdx];
129
130 unsigned int neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 0);
131 if (scvf.center()[normDim] < ParentType::cellCenter(elementIdx)[normDim])
132 neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 1);
133
134 storedViscosityTildeGradient_[elementIdx][normDim]
135 = (storedViscosityTilde(neighborIndex) - dirichletViscosityTilde)
136 / (ParentType::cellCenter(neighborIndex)[normDim] - scvf.center()[normDim]);
137 }
138 }
139 }
140 }
141
143 {
144 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(), "RANS.UseStoredEddyViscosity", false);
145 return useStoredEddyViscosity;
146 }
147
148 Scalar storedDynamicEddyViscosity(const int elementIdx) const
149 { return storedDynamicEddyViscosity_[elementIdx]; }
150
151 Scalar storedViscosityTilde(const int elementIdx) const
152 { return storedViscosityTilde_[elementIdx]; }
153
154 DimVector storedViscosityTildeGradient(const int elementIdx) const
155 { return storedViscosityTildeGradient_[elementIdx]; }
156
157private:
158 std::vector<Scalar> storedDynamicEddyViscosity_;
159 std::vector<Scalar> storedViscosityTilde_;
160 std::vector<DimVector> storedViscosityTildeGradient_;
161
163 Implementation &asImp_()
164 { return *static_cast<Implementation *>(this); }
165
167 const Implementation &asImp_() const
168 { return *static_cast<const Implementation *>(this); }
169};
170
171} // end namespace Dumux
172
173#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
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/oneeq/problem.hh:63
bool useStoredEddyViscosity() const
Definition: freeflow/rans/oneeq/problem.hh:142
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/oneeq/problem.hh:79
Scalar storedViscosityTilde(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:151
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/oneeq/problem.hh:56
DimVector storedViscosityTildeGradient(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:154
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:148
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.