12#ifndef DUMUX_ONEEQ_PROBLEM_HH
13#define DUMUX_ONEEQ_PROBLEM_HH
33template<
class TypeTag>
40 using Grid =
typename GridView::Grid;
43 using DimVector = Dune::FieldVector<Scalar, Grid::dimension>;
47 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
56 RANSProblemImpl(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
65 ParentType::updateStaticWallProperties();
68 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
69 storedViscosityTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
70 storedViscosityTildeGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
78 template<
class SolutionVector>
81 ParentType::updateDynamicWallProperties(curSol);
83 auto fvGeometry =
localView(this->gridGeometry());
84 for (
const auto& element : elements(this->gridGeometry().gridView()))
86 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
87 fvGeometry.bindElement(element);
89 for (
auto&& scv : scvs(fvGeometry))
91 const int dofIdx = scv.dofIndex();
92 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
93 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
94 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
96 storedViscosityTilde_[elementIdx] = elemSol[0][Indices::viscosityTildeIdx];
98 VolumeVariables volVars;
99 volVars.update(elemSol, asImp_(), element, scv);
100 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
105 for (
const auto& element : elements(this->gridGeometry().gridView()))
107 const unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
108 fvGeometry.bindElement(element);
110 for (
unsigned int axisIdx = 0; axisIdx < Grid::dimension; ++axisIdx)
112 const unsigned int neighborIndex0 = ParentType::neighborIndex(elementIdx, axisIdx, 0);
113 const unsigned int neighborIndex1 = ParentType::neighborIndex(elementIdx, axisIdx, 1);
116 storedViscosityTildeGradient_[elementIdx][axisIdx]
117 = (storedViscosityTilde(neighborIndex1) - storedViscosityTilde(neighborIndex0))
118 / (ParentType::cellCenter(neighborIndex1)[axisIdx] - ParentType::cellCenter(neighborIndex0)[axisIdx]);
122 for (
auto&& scvf : scvfs(fvGeometry))
124 const unsigned int normDim = scvf.directionIndex();
125 if (scvf.boundary() && asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::viscosityTildeIdx))
128 Scalar dirichletViscosityTilde = asImp_().dirichlet(element, scvf)[Indices::viscosityTildeIdx];
130 unsigned int neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 0);
131 if (scvf.center()[normDim] < ParentType::cellCenter(elementIdx)[normDim])
132 neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 1);
134 storedViscosityTildeGradient_[elementIdx][normDim]
135 = (storedViscosityTilde(neighborIndex) - dirichletViscosityTilde)
136 / (ParentType::cellCenter(neighborIndex)[normDim] - scvf.center()[normDim]);
144 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
145 return useStoredEddyViscosity;
149 {
return storedDynamicEddyViscosity_[elementIdx]; }
152 {
return storedViscosityTilde_[elementIdx]; }
155 {
return storedViscosityTildeGradient_[elementIdx]; }
158 std::vector<Scalar> storedDynamicEddyViscosity_;
159 std::vector<Scalar> storedViscosityTilde_;
160 std::vector<DimVector> storedViscosityTildeGradient_;
163 Implementation &asImp_()
164 {
return *
static_cast<Implementation *
>(
this); }
167 const Implementation &asImp_()
const
168 {
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
void updateStaticWallProperties()
Correct size of the static (solution independent) wall variables.
Definition: freeflow/rans/oneeq/problem.hh:63
bool useStoredEddyViscosity() const
Definition: freeflow/rans/oneeq/problem.hh:142
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: freeflow/rans/oneeq/problem.hh:79
Scalar storedViscosityTilde(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:151
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:56
DimVector storedViscosityTildeGradient(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:154
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:148
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
Free function to get the local view of a grid cache object.
The available discretization methods in Dumux.
Adaption of the fully implicit scheme to the tracer transport model.
The local element solution class for staggered methods.
Base class for all staggered fv problems.
The available free flow turbulence models in Dumux.