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 =
"")
70 : ParentType(gridGeometry, 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 for (
const auto& element : elements(this->gridGeometry().gridView()))
97 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
99 auto fvGeometry =
localView(this->gridGeometry());
100 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);
121 for (
unsigned int dimIdx = 0; dimIdx < Grid::dimension; ++dimIdx)
123 const unsigned int neighborIndex0 = ParentType::neighborIndex(elementIdx, dimIdx, 0);
124 const unsigned int neighborIndex1 = ParentType::neighborIndex(elementIdx, dimIdx, 1);
127 storedViscosityTildeGradient_[elementIdx][dimIdx]
128 = (storedViscosityTilde(neighborIndex1) - storedViscosityTilde(neighborIndex0))
129 / (ParentType::cellCenter(neighborIndex1)[dimIdx] - ParentType::cellCenter(neighborIndex0)[dimIdx]);
133 auto fvGeometry =
localView(this->gridGeometry());
134 fvGeometry.bindElement(element);
135 for (
auto&& scvf : scvfs(fvGeometry))
137 const unsigned int normDim = scvf.directionIndex();
138 if (scvf.boundary() && asImp_().boundaryTypes(element, scvf).isDirichlet(Indices::viscosityTildeIdx))
141 Scalar dirichletViscosityTilde = asImp_().dirichlet(element, scvf)[Indices::viscosityTildeIdx];
143 unsigned int neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 0);
144 if (scvf.center()[normDim] < ParentType::cellCenter(elementIdx)[normDim])
145 neighborIndex = ParentType::neighborIndex(elementIdx, normDim, 1);
147 storedViscosityTildeGradient_[elementIdx][normDim]
148 = (storedViscosityTilde(neighborIndex) - dirichletViscosityTilde)
149 / (ParentType::cellCenter(neighborIndex)[normDim] - scvf.center()[normDim]);
157 static const bool useStoredEddyViscosity = getParamFromGroup<bool>(this->paramGroup(),
"RANS.UseStoredEddyViscosity",
false);
158 return useStoredEddyViscosity;
162 {
return storedDynamicEddyViscosity_[elementIdx]; }
165 {
return storedViscosityTilde_[elementIdx]; }
168 {
return storedViscosityTildeGradient_[elementIdx]; }
171 std::vector<Scalar> storedDynamicEddyViscosity_;
172 std::vector<Scalar> storedViscosityTilde_;
173 std::vector<DimVector> storedViscosityTildeGradient_;
176 Implementation &asImp_()
177 {
return *
static_cast<Implementation *
>(
this); }
180 const Implementation &asImp_()
const
181 {
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
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:155
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:164
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:167
Scalar storedDynamicEddyViscosity(const int elementIdx) const
Definition: freeflow/rans/oneeq/problem.hh:161
forward declare
Definition: freeflow/rans/problem.hh:43
Reynolds-Averaged Navier-Stokes problem base class.
Definition: freeflow/rans/problem.hh:59
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.