24#ifndef DUMUX_NAVIERSTOKES_MOMENTUM_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_MOMENTUM_FLUXVARIABLES_HH
47template<
class TypeTag>
52 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
53 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
54 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
56 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
57 using ElementFluxVariablesCache =
typename GridFluxVariablesCache::LocalView;
58 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
60 using GridGeometry =
typename GridVariables::GridGeometry;
61 using FVElementGeometry =
typename GridGeometry::LocalView;
62 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
63 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
68 using GridView =
typename GridGeometry::GridView;
69 using Element =
typename GridView::template Codim<0>::Entity;
72 using Indices =
typename ModelTraits::Indices;
75 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
85 const SubControlVolumeFace&
scvFace,
99 ==
static_cast<std::size_t
>(GridView::dimension),
100 "Expects problem.dirichlet to return an array with as many entries as dimensions."
105 {
return *problemPtr_; }
108 {
return *elementPtr_; }
111 {
return *scvFacePtr_; }
114 {
return *fvGeometryPtr_; }
117 {
return *elemVolVarsPtr_; }
120 {
return *elemFluxVarsCachePtr_; }
123 {
return *elemBcTypesPtr_; }
130 if (!this->
problem().enableInertiaTerms())
131 return NumEqVector(0.0);
133 if (this->
scvFace().isFrontal())
144 if (this->
scvFace().isFrontal())
170 const auto& scvf = this->
scvFace();
171 assert(scvf.isFrontal());
173 NumEqVector result(0.0);
183 static const bool enableUnsymmetrizedVelocityGradient
184 = getParamFromGroup<bool>(this->
problem().paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
185 static const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0;
188 result -= factor * mu * velGradII * Extrusion::area(
fvGeometry, scvf) *
elemVolVars[scvf.insideScvIdx()].extrusionFactor() * scvf.directionSign();
190 static const bool enableDilatationTerm = getParamFromGroup<bool>(this->
problem().paramGroup(),
"FreeFlow.EnableDilatationTerm",
false);
191 if (enableDilatationTerm)
193 Scalar divergence = velGradII;
196 const auto otherFrontalScvf = *(scvfs(
fvGeometry, scv).begin());
197 assert(otherFrontalScvf.isFrontal() && !otherFrontalScvf.boundary());
198 if (otherFrontalScvf.index() != scvf.index())
202 result += 2.0/3.0 * mu * divergence * scvf.directionSign() * Extrusion::area(
fvGeometry, scvf) *
elemVolVars[scvf.insideScvIdx()].extrusionFactor();
232 const auto& scvf = this->
scvFace();
233 assert(scvf.isLateral());
235 NumEqVector result(0.0);
239 const auto& scv =
fvGeometry.scv(scvf.insideScvIdx());
241 static const bool enableUnsymmetrizedVelocityGradient
242 = getParamFromGroup<bool>(
problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
250 GlobalPosition gradVn(0.0);
251 gradV.mv(scvf.unitOuterNormal(), gradVn);
252 const Scalar velocityGrad_ij = gradVn[scv.dofAxis()];
253 result -= mu * velocityGrad_ij;
256 if (!enableUnsymmetrizedVelocityGradient)
258 GlobalPosition gradVTransposedN(0.0);
259 gradV.mtv(scvf.unitOuterNormal(), gradVTransposedN);
260 const Scalar velocityGrad_ji = gradVTransposedN[scv.dofAxis()];
261 result -= mu * velocityGrad_ji;
265 return result * Extrusion::area(
fvGeometry, scvf) *
elemVolVars[scvf.insideScvIdx()].extrusionFactor();
286 NumEqVector result(0.0);
287 const auto& scvf = this->
scvFace();
288 if (scvf.isLateral() || scvf.boundary())
304 result -= referencePressure*scvf.area();
307 result *= scvf.directionSign();
331 NumEqVector flux(0.0);
332 const auto& scvf = this->
scvFace();
333 assert(scvf.isFrontal());
337 const auto velocitySelf =
elemVolVars[scvf.insideScvIdx()].velocity();
338 const auto velocityOpposite =
elemVolVars[scvf.outsideScvIdx()].velocity();
341 const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5;
343 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
346 static const auto upwindWeight = getParamFromGroup<Scalar>(
problem.paramGroup(),
"Flux.UpwindWeight");
347 const Scalar transportedMomentum = selfIsUpstream ? (upwindWeight * velocitySelf + (1.0 - upwindWeight) * velocityOpposite) *
density
348 : (upwindWeight * velocityOpposite + (1.0 - upwindWeight) * velocitySelf) *
density;
350 return transportingVelocity * transportedMomentum * scvf.directionSign() * Extrusion::area(this->
fvGeometry(), scvf) * extrusionFactor_(
elemVolVars, scvf);
377 NumEqVector flux(0.0);
378 const auto& scvf = this->
scvFace();
379 assert(scvf.isLateral());
386 const Scalar transportingVelocity = [&]()
388 const auto& orthogonalScvf =
fvGeometry.lateralOrthogonalScvf(scvf);
389 const Scalar innerTransportingVelocity =
elemVolVars[orthogonalScvf.insideScvIdx()].velocity();
391 static const bool useOldScheme = getParam<bool>(
"FreeFlow.UseOldTransportingVelocity",
true);
395 if (scvf.boundary() &&
fvGeometry.scv(scvf.insideScvIdx()).boundary())
397 if (this->
elemBcTypes()[scvf.localIndex()].isDirichlet(scvf.normalAxis()))
398 return problem.dirichlet(this->
element(), scvf)[scvf.normalAxis()];
401 return innerTransportingVelocity;
407 if (this->
elemBcTypes()[scvf.localIndex()].isDirichlet(scvf.normalAxis()))
408 return 0.5*(
problem.dirichlet(this->
element(), scvf)[scvf.normalAxis()] + innerTransportingVelocity);
411 if (orthogonalScvf.boundary())
413 if (this->
elemBcTypes()[orthogonalScvf.localIndex()].isDirichlet(scvf.normalAxis()))
414 return 0.5*(
problem.dirichlet(this->element(), scvf)[scvf.normalAxis()] + innerTransportingVelocity);
416 return innerTransportingVelocity;
420 const auto insideVolume =
fvGeometry.scv(orthogonalScvf.insideScvIdx()).
volume();
421 const auto outsideVolume =
fvGeometry.scv(orthogonalScvf.outsideScvIdx()).
volume();
422 const auto outerTransportingVelocity =
elemVolVars[orthogonalScvf.outsideScvIdx()].velocity();
423 return (insideVolume*innerTransportingVelocity + outsideVolume*outerTransportingVelocity) / (insideVolume + outsideVolume);
426 const Scalar transportedMomentum = [&]()
428 const auto& insideScv =
fvGeometry.scv(scvf.insideScvIdx());
430 auto getDirichletMomentumFlux = [&]()
438 if (!this->
elemBcTypes()[scvf.localIndex()].isDirichlet(insideScv.dofAxis()))
439 DUNE_THROW(Dune::InvalidStateException,
"Neither Dirichlet nor Neumann BC set at " << scvf.ipGlobal());
441 return getDirichletMomentumFlux();
445 if (
fvGeometry.scvfIntegrationPointInConcaveCorner(scvf))
448 const auto& outsideScvfWithSameIntegrationPoint =
fvGeometry.outsideScvfWithSameIntegrationPoint(scvf);
449 if (!this->
problem().boundaryTypes(this->
element(), outsideScvfWithSameIntegrationPoint).isDirichlet(insideScv.dofAxis()))
450 DUNE_THROW(Dune::InvalidStateException,
"Neither Dirichlet nor Neumann BC set at " << outsideScvfWithSameIntegrationPoint.ipGlobal());
452 return getDirichletMomentumFlux();
456 const bool selfIsUpstream = scvf.directionSign() ==
sign(transportingVelocity);
458 const auto innerVelocity =
elemVolVars[scvf.insideScvIdx()].velocity();
459 const auto outerVelocity =
elemVolVars[scvf.outsideScvIdx()].velocity();
463 const auto insideMomentum = innerVelocity * rho.first;
464 const auto outsideMomentum = outerVelocity * rho.second;
467 static const auto upwindWeight = getParamFromGroup<Scalar>(
problem.paramGroup(),
"Flux.UpwindWeight");
469 return selfIsUpstream ? (upwindWeight * insideMomentum + (1.0 - upwindWeight) * outsideMomentum)
470 : (upwindWeight * outsideMomentum + (1.0 - upwindWeight) * insideMomentum);
473 return transportingVelocity * transportedMomentum * scvf.directionSign() * Extrusion::area(
fvGeometry, scvf) * extrusionFactor_(
elemVolVars, scvf);
478 template<
class ElementVolumeVariables,
class SubControlVolumeFace>
479 Scalar extrusionFactor_(
const ElementVolumeVariables&
elemVolVars,
const SubControlVolumeFace& scvf)
const
481 const auto& insideVolVars =
elemVolVars[scvf.insideScvIdx()];
482 const auto& outsideVolVars =
elemVolVars[scvf.outsideScvIdx()];
483 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
487 const Problem* problemPtr_;
488 const Element* elementPtr_;
489 const FVElementGeometry* fvGeometryPtr_;
490 const SubControlVolumeFace* scvFacePtr_;
491 const ElementVolumeVariables* elemVolVarsPtr_;
492 const ElementFluxVariablesCache* elemFluxVarsCachePtr_;
493 const ElementBoundaryTypes* elemBcTypesPtr_;
Some exceptions thrown in DuMux
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
constexpr Scalar harmonicMean(Scalar x, Scalar y, Scalar wx=1.0, Scalar wy=1.0) noexcept
Calculate the (weighted) harmonic mean of two scalar values.
Definition: math.hh:69
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
The flux variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:49
NumEqVector pressureContribution() const
Returns the frontal pressure contribution.
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:284
NumEqVector diffusiveMomentumFlux() const
Returns the diffusive momentum flux due to viscous forces.
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:142
NumEqVector frontalAdvectiveMomentumFlux() const
Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the ce...
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:329
NumEqVector lateralAdvectiveMomentumFlux() const
Returns the advective momentum flux over the staggered face perpendicular to the scvf where the veloc...
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:375
NumEqVector frontalDiffusiveMomentumFlux() const
Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the ce...
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:168
const ElementBoundaryTypes & elemBcTypes() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:122
const Problem & problem() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:104
const SubControlVolumeFace & scvFace() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:110
const FVElementGeometry & fvGeometry() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:113
NumEqVector advectiveMomentumFlux() const
Returns the diffusive momentum flux due to viscous forces.
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:128
NavierStokesMomentumFluxVariables(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvFace, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const ElementBoundaryTypes &elemBcTypes)
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:82
const ElementVolumeVariables & elemVolVars() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:116
const ElementFluxVariablesCache & elemFluxVarsCache() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:119
const Element & element() const
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:107
NumEqVector lateralDiffusiveMomentumFlux() const
Returns the diffusive momentum flux over the staggered face perpendicular to the scvf where the veloc...
Definition: freeflow/navierstokes/momentum/fluxvariables.hh:230
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:58
static auto velocityGradient(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars, bool fullGradient=false)
Definition: momentum/velocitygradients.hh:62
static auto velocityGradII(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const ElemVolVars &elemVolVars)
Returns the in-axis velocity gradient.
Definition: momentum/velocitygradients.hh:158
Declares all properties used in Dumux.