24#ifndef DUMUX_NAVIERSTOKES_PROBLEM_HH
25#define DUMUX_NAVIERSTOKES_PROBLEM_HH
27#include <dune/common/exceptions.hh>
28#include <dune/common/typetraits.hh>
38template<
class TypeTag>
45template<
class TypeTag>
58template<
class TypeTag>
65 using GridView =
typename GridGeometry::GridView;
66 using Element =
typename GridView::template Codim<0>::Entity;
69 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
70 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
71 using FaceVariables =
typename GridFaceVariables::FaceVariables;
72 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
73 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
76 using FVElementGeometry =
typename GridGeometry::LocalView;
77 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
78 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
83 dim = GridView::dimension,
84 dimWorld = GridView::dimensionworld
87 using GlobalPosition =
typename SubControlVolumeFace::GlobalPosition;
88 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
89 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
97 NavierStokesProblem(std::shared_ptr<const GridGeometry> gridGeometry,
const std::string& paramGroup =
"")
98 : ParentType(gridGeometry, paramGroup)
101 if (getParamFromGroup<bool>(paramGroup,
"Problem.EnableGravity"))
102 gravity_[dim-1] = -9.81;
104 enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup,
"Problem.EnableInertiaTerms");
116 {
return asImp_().temperature(); }
124 { DUNE_THROW(Dune::NotImplemented,
"temperature() method not implemented by the actual problem"); }
139 {
return enableInertiaTerms_; }
142 template <
class SolutionVector,
class G = Gr
idGeometry>
143 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, void>::type
145 const SubControlVolumeFace& scvf,
146 const PrimaryVariables& initSol)
const
148 sol[GridGeometry::faceIdx()][scvf.dofIndex()][0] = initSol[Indices::velocity(scvf.directionIndex())];
167 const Scalar factor = 8.0)
const
169 static_assert(dim == 2,
"Pseudo 3D wall friction may only be used in 2D");
170 return -factor * velocity *
viscosity / (height*height);
174 template <
class ElementVolumeVariables,
class ElementFaceVariables,
class G = Gr
idGeometry>
175 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, Scalar>::type
177 const ElementVolumeVariables& elemVolVars,
178 const ElementFaceVariables& elemFaceVars,
180 const Scalar factor = 8.0)
const
182 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
183 const Scalar
viscosity = elemVolVars[scvf.insideScvIdx()].effectiveViscosity();
192 Scalar
permeability(
const Element& element,
const SubControlVolumeFace& scvf)
const
194 DUNE_THROW(Dune::NotImplemented,
"When using the Beavers-Joseph-Saffman boundary condition, the permeability must be returned in the acutal problem");
202 Scalar
alphaBJ(
const SubControlVolumeFace& scvf)
const
204 DUNE_THROW(Dune::NotImplemented,
"When using the Beavers-Joseph-Saffman boundary condition, the alpha value must be returned in the acutal problem");
210 Scalar
betaBJ(
const Element& element,
const SubControlVolumeFace& scvf,
const GlobalPosition& tangentialVector)
const
212 const Scalar interfacePermeability = interfacePermeability_(element, scvf, tangentialVector);
214 return asImp_().alphaBJ(scvf) / sqrt(interfacePermeability);
220 [[deprecated(
"Use betaBJ with tangential vector instead. Will be removed after 3.3")]]
221 Scalar
betaBJ(
const Element& element,
const SubControlVolumeFace& scvf)
const
224 return asImp_().alphaBJ(scvf) / sqrt(asImp_().
permeability(element, scvf));
232 return VelocityVector(0.0);
240 const SubControlVolume& scv,
241 const SubControlVolumeFace& ownScvf,
242 const SubControlVolumeFace& faceOnPorousBoundary,
243 const Scalar velocitySelf,
244 const Scalar tangentialVelocityGradient)
const
247 GlobalPosition orientation = ownScvf.unitOuterNormal();
248 orientation[ownScvf.directionIndex()] = 1.0;
252 const Scalar
betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary, orientation);
253 const Scalar distanceNormalToBoundary = (faceOnPorousBoundary.center() - scv.center()).two_norm();
255 return (tangentialVelocityGradient*distanceNormalToBoundary
256 + asImp_().porousMediumVelocity(element, faceOnPorousBoundary) * orientation *
betaBJ * distanceNormalToBoundary
257 + velocitySelf) / (
betaBJ*distanceNormalToBoundary + 1.0);
262 Scalar interfacePermeability_(
const Element& element,
const SubControlVolumeFace& scvf,
const GlobalPosition& tangentialVector)
const
264 const auto& K = asImp_().permeability(element, scvf);
267 if constexpr (Dune::IsNumber<std::decay_t<
decltype(K)>>::value)
270 return vtmv(tangentialVector, K, tangentialVector);
274 Implementation &asImp_()
275 {
return *
static_cast<Implementation *
>(
this); }
278 const Implementation &asImp_()
const
279 {
return *
static_cast<const Implementation *
>(
this); }
281 GravityVector gravity_;
282 bool enableInertiaTerms_;
Base class for all staggered fv problems.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:849
typename NavierStokesParentProblemImpl< TypeTag, GetPropType< TypeTag, Properties::GridGeometry >::discMethod >::type NavierStokesParentProblem
The actual NavierStokesParentProblem.
Definition: freeflow/navierstokes/problem.hh:48
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
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
Base class for all staggered finite-volume problems.
Definition: staggeredfvproblem.hh:47
The implementation is specialized for the different discretizations.
Definition: freeflow/navierstokes/problem.hh:36
Navier-Stokes problem base class.
Definition: freeflow/navierstokes/problem.hh:60
Scalar temperature() const
Returns the temperature within the domain.
Definition: freeflow/navierstokes/problem.hh:123
Scalar pseudo3DWallFriction(const Scalar velocity, const Scalar viscosity, const Scalar height, const Scalar factor=8.0) const
An additional drag term can be included as source term for the momentum balance to mimic 3D flow beha...
Definition: freeflow/navierstokes/problem.hh:164
Scalar permeability(const Element &element, const SubControlVolumeFace &scvf) const
Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boun...
Definition: freeflow/navierstokes/problem.hh:192
VelocityVector porousMediumVelocity(const Element &element, const SubControlVolumeFace &scvf) const
Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
Definition: freeflow/navierstokes/problem.hh:230
std::enable_if< G::discMethod==DiscretizationMethod::staggered, void >::type applyInitialFaceSolution(SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const
Applys the initial face solution (velocities on the faces). Specialization for staggered grid discret...
Definition: freeflow/navierstokes/problem.hh:144
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: freeflow/navierstokes/problem.hh:202
const GravityVector & gravity() const
Returns the acceleration due to gravity.
Definition: freeflow/navierstokes/problem.hh:132
std::enable_if< G::discMethod==DiscretizationMethod::staggered, Scalar >::type pseudo3DWallFriction(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const Scalar height, const Scalar factor=8.0) const
Convenience function for staggered grid implementation.
Definition: freeflow/navierstokes/problem.hh:176
Scalar betaBJ(const Element &element, const SubControlVolumeFace &scvf, const GlobalPosition &tangentialVector) const
Returns the beta value which is the alpha value divided by the square root of the (scalar-valued) int...
Definition: freeflow/navierstokes/problem.hh:210
const Scalar beaversJosephVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &ownScvf, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf, const Scalar tangentialVelocityGradient) const
Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition.
Definition: freeflow/navierstokes/problem.hh:239
Scalar betaBJ(const Element &element, const SubControlVolumeFace &scvf) const
Returns the beta value which is the alpha value divided by the square root of the (scalar-valued) int...
Definition: freeflow/navierstokes/problem.hh:221
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature at a given global position.
Definition: freeflow/navierstokes/problem.hh:115
bool enableInertiaTerms() const
Returns whether interia terms should be considered.
Definition: freeflow/navierstokes/problem.hh:138
NavierStokesProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string ¶mGroup="")
The constructor.
Definition: freeflow/navierstokes/problem.hh:97
Declares all properties used in Dumux.