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)
72 useStoredEddyViscosity_ = getParamFromGroup<bool>(this->paramGroup(),
73 "RANS.UseStoredEddyViscosity",
false);
81 ParentType::updateStaticWallProperties();
84 storedDynamicEddyViscosity_.resize(this->gridGeometry().elementMapper().size(), 0.0);
85 storedViscosityTilde_.resize(this->gridGeometry().elementMapper().size(), 0.0);
86 storedViscosityTildeGradient_.resize(this->gridGeometry().elementMapper().size(), DimVector(0.0));
96 ParentType::updateDynamicWallProperties(curSol);
98 for (
const auto& element : elements(this->gridGeometry().gridView()))
100 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
102 auto fvGeometry =
localView(this->gridGeometry());
103 fvGeometry.bindElement(element);
104 for (
auto&& scv : scvs(fvGeometry))
106 const int dofIdx = scv.dofIndex();
107 const auto& cellCenterPriVars = curSol[GridGeometry::cellCenterIdx()][dofIdx];
108 PrimaryVariables priVars = makePriVarsFromCellCenterPriVars<PrimaryVariables>(cellCenterPriVars);
109 auto elemSol = elementSolution<typename GridGeometry::LocalView>(std::move(priVars));
111 storedViscosityTilde_[elementIdx] = elemSol[0][Indices::viscosityTildeIdx];
113 VolumeVariables volVars;
114 volVars.update(elemSol, asImp_(), element, scv);
115 storedDynamicEddyViscosity_[elementIdx] = volVars.calculateEddyViscosity();
120 for (
const auto& element : elements(this->gridGeometry().gridView()))
122 unsigned int elementIdx = this->gridGeometry().elementMapper().index(element);
124 for (
unsigned int dimIdx = 0; dimIdx < Grid::dimension; ++dimIdx)
126 storedViscosityTildeGradient_[elementIdx][dimIdx]
127 = (storedViscosityTilde_[ParentType::neighborIdx_[elementIdx][dimIdx][1]]
128 - storedViscosityTilde_[ParentType::neighborIdx_[elementIdx][dimIdx][0]])
129 / (ParentType::cellCenter_[ParentType::neighborIdx_[elementIdx][dimIdx][1]][dimIdx]
130 - ParentType::cellCenter_[ParentType::neighborIdx_[elementIdx][dimIdx][0]][dimIdx]);
133 auto fvGeometry =
localView(this->gridGeometry());
134 fvGeometry.bindElement(element);
135 for (
auto&& scvf : scvfs(fvGeometry))
137 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 neighborIdx = ParentType::neighborIdx_[elementIdx][normDim][0];
144 if (scvf.center()[normDim] < ParentType::cellCenter_[elementIdx][normDim])
145 neighborIdx = ParentType::neighborIdx_[elementIdx][normDim][1];
147 storedViscosityTildeGradient_[elementIdx][normDim]
148 = (storedViscosityTilde_[neighborIdx] - dirichletViscosityTilde)
149 / (ParentType::cellCenter_[neighborIdx][normDim] - scvf.center()[normDim]);
163 Implementation &asImp_()
164 {
return *
static_cast<Implementation *
>(
this); }
167 const Implementation &asImp_()
const
168 {
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:37
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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: dumux/freeflow/rans/oneeq/problem.hh:79
void updateDynamicWallProperties(const SolutionVector &curSol)
Update the dynamic (solution dependent) relations to the walls.
Definition: dumux/freeflow/rans/oneeq/problem.hh:94
std::vector< Scalar > storedDynamicEddyViscosity_
Definition: dumux/freeflow/rans/oneeq/problem.hh:156
RANSProblemImpl(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor sets the gravity, if desired by the user.
Definition: dumux/freeflow/rans/oneeq/problem.hh:69
bool useStoredEddyViscosity_
Definition: dumux/freeflow/rans/oneeq/problem.hh:159
std::vector< DimVector > storedViscosityTildeGradient_
Definition: dumux/freeflow/rans/oneeq/problem.hh:158
std::vector< Scalar > storedViscosityTilde_
Definition: dumux/freeflow/rans/oneeq/problem.hh:157
forward declare
Definition: dumux/freeflow/rans/problem.hh:41
Reynolds-Averaged Navier-Stokes problem base class.
Definition: dumux/freeflow/rans/problem.hh:57
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.