12#ifndef DUMUX_SST_PROBLEM_HH
13#define DUMUX_SST_PROBLEM_HH
33template<
class TypeTag>
48 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
49 using DimVector =
typename Element::Geometry::GlobalCoordinate;
52 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
55 sstModelVersion_ =
sstModelFromString(getParamFromGroup<std::string>(paramGroup,
"RANS.SSTModelVersion",
"SST"));
63 ParentType::updateStaticWallProperties();
65 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
66 storedDissipation_.resize(this->gridGeometry().elementMapper().size(), 0.0);
67 storedDissipationGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
68 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
69 storedTurbulentKineticEnergyGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
77 template<
class SolutionVector>
80 ParentType::updateDynamicWallProperties(curSol);
82 for (
const auto& element : elements(this->gridGeometry().gridView()))
84 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
86 auto fvGeometry =
localView(this->gridGeometry());
87 fvGeometry.bindElement(element);
88 for (
auto&& scv : scvs(fvGeometry))
90 const int dofIdx = scv.dofIndex();
91 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
92 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
93 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
95 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
96 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
98 VolumeVariables volVars;
99 volVars.update(elemSol, asImp_(), element, scv);
100 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*
this);
105 for (
const auto& element : elements(this->gridGeometry().gridView()))
107 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
109 for (
unsigned int axisIdx = 0; axisIdx < DimVector::dimension; ++axisIdx)
111 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
112 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
115 storedTurbulentKineticEnergyGradient_[elementIdx][axisIdx]
116 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
117 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
119 storedDissipationGradient_[elementIdx][axisIdx]
120 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
121 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
128 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
129 return useStoredEddyViscosity;
133 {
return storedDynamicEddyViscosity_[elementIdx]; }
136 {
return storedTurbulentKineticEnergy_[elementIdx]; }
139 {
return storedDissipation_[elementIdx]; }
142 {
return storedTurbulentKineticEnergyGradient_[elementIdx]; }
145 {
return storedDissipationGradient_[elementIdx]; }
148 {
return sstModelVersion_; }
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_;
160 Implementation &asImp_()
161 {
return *
static_cast<Implementation *
>(
this); }
164 const Implementation &asImp_()
const
165 {
return *
static_cast<const Implementation *
>(
this); }
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
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:144
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/sst/problem.hh:126
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/sst/problem.hh:61
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:132
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:135
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:141
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
Definition: freeflow/rans/twoeq/sst/problem.hh:52
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:138
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/sst/problem.hh:78
SSTModel sstModelVersion() const
Definition: freeflow/rans/twoeq/sst/problem.hh:147
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
SSTModel
The available variations of the SST Turbulence Model.
Definition: turbulencemodel.hh:63
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
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:82
The local element solution class for staggered methods.
Base class for all staggered fv problems.
The available free flow turbulence models in Dumux.