24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH
52template<
class TypeTag,
class DiscretizationMethod>
53class NavierStokesFluxVariablesImpl;
55template<
class TypeTag>
58 typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView,
59 typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView,
60 typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView>
64 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
65 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
66 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
68 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
69 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
71 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
72 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
73 using FaceVariables =
typename GridFaceVariables::FaceVariables;
77 using FVElementGeometry =
typename GridGeometry::LocalView;
78 using GridView =
typename GridGeometry::GridView;
79 using Element =
typename GridView::template Codim<0>::Entity;
82 using Indices =
typename ModelTraits::Indices;
83 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
84 using FaceFrontalSubControlVolumeFace =
typename GridGeometry::Traits::FaceFrontalSubControlVolumeFace;
85 using FaceLateralSubControlVolumeFace =
typename GridGeometry::Traits::FaceLateralSubControlVolumeFace;
92 static constexpr bool normalizePressure = getPropValue<TypeTag, Properties::NormalizePressure>();
94 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
96 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
97 static constexpr auto faceIdx = GridGeometry::faceIdx();
99 static constexpr int upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
100 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
114 template<
class UpwindTerm>
115 [[deprecated(
"Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
117 const ElementVolumeVariables& elemVolVars,
118 const ElementFaceVariables& elemFaceVars,
119 const SubControlVolumeFace &scvf,
120 UpwindTerm upwindTerm)
122 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
123 const bool insideIsUpstream = scvf.directionSign() ==
sign(velocity);
124 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
126 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
127 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
129 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
130 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
132 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
133 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
134 * velocity * Extrusion::area(scvf) * scvf.directionSign();
136 return flux * extrusionFactor_(elemVolVars, scvf);
148 template<
class UpwindTerm>
150 const FVElementGeometry& fvGeometry,
151 const ElementVolumeVariables& elemVolVars,
152 const ElementFaceVariables& elemFaceVars,
153 const SubControlVolumeFace &scvf,
154 UpwindTerm upwindTerm)
156 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
157 const bool insideIsUpstream = scvf.directionSign() ==
sign(velocity);
158 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
160 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
161 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
163 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
164 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
166 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
167 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
168 * velocity * Extrusion::area(fvGeometry, scvf) * scvf.directionSign();
170 return flux * extrusionFactor_(elemVolVars, scvf);
189 const Element& element,
190 const FVElementGeometry& fvGeometry,
191 const ElementVolumeVariables& elemVolVars,
192 const ElementFaceVariables& elemFaceVars,
193 const SubControlVolumeFace& scvf,
194 const FluxVariablesCache& fluxVarsCache)
197 auto upwindTerm = [](
const auto& volVars) {
return volVars.density(); };
200 CellCenterPrimaryVariables result(0.0);
201 result[Indices::conti0EqIdx - ModelTraits::dim()] = advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTerm);
210 const Element& element,
211 const SubControlVolumeFace& scvf,
212 const FVElementGeometry& fvGeometry,
213 const ElementVolumeVariables& elemVolVars,
214 const ElementFaceVariables& elemFaceVars,
215 const GridFluxVariablesCache& gridFluxVarsCache)
217 return computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache) +
218 computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache);
239 const Element& element,
240 const SubControlVolumeFace& scvf,
241 const FVElementGeometry& fvGeometry,
242 const ElementVolumeVariables& elemVolVars,
243 const ElementFaceVariables& elemFaceVars,
244 const GridFluxVariablesCache& gridFluxVarsCache)
246 FacePrimaryVariables frontalFlux(0.0);
249 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
250 const Scalar velocityOpposite = elemFaceVars[scvf].velocityOpposite();
251 const auto& faceVars = elemFaceVars[scvf];
254 if (problem.enableInertiaTerms())
257 const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5;
258 const bool selfIsUpstream = scvf.directionSign() !=
sign(transportingVelocity);
262 * transportingVelocity * -1.0 * scvf.directionSign();
267 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
270 const Scalar velocityGrad_ii = VelocityGradients::velocityGradII(scvf, faceVars) * scvf.directionSign();
272 static const bool enableUnsymmetrizedVelocityGradient
273 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
274 const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0;
275 frontalFlux -= factor * insideVolVars.effectiveViscosity() * velocityGrad_ii;
281 const Scalar
pressure = normalizePressure ?
282 insideVolVars.pressure() - problem.initial(scvf)[Indices::pressureIdx]
283 : insideVolVars.pressure();
287 frontalFlux +=
pressure * -1.0 * scvf.directionSign();
291 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
292 FaceFrontalSubControlVolumeFace frontalFace(scv.center(), scvf.area());
293 return frontalFlux * Extrusion::area(fvGeometry, frontalFace) * insideVolVars.extrusionFactor();
314 const Element& element,
315 const SubControlVolumeFace& scvf,
316 const FVElementGeometry& fvGeometry,
317 const ElementVolumeVariables& elemVolVars,
318 const ElementFaceVariables& elemFaceVars,
319 const GridFluxVariablesCache& gridFluxVarsCache)
321 FacePrimaryVariables lateralFlux(0.0);
322 const auto& faceVars = elemFaceVars[scvf];
323 const std::size_t numSubFaces = scvf.pairData().size();
327 std::optional<BoundaryTypes> currentScvfBoundaryTypes;
329 currentScvfBoundaryTypes.emplace(problem.boundaryTypes(element, scvf));
332 for (
int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
334 const auto eIdx = scvf.insideScvIdx();
336 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
339 std::optional<BoundaryTypes> lateralFaceBoundaryTypes;
343 if (lateralFace.boundary())
356 lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralFace));
362 if (currentScvfBoundaryTypes)
365 if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralFace.directionIndex())))
368 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
369 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
371 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())]
372 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf) * lateralFace.directionSign();
381 if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())) && lateralFaceBoundaryTypes &&
382 lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
384 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
385 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
387 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
388 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf) * lateralFace.directionSign();
396 if (lateralFace.boundary())
400 if (lateralFaceBoundaryTypes->isSymmetry())
404 if (lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
408 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
409 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
411 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
412 * elemVolVars[lateralFace.insideScvIdx()].extrusionFactor() * Extrusion::area(fvGeometry, lateralScvf);
418 if (incorporateWallFunction_(lateralFlux, problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, localSubFaceIdx))
423 if (lateralFaceBoundaryTypes)
425 std::bitset<3> admittableBcTypes;
426 admittableBcTypes.set(0, lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
427 admittableBcTypes.set(1, lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())));
428 admittableBcTypes.set(2, lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())));
429 if (admittableBcTypes.count() != 1)
431 DUNE_THROW(Dune::InvalidStateException,
"Invalid boundary conditions for lateral scvf "
432 "for the momentum equations at global position " << lateralStaggeredFaceCenter_(scvf, localSubFaceIdx)
433 <<
", current scvf global position " << scvf.center());
438 if (problem.enableInertiaTerms())
439 lateralFlux += computeAdvectivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
440 scvf, elemVolVars, elemFaceVars,
442 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
445 lateralFlux += computeDiffusivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
446 scvf, elemVolVars, faceVars,
447 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
469 [[deprecated(
"Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
471 const Element& element,
472 const SubControlVolumeFace& scvf,
473 const ElementVolumeVariables& elemVolVars,
474 const ElementFaceVariables& elemFaceVars)
const
476 FacePrimaryVariables inOrOutflow(0.0);
477 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
480 if (problem.enableInertiaTerms())
482 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
483 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
484 const auto& upVolVars = (scvf.directionSign() ==
sign(velocitySelf)) ?
485 insideVolVars : outsideVolVars;
487 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
491 const Scalar boundaryPressure = normalizePressure
492 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
493 problem.initial(scvf)[Indices::pressureIdx])
494 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
495 inOrOutflow += boundaryPressure;
498 return inOrOutflow * scvf.directionSign() * Extrusion::area(scvf) * insideVolVars.extrusionFactor();
518 const SubControlVolumeFace& scvf,
519 const FVElementGeometry& fvGeometry,
520 const ElementVolumeVariables& elemVolVars,
521 const ElementFaceVariables& elemFaceVars)
const
523 FacePrimaryVariables inOrOutflow(0.0);
524 const auto& element = fvGeometry.element();
525 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
528 if (problem.enableInertiaTerms())
530 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
531 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
532 const auto& upVolVars = (scvf.directionSign() ==
sign(velocitySelf)) ?
533 insideVolVars : outsideVolVars;
535 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
539 const Scalar boundaryPressure = normalizePressure
540 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
541 problem.initial(scvf)[Indices::pressureIdx])
542 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
543 inOrOutflow += boundaryPressure;
546 return inOrOutflow * scvf.directionSign() * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
572 FacePrimaryVariables computeAdvectivePartOfLateralMomentumFlux_(
const Problem& problem,
573 const FVElementGeometry& fvGeometry,
574 const Element& element,
575 const SubControlVolumeFace& scvf,
576 const ElementVolumeVariables& elemVolVars,
577 const ElementFaceVariables& elemFaceVars,
578 const GridFluxVariablesCache& gridFluxVarsCache,
579 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
580 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
581 const int localSubFaceIdx)
583 const auto eIdx = scvf.insideScvIdx();
584 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
588 const Scalar transportingVelocity = [&]()
590 const auto& faceVars = elemFaceVars[scvf];
591 if (!scvf.boundary())
592 return faceVars.velocityLateralInside(localSubFaceIdx);
596 const auto bcTypes = problem.boundaryTypes(element, scvf);
598 if (bcTypes.isDirichlet(Indices::velocity(lateralFace.directionIndex())))
602 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
604 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())];
606 else if (bcTypes.isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
608 return VelocityGradients::beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
609 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
612 return faceVars.velocityLateralInside(localSubFaceIdx);
616 const bool selfIsUpstream = lateralFace.directionSign() ==
sign(transportingVelocity);
617 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
618 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
619 return upwindHelper.computeUpwindLateralMomentum(selfIsUpstream, lateralFace, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes)
620 * transportingVelocity * lateralFace.directionSign() * Extrusion::area(fvGeometry, lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
645 FacePrimaryVariables computeDiffusivePartOfLateralMomentumFlux_(
const Problem& problem,
646 const FVElementGeometry& fvGeometry,
647 const Element& element,
648 const SubControlVolumeFace& scvf,
649 const ElementVolumeVariables& elemVolVars,
650 const FaceVariables& faceVars,
651 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
652 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
653 const int localSubFaceIdx)
655 const auto eIdx = scvf.insideScvIdx();
656 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
658 FacePrimaryVariables lateralDiffusiveFlux(0.0);
660 static const bool enableUnsymmetrizedVelocityGradient
661 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
666 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
667 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
670 const Scalar muAvg = lateralFace.boundary()
671 ? insideVolVars.effectiveViscosity()
672 : (insideVolVars.effectiveViscosity() + outsideVolVars.effectiveViscosity()) * 0.5;
675 if (!enableUnsymmetrizedVelocityGradient)
677 if (!scvf.boundary() ||
678 currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralFace.directionIndex())) ||
679 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
681 const Scalar velocityGrad_ji = VelocityGradients::velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
683 lateralDiffusiveFlux -= muAvg * velocityGrad_ji * lateralFace.directionSign();
690 if (!lateralFace.boundary() || !lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx))
692 const Scalar velocityGrad_ij = VelocityGradients::velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
693 lateralDiffusiveFlux -= muAvg * velocityGrad_ij * lateralFace.directionSign();
697 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
698 return lateralDiffusiveFlux * Extrusion::area(fvGeometry, lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
715 const GlobalPosition& lateralStaggeredFaceCenter_(
const SubControlVolumeFace& scvf,
const int localSubFaceIdx)
const
717 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
734 GlobalPosition lateralStaggeredSCVFCenter_(
const SubControlVolumeFace& lateralFace,
735 const SubControlVolumeFace& currentFace,
736 const int localSubFaceIdx)
const
738 auto pos = lateralStaggeredFaceCenter_(currentFace, localSubFaceIdx) + lateralFace.center();
744 static Scalar extrusionFactor_(
const ElementVolumeVariables& elemVolVars,
const SubControlVolumeFace& scvf)
746 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
747 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
748 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
752 template<
class ...Args,
bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<!turbulenceModel, int> = 0>
753 bool incorporateWallFunction_(Args&&... args)
const
757 template<
bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<turbulenceModel,
int> = 0>
758 bool incorporateWallFunction_(FacePrimaryVariables& lateralFlux,
759 const Problem& problem,
760 const Element& element,
761 const FVElementGeometry& fvGeometry,
762 const SubControlVolumeFace& scvf,
763 const ElementVolumeVariables& elemVolVars,
764 const ElementFaceVariables& elemFaceVars,
765 const std::size_t localSubFaceIdx)
const
767 const auto eIdx = scvf.insideScvIdx();
768 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
770 if (problem.useWallFunction(element, lateralFace, Indices::velocity(scvf.directionIndex())))
772 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
773 const auto lateralBoundaryFace =
makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter_(scvf, localSubFaceIdx));
774 lateralFlux += problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
775 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf);
Some exceptions thrown in DuMux
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.
Base class for the flux variables living on a sub control volume face.
SubControlVolumeFace makeStaggeredBoundaryFace(const SubControlVolumeFace &scvf, const typename SubControlVolumeFace::GlobalPosition &newCenter)
Helper function to turn a given cell scvface into a fake boundary face.
Definition: discretization/staggered/freeflow/subcontrolvolumeface.hh:104
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
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
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:46
Base class for the flux variables living on a sub control volume face.
Definition: fluxvariablesbase.hh:45
The flux variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: freeflow/navierstokes/fluxvariables.hh:35
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:58
static Scalar advectiveFluxForCellCenter(const Problem &problem, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, UpwindTerm upwindTerm)
Returns the advective flux over a sub control volume face.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:149
FacePrimaryVariables computeLateralMomentumFlux(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
Returns the momentum flux over the staggered faces perpendicular to the scvf where the velocity dof o...
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:313
FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Returns the momentum flux over an inflow or outflow boundary face.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:517
FacePrimaryVariables computeFrontalMomentumFlux(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
Returns the frontal part of the momentum flux. This treats the flux over the staggered face at the ce...
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:238
static Scalar advectiveFluxForCellCenter(const Problem &problem, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, UpwindTerm upwindTerm)
Returns the advective flux over a sub control volume face.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:116
FacePrimaryVariables computeMomentumFlux(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const GridFluxVariablesCache &gridFluxVarsCache)
Returns the momentum flux over all staggered faces.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:209
FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem &problem, const Element &element, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Returns the momentum flux over an inflow or outflow boundary face.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:470
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:104
CellCenterPrimaryVariables computeMassFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf, const FluxVariablesCache &fluxVarsCache)
Computes the flux for the cell center residual (mass balance).
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:188
The upwinding variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: staggeredupwindhelper.hh:48
FacePrimaryVariables computeUpwindFrontalMomentum(const bool selfIsUpstream) const
Returns the momentum in the frontal direction.
Definition: staggeredupwindhelper.hh:115
Declares all properties used in Dumux.
Type traits for problem classes.