24#ifndef DUMUX_ONEEQ_PROBLEM_HH
25#define DUMUX_ONEEQ_PROBLEM_HH
45template<
class TypeTag>
52 using Grid =
typename GridView::Grid;
55 using DimVector = Dune::FieldVector<Scalar, Grid::dimension>;
59 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
69 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
78 ParentType::updateStaticWallProperties();
81 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
82 storedViscosityTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
83 storedViscosityTildeGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
93 ParentType::updateDynamicWallProperties(curSol);
95 auto fvGeometry =
localView(this->gridGeometry());
96 for (
const auto& element : elements(this->gridGeometry().gridView()))
98 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
99 fvGeometry.bindElement(element);
101 for (
auto&& scv : scvs(fvGeometry))
103 const int dofIdx = scv.dofIndex();
104 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
105 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
106 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
108 storedViscosityTilde_[elementIdx] = elemSol[0][Indices::viscosityTildeIdx];
110 VolumeVariables volVars;
111 volVars.update(elemSol, asImp_(), element, scv);
112 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
117 for (
const auto& element : elements(this->gridGeometry().gridView()))
119 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
120 fvGeometry.bindElement(element);
122 for (
unsigned int axisIdx = 0; axisIdx < Grid::dimension; ++axisIdx)
124 const unsigned int neighborIndex0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
125 const unsigned int neighborIndex1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
128 storedViscosityTildeGradient_[elementIdx][axisIdx]
129 = (storedViscosityTilde(neighborIndex1) - storedViscosityTilde(neighborIndex0))
130 / (ParentType::cellCenter(neighborIndex1)[axisIdx] - ParentType::cellCenter(neighborIndex0)[axisIdx]);
134 for (
auto&& scvf : scvfs(fvGeometry))
136 const unsigned int normDim = scvf.directionIndex();
137 if (scvf.boundary() && asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::viscosityTildeIdx))
140 Scalar dirichletViscosityTilde = asImp_().dirichlet(element, scvf)[Indices::viscosityTildeIdx];
142 unsigned int neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 0);
143 if (scvf.center()[normDim] < ParentType::cellCenter(elementIdx)[normDim])
144 neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 1);
146 storedViscosityTildeGradient_[elementIdx][normDim]
147 = (storedViscosityTilde(neighborIndex) - dirichletViscosityTilde)
148 / (ParentType::cellCenter(neighborIndex)[normDim] - scvf.center()[normDim]);
156 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
157 return useStoredEddyViscosity;
161 {
return storedDynamicEddyViscosity_[elementIdx]; }
164 {
return storedViscosityTilde_[elementIdx]; }
167 {
return storedViscosityTildeGradient_[elementIdx]; }
170 std::vector<Scalar> storedDynamicEddyViscosity_;
171 std::vector<Scalar> storedViscosityTilde_;
172 std::vector<DimVector> storedViscosityTildeGradient_;
175 Implementation &asImp_()
176 {
return *
static_cast<Implementation *
>(
this); }
179 const Implementation &asImp_()
const
180 {
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
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
Definition: propertysystem.hh:150
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/oneeq/problem.hh:76
bool useStoredEddyViscosity() const
Definition: freeflow/rans/oneeq/problem.hh:154
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/oneeq/problem.hh:91
Scalar storedViscosityTilde(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:163
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor sets the gravity, if desired by the user.
Definition: freeflow/rans/oneeq/problem.hh:69
DimVector storedViscosityTildeGradient(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:166
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:160
forward declare
Definition: freeflow/rans/problem.hh:44
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:60
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.