24#ifndef DUMUX_SST_PROBLEM_HH
25#define DUMUX_SST_PROBLEM_HH
45template<
class TypeTag>
60 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
61 using DimVector =
typename Element::Geometry::GlobalCoordinate;
64 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
67 sstModelVersion_ =
sstModelFromString(getParamFromGroup<std::string>(paramGroup,
"RANS.SSTModelVersion",
"SST"));
75 ParentType::updateStaticWallProperties();
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));
89 template<
class SolutionVector>
92 ParentType::updateDynamicWallProperties(curSol);
94 for (
const auto& element : elements(this->gridGeometry().gridView()))
96 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
98 auto fvGeometry =
localView(this->gridGeometry());
99 fvGeometry.bindElement(element);
100 for (
auto&& scv : scvs(fvGeometry))
102 const int dofIdx = scv.dofIndex();
103 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
104 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
105 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
107 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
108 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
110 VolumeVariables volVars;
111 volVars.update(elemSol, asImp_(), element, scv);
112 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*
this);
117 for (
const auto& element : elements(this->gridGeometry().gridView()))
119 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
121 for (
unsigned int axisIdx = 0; axisIdx < DimVector::dimension; ++axisIdx)
123 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
124 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
127 storedTurbulentKineticEnergyGradient_[elementIdx][axisIdx]
128 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
129 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
131 storedDissipationGradient_[elementIdx][axisIdx]
132 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
133 / (ParentType::cellCenter(neighborIdx1)[axisIdx] - ParentType::cellCenter(neighborIdx0)[axisIdx]);
140 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
141 return useStoredEddyViscosity;
145 {
return storedDynamicEddyViscosity_[elementIdx]; }
148 {
return storedTurbulentKineticEnergy_[elementIdx]; }
151 {
return storedDissipation_[elementIdx]; }
154 {
return storedTurbulentKineticEnergyGradient_[elementIdx]; }
157 {
return storedDissipationGradient_[elementIdx]; }
160 {
return sstModelVersion_; }
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_;
172 Implementation &asImp_()
173 {
return *
static_cast<Implementation *
>(
this); }
176 const Implementation &asImp_()
const
177 {
return *
static_cast<const Implementation *
>(
this); }
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:38
SSTModel
The available variations of the SST Turbulence Model.
Definition: turbulencemodel.hh:75
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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:94
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:55
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:156
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/sst/problem.hh:138
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/sst/problem.hh:73
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:144
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:147
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:153
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
Definition: freeflow/rans/twoeq/sst/problem.hh:64
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/sst/problem.hh:150
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/sst/problem.hh:90
SSTModel sstModelVersion() const
Definition: freeflow/rans/twoeq/sst/problem.hh:159
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.