24#ifndef DUMUX_LOWREKEPSILON_PROBLEM_HH
25#define DUMUX_LOWREKEPSILON_PROBLEM_HH
45template<
class TypeTag>
62 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
71 ParentType::updateStaticWallProperties();
74 storedDissipationTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
75 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
76 storedTurbulentKineticEnergy_.resize(this->gridGeometry().elementMapper().size(), 0.0);
84 template<
class SolutionVector>
87 ParentType::updateDynamicWallProperties(curSol);
89 auto fvGeometry =
localView(this->gridGeometry());
90 for (
const auto& element : elements(this->gridGeometry().gridView()))
92 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
93 fvGeometry.bindElement(element);
95 for (
auto&& scv : scvs(fvGeometry))
97 const int dofIdx = scv.dofIndex();
98 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
99 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
100 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
102 storedDissipationTilde_[elementIdx] = elemSol[0][Indices::dissipationEqIdx];
103 storedTurbulentKineticEnergy_[elementIdx] = elemSol[0][Indices::turbulentKineticEnergyEqIdx];
105 VolumeVariables volVars;
106 volVars.update(elemSol, asImp_(), element, scv);
107 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
114 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
true);
115 return useStoredEddyViscosity;
119 {
return storedDissipationTilde_[elementIdx]; }
122 {
return storedDynamicEddyViscosity_[elementIdx]; }
125 {
return storedTurbulentKineticEnergy_[elementIdx]; }
128 std::vector<Scalar> storedDissipationTilde_;
129 std::vector<Scalar> storedDynamicEddyViscosity_;
130 std::vector<Scalar> storedTurbulentKineticEnergy_;
133 Implementation &asImp_()
134 {
return *
static_cast<Implementation *
>(
this); }
137 const Implementation &asImp_()
const
138 {
return *
static_cast<const Implementation *
>(
this); }
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.
Base class for all staggered fv problems.
TurbulenceModel
The available free flow turbulence models in Dumux.
Definition: turbulencemodel.hh:38
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
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
Scalar storedDissipationTilde(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:118
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:85
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:69
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:62
Scalar storedTurbulentKineticEnergy(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:124
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:121
bool useStoredEddyViscosity() const
Definition: freeflow/rans/twoeq/lowrekepsilon/problem.hh:112
The local element solution class for staggered methods.
Declares all properties used in Dumux.
Adaption of the fully implicit scheme to the tracer transport model.