24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH
52template<
class TypeTag, DiscretizationMethod discMethod>
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>
116 const ElementVolumeVariables& elemVolVars,
117 const ElementFaceVariables& elemFaceVars,
118 const SubControlVolumeFace &scvf,
119 UpwindTerm upwindTerm)
121 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
122 const bool insideIsUpstream = scvf.directionSign() ==
sign(velocity);
123 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(),
"Flux.UpwindWeight");
125 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
126 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
128 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
129 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
131 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
132 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
133 * velocity * Extrusion::area(scvf) * scvf.directionSign();
135 return flux * extrusionFactor_(elemVolVars, scvf);
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const ElementFaceVariables& elemFaceVars,
158 const SubControlVolumeFace& scvf,
159 const FluxVariablesCache& fluxVarsCache)
162 auto upwindTerm = [](
const auto& volVars) {
return volVars.density(); };
165 CellCenterPrimaryVariables result(0.0);
166 result[Indices::conti0EqIdx - ModelTraits::dim()] = advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTerm);
175 const Element& element,
176 const SubControlVolumeFace& scvf,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const ElementFaceVariables& elemFaceVars,
180 const GridFluxVariablesCache& gridFluxVarsCache)
182 return computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache) +
183 computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache);
204 const Element& element,
205 const SubControlVolumeFace& scvf,
206 const FVElementGeometry& fvGeometry,
207 const ElementVolumeVariables& elemVolVars,
208 const ElementFaceVariables& elemFaceVars,
209 const GridFluxVariablesCache& gridFluxVarsCache)
211 FacePrimaryVariables frontalFlux(0.0);
214 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
215 const Scalar velocityOpposite = elemFaceVars[scvf].velocityOpposite();
216 const auto& faceVars = elemFaceVars[scvf];
219 if (problem.enableInertiaTerms())
222 const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5;
223 const bool selfIsUpstream = scvf.directionSign() !=
sign(transportingVelocity);
227 * transportingVelocity * -1.0 * scvf.directionSign();
232 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
235 const Scalar velocityGrad_ii = VelocityGradients::velocityGradII(scvf, faceVars) * scvf.directionSign();
237 static const bool enableUnsymmetrizedVelocityGradient
238 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
239 const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0;
240 frontalFlux -= factor * insideVolVars.effectiveViscosity() * velocityGrad_ii;
246 const Scalar
pressure = normalizePressure ?
247 insideVolVars.pressure() - problem.initial(scvf)[Indices::pressureIdx]
248 : insideVolVars.pressure();
252 frontalFlux +=
pressure * -1.0 * scvf.directionSign();
256 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
257 FaceFrontalSubControlVolumeFace frontalFace(scv.center(), scvf.area());
258 return frontalFlux * Extrusion::area(frontalFace) * insideVolVars.extrusionFactor();
279 const Element& element,
280 const SubControlVolumeFace& scvf,
281 const FVElementGeometry& fvGeometry,
282 const ElementVolumeVariables& elemVolVars,
283 const ElementFaceVariables& elemFaceVars,
284 const GridFluxVariablesCache& gridFluxVarsCache)
286 FacePrimaryVariables lateralFlux(0.0);
287 const auto& faceVars = elemFaceVars[scvf];
288 const std::size_t numSubFaces = scvf.pairData().size();
292 std::optional<BoundaryTypes> currentScvfBoundaryTypes;
294 currentScvfBoundaryTypes.emplace(problem.boundaryTypes(element, scvf));
297 for (
int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
299 const auto eIdx = scvf.insideScvIdx();
301 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
304 std::optional<BoundaryTypes> lateralFaceBoundaryTypes;
308 if (lateralFace.boundary())
321 lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralFace));
327 if (currentScvfBoundaryTypes)
330 if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralFace.directionIndex())))
333 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
334 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
336 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())]
337 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf) * lateralFace.directionSign();
346 if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())) && lateralFaceBoundaryTypes &&
347 lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
349 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
350 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
352 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
353 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf) * lateralFace.directionSign();
361 if (lateralFace.boundary())
365 if (lateralFaceBoundaryTypes->isSymmetry())
369 if (lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
373 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
374 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
376 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
377 * elemVolVars[lateralFace.insideScvIdx()].extrusionFactor() * Extrusion::area(lateralScvf);
383 if (incorporateWallFunction_(lateralFlux, problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, localSubFaceIdx))
388 if (lateralFaceBoundaryTypes)
390 std::bitset<3> admittableBcTypes;
391 admittableBcTypes.set(0, lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
392 admittableBcTypes.set(1, lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())));
393 admittableBcTypes.set(2, lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())));
394 if (admittableBcTypes.count() != 1)
396 DUNE_THROW(Dune::InvalidStateException,
"Invalid boundary conditions for lateral scvf "
397 "for the momentum equations at global position " << lateralStaggeredFaceCenter_(scvf, localSubFaceIdx)
398 <<
", current scvf global position " << scvf.center());
403 if (problem.enableInertiaTerms())
404 lateralFlux += computeAdvectivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
405 scvf, elemVolVars, elemFaceVars,
407 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
410 lateralFlux += computeDiffusivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
411 scvf, elemVolVars, faceVars,
412 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
435 const Element& element,
436 const SubControlVolumeFace& scvf,
437 const ElementVolumeVariables& elemVolVars,
438 const ElementFaceVariables& elemFaceVars)
const
440 FacePrimaryVariables inOrOutflow(0.0);
441 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
444 if (problem.enableInertiaTerms())
446 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
447 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
448 const auto& upVolVars = (scvf.directionSign() ==
sign(velocitySelf)) ?
449 insideVolVars : outsideVolVars;
451 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
455 const Scalar boundaryPressure = normalizePressure
456 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
457 problem.initial(scvf)[Indices::pressureIdx])
458 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
459 inOrOutflow += boundaryPressure;
462 return inOrOutflow * scvf.directionSign() * Extrusion::area(scvf) * insideVolVars.extrusionFactor();
488 FacePrimaryVariables computeAdvectivePartOfLateralMomentumFlux_(
const Problem& problem,
489 const FVElementGeometry& fvGeometry,
490 const Element& element,
491 const SubControlVolumeFace& scvf,
492 const ElementVolumeVariables& elemVolVars,
493 const ElementFaceVariables& elemFaceVars,
494 const GridFluxVariablesCache& gridFluxVarsCache,
495 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
496 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
497 const int localSubFaceIdx)
499 const auto eIdx = scvf.insideScvIdx();
500 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
504 const Scalar transportingVelocity = [&]()
506 const auto& faceVars = elemFaceVars[scvf];
507 if (!scvf.boundary())
508 return faceVars.velocityLateralInside(localSubFaceIdx);
512 const auto bcTypes = problem.boundaryTypes(element, scvf);
514 if (bcTypes.isDirichlet(Indices::velocity(lateralFace.directionIndex())))
518 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
520 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())];
522 else if (bcTypes.isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
524 return VelocityGradients::beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
525 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
528 return faceVars.velocityLateralInside(localSubFaceIdx);
532 const bool selfIsUpstream = lateralFace.directionSign() ==
sign(transportingVelocity);
533 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
534 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
535 return upwindHelper.computeUpwindLateralMomentum(selfIsUpstream, lateralFace, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes)
536 * transportingVelocity * lateralFace.directionSign() * Extrusion::area(lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
561 FacePrimaryVariables computeDiffusivePartOfLateralMomentumFlux_(
const Problem& problem,
562 const FVElementGeometry& fvGeometry,
563 const Element& element,
564 const SubControlVolumeFace& scvf,
565 const ElementVolumeVariables& elemVolVars,
566 const FaceVariables& faceVars,
567 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
568 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
569 const int localSubFaceIdx)
571 const auto eIdx = scvf.insideScvIdx();
572 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
574 FacePrimaryVariables lateralDiffusiveFlux(0.0);
576 static const bool enableUnsymmetrizedVelocityGradient
577 = getParamFromGroup<bool>(problem.paramGroup(),
"FreeFlow.EnableUnsymmetrizedVelocityGradient",
false);
582 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
583 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
586 const Scalar muAvg = lateralFace.boundary()
587 ? insideVolVars.effectiveViscosity()
588 : (insideVolVars.effectiveViscosity() + outsideVolVars.effectiveViscosity()) * 0.5;
591 if (!enableUnsymmetrizedVelocityGradient)
593 if (!scvf.boundary() ||
594 currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralFace.directionIndex())) ||
595 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
597 const Scalar velocityGrad_ji = VelocityGradients::velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
599 lateralDiffusiveFlux -= muAvg * velocityGrad_ji * lateralFace.directionSign();
606 if (!lateralFace.boundary() || !lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx))
608 const Scalar velocityGrad_ij = VelocityGradients::velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
609 lateralDiffusiveFlux -= muAvg * velocityGrad_ij * lateralFace.directionSign();
613 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
614 return lateralDiffusiveFlux * Extrusion::area(lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
631 const GlobalPosition& lateralStaggeredFaceCenter_(
const SubControlVolumeFace& scvf,
const int localSubFaceIdx)
const
633 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
650 GlobalPosition lateralStaggeredSCVFCenter_(
const SubControlVolumeFace& lateralFace,
651 const SubControlVolumeFace& currentFace,
652 const int localSubFaceIdx)
const
654 auto pos = lateralStaggeredFaceCenter_(currentFace, localSubFaceIdx) + lateralFace.center();
660 static Scalar extrusionFactor_(
const ElementVolumeVariables& elemVolVars,
const SubControlVolumeFace& scvf)
662 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
663 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
664 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
668 template<
class ...Args,
bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<!turbulenceModel, int> = 0>
669 bool incorporateWallFunction_(Args&&... args)
const
673 template<
bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<turbulenceModel,
int> = 0>
674 bool incorporateWallFunction_(FacePrimaryVariables& lateralFlux,
675 const Problem& problem,
676 const Element& element,
677 const FVElementGeometry& fvGeometry,
678 const SubControlVolumeFace& scvf,
679 const ElementVolumeVariables& elemVolVars,
680 const ElementFaceVariables& elemFaceVars,
681 const std::size_t localSubFaceIdx)
const
683 const auto eIdx = scvf.insideScvIdx();
684 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
686 if (problem.useWallFunction(element, lateralFace, Indices::velocity(scvf.directionIndex())))
688 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
689 const auto lateralBoundaryFace =
makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter_(scvf, localSubFaceIdx));
690 lateralFlux += problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
691 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf);
Some exceptions thrown in DuMux
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Helper classes to compute the integration elements.
Base class for the flux variables living on a sub control volume face.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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:68
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:627
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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 pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
typename Detail::template ProblemTraits< Problem, GridGeometry::discMethod >::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
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:278
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:203
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:115
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:174
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:104
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:434
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:153
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 directon.
Definition: staggeredupwindhelper.hh:115
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: velocitygradients.hh:39
Declares all properties used in Dumux.
Type traits for problem classes.