24#ifndef DUMUX_KOMEGA_PROBLEM_HH
25#define DUMUX_KOMEGA_PROBLEM_HH
45template<
class TypeTag>
61 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
62 using DimVector =
typename Element::Geometry::GlobalCoordinate;
65 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
66 : ParentType(gridGeometry, paramGroup)
74 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));
91 ParentType::updateDynamicWallProperties(curSol);
93 for (
const auto& element : elements(this->gridGeometry().gridView()))
95 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
97 auto fvGeometry =
localView(this->gridGeometry());
98 fvGeometry.bindElement(element);
99 for (
auto&& scv : scvs(fvGeometry))
101 const int dofIdx = scv.dofIndex();
102 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
103 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
104 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
106 storedDissipation_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
107 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
109 VolumeVariables volVars;
110 volVars.update(elemSol, asImp_(), element, scv);
111 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity(*
this);
116 for (
const auto& element : elements(this->gridGeometry().gridView()))
118 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
120 for (
unsigned int dimIdx = 0; dimIdx < DimVector::dimension; ++dimIdx)
122 const unsigned neighborIdx0 = ParentType::neighborIndex(elementIdx, dimIdx, 0);
123 const unsigned neighborIdx1 = ParentType::neighborIndex(elementIdx, dimIdx, 1);
126 storedTurbulentKineticEnergyGradient_[elementIdx][dimIdx]
127 = (storedTurbulentKineticEnergy(neighborIdx1) - storedTurbulentKineticEnergy(neighborIdx0))
128 / (ParentType::cellCenter(neighborIdx1)[dimIdx] - ParentType::cellCenter(neighborIdx0)[dimIdx]);
130 storedDissipationGradient_[elementIdx][dimIdx]
131 = (storedDissipation(neighborIdx1) - storedDissipation(neighborIdx0))
132 / (ParentType::cellCenter(neighborIdx1)[dimIdx] - ParentType::cellCenter(neighborIdx0)[dimIdx]);
143 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
144 return useStoredEddyViscosity;
148 {
return storedDynamicEddyViscosity_[elementIdx]; }
151 {
return storedTurbulentKineticEnergy_[elementIdx]; }
154 {
return storedDissipation_[elementIdx]; }
157 {
return storedTurbulentKineticEnergyGradient_[elementIdx]; }
160 {
return storedDissipationGradient_[elementIdx]; }
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_;
170 Implementation &asImp_()
171 {
return *
static_cast<Implementation *
>(
this); }
174 const Implementation &asImp_()
const
175 {
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:37
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:150
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/komega/problem.hh:89
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/komega/problem.hh:141
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:147
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/komega/problem.hh:72
DimVector storedTurbulentKineticEnergyGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:156
DimVector storedDissipationGradient(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:159
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
Definition: freeflow/rans/twoeq/komega/problem.hh:65
Scalar storedDissipation(const int elementIdx) const
Definition: freeflow/rans/twoeq/komega/problem.hh:153
const Scalar betaOmega() const
Returns the constant.
Definition: freeflow/rans/twoeq/komega/problem.hh:138
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.