24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
45template<
class TypeTag,
int upwindSchemeOrder>
50 using GridVolumeVariables =
typename GridVariables::GridVolumeVariables;
51 using ElementVolumeVariables =
typename GridVolumeVariables::LocalView;
52 using VolumeVariables =
typename GridVolumeVariables::VolumeVariables;
54 using GridFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache;
55 using FluxVariablesCache =
typename GridFluxVariablesCache::FluxVariablesCache;
57 using GridFaceVariables =
typename GridVariables::GridFaceVariables;
58 using ElementFaceVariables =
typename GridFaceVariables::LocalView;
59 using FaceVariables =
typename GridFaceVariables::FaceVariables;
63 using FVElementGeometry =
typename GridGeometry::LocalView;
64 using GridView =
typename GridGeometry::GridView;
65 using Element =
typename GridView::template Codim<0>::Entity;
68 using Indices =
typename ModelTraits::Indices;
69 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
70 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
76 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
77 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
89 const ElementFaceVariables& elemFaceVars,
90 const ElementVolumeVariables& elemVolVars,
91 const GridFluxVariablesCache& gridFluxVarsCache,
92 const Scalar transportingVelocity)
94 const bool canHigherOrder = canFrontalSecondOrder_(scvf, transportingVelocity, std::integral_constant<bool, useHigherOrder>{});
98 const auto upwindingMomenta = getFrontalUpwindingMomenta_(scvf, elemFaceVars, elemVolVars[scvf.insideScvIdx()].density(),
99 transportingVelocity, std::integral_constant<bool, useHigherOrder>{});
100 return doFrontalMomentumUpwinding_(scvf, upwindingMomenta, transportingVelocity, gridFluxVarsCache);
104 const auto upwindingMomenta = getFrontalUpwindingMomenta_(scvf, elemFaceVars, elemVolVars[scvf.insideScvIdx()].density(),
105 transportingVelocity, std::integral_constant<bool, false>{});
106 return doFrontalMomentumUpwinding_(scvf, upwindingMomenta, transportingVelocity, gridFluxVarsCache);
120 const FVElementGeometry& fvGeometry,
121 const Element& element,
122 const SubControlVolumeFace& scvf,
123 const ElementVolumeVariables& elemVolVars,
124 const FaceVariables& faceVars,
125 const GridFluxVariablesCache& gridFluxVarsCache,
126 const int localSubFaceIdx,
127 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
128 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes)
130 const auto eIdx = scvf.insideScvIdx();
131 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
135 const Scalar transportingVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
138 const bool selfIsUpstream = ( lateralFace.directionSign() ==
sign(transportingVelocity) );
140 const bool canHigherOrder = canLateralSecondOrder_(scvf, selfIsUpstream, localSubFaceIdx, std::integral_constant<bool, useHigherOrder>{});
144 const auto parallelUpwindingMomenta = getLateralUpwindingMomenta_(problem, fvGeometry, element, scvf, elemVolVars, faceVars,
145 transportingVelocity, localSubFaceIdx, currentScvfBoundaryTypes,
146 lateralFaceBoundaryTypes, std::integral_constant<bool, useHigherOrder>{});
147 return doLateralMomentumUpwinding_(fvGeometry, scvf, parallelUpwindingMomenta, transportingVelocity, localSubFaceIdx, gridFluxVarsCache);
151 const auto parallelUpwindingMomenta = getLateralUpwindingMomenta_(problem, fvGeometry, element, scvf, elemVolVars, faceVars,
152 transportingVelocity, localSubFaceIdx, currentScvfBoundaryTypes,
153 lateralFaceBoundaryTypes, std::integral_constant<bool, false>{});
154 return doLateralMomentumUpwinding_(fvGeometry, scvf, parallelUpwindingMomenta, transportingVelocity, localSubFaceIdx, gridFluxVarsCache);
169 static bool canFrontalSecondOrder_(
const SubControlVolumeFace& ownScvf,
170 const Scalar transportingVelocity,
187 static bool canFrontalSecondOrder_(
const SubControlVolumeFace& ownScvf,
188 const Scalar transportingVelocity,
192 const bool selfIsUpstream = ownScvf.directionSign() !=
sign(transportingVelocity);
193 return selfIsUpstream ? ownScvf.hasForwardNeighbor(0) : ownScvf.hasBackwardNeighbor(0);
205 static std::array<Scalar, 2> getFrontalUpwindingMomenta_(
const SubControlVolumeFace& scvf,
206 const ElementFaceVariables& elemFaceVars,
208 const Scalar transportingVelocity,
211 const Scalar momentumSelf = elemFaceVars[scvf].velocitySelf() *
density;
212 const Scalar momentumOpposite = elemFaceVars[scvf].velocityOpposite() *
density;
213 const bool selfIsUpstream = scvf.directionSign() !=
sign(transportingVelocity);
215 return selfIsUpstream ? std::array<Scalar, 2>{momentumOpposite, momentumSelf}
216 : std::array<Scalar, 2>{momentumSelf, momentumOpposite};
228 static std::array<Scalar, 3> getFrontalUpwindingMomenta_(
const SubControlVolumeFace& scvf,
229 const ElementFaceVariables& elemFaceVars,
231 const Scalar transportingVelocity,
234 const bool selfIsUpstream = scvf.directionSign() !=
sign(transportingVelocity);
235 std::array<Scalar, 3> momenta;
239 momenta[0] = elemFaceVars[scvf].velocityOpposite() *
density;
240 momenta[1] = elemFaceVars[scvf].velocitySelf() *
density;
241 momenta[2] = elemFaceVars[scvf].velocityForward(0) *
density;
245 momenta[0] = elemFaceVars[scvf].velocitySelf() *
density;
246 momenta[1] = elemFaceVars[scvf].velocityOpposite() *
density;
247 momenta[2] = elemFaceVars[scvf].velocityBackward(0) *
density;
261 static Scalar doFrontalMomentumUpwinding_(
const SubControlVolumeFace& scvf,
262 const std::array<Scalar, 2>& momenta,
263 const Scalar transportingVelocity,
264 const GridFluxVariablesCache& gridFluxVarsCache)
266 const auto& upwindScheme = gridFluxVarsCache.staggeredUpwindMethods();
267 return upwindScheme.upwind(momenta[0], momenta[1]);
278 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
279 static Scalar doFrontalMomentumUpwinding_(
const SubControlVolumeFace& scvf,
280 const std::array<Scalar, 3>& momenta,
281 const Scalar transportingVelocity,
282 const GridFluxVariablesCache& gridFluxVarsCache)
284 const auto& upwindScheme = gridFluxVarsCache.staggeredUpwindMethods();
285 const bool selfIsUpstream = scvf.directionSign() !=
sign(transportingVelocity);
286 std::array<Scalar,3> distances = getFrontalDistances_(scvf, selfIsUpstream);
287 return upwindScheme.tvd(momenta, distances, selfIsUpstream, upwindScheme.tvdApproach());
296 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
297 static std::array<Scalar, 3> getFrontalDistances_(
const SubControlVolumeFace& ownScvf,
298 const bool selfIsUpstream)
302 std::array<Scalar, 3> distances;
306 distances[0] = ownScvf.selfToOppositeDistance();
307 distances[1] = ownScvf.axisData().inAxisForwardDistances[0];
308 distances[2] = 0.5 * (ownScvf.axisData().selfToOppositeDistance + ownScvf.axisData().inAxisBackwardDistances[0]);
312 distances[0] = ownScvf.selfToOppositeDistance();
313 distances[1] = ownScvf.axisData().inAxisBackwardDistances[0];
314 distances[2] = 0.5 * (ownScvf.axisData().selfToOppositeDistance + ownScvf.axisData().inAxisForwardDistances[0]);
330 static bool canLateralSecondOrder_(
const SubControlVolumeFace& ownScvf,
331 const bool selfIsUpstream,
332 const int localSubFaceIdx,
344 return ownScvf.hasParallelNeighbor(localSubFaceIdx, 0);
352 static bool canLateralSecondOrder_(
const SubControlVolumeFace& ownScvf,
353 const bool selfIsUpstream,
354 const int localSubFaceIdx,
364 static std::array<Scalar, 3> getLateralUpwindingMomenta_(
const Problem& problem,
365 const FVElementGeometry& fvGeometry,
366 const Element& element,
367 const SubControlVolumeFace& ownScvf,
368 const ElementVolumeVariables& elemVolVars,
369 const FaceVariables& faceVars,
370 const Scalar transportingVelocity,
371 const int localSubFaceIdx,
372 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
373 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
376 const SubControlVolumeFace& lateralFace = fvGeometry.scvf(ownScvf.insideScvIdx(), ownScvf.pairData(localSubFaceIdx).localLateralFaceIdx);
379 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
380 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
384 if (!ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
386 const Scalar boundaryMomentum = getParallelVelocityFromBoundary_(problem, element, fvGeometry, ownScvf,
387 faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
388 localSubFaceIdx) * insideVolVars.density();
390 return std::array<Scalar, 3>{boundaryMomentum, boundaryMomentum, boundaryMomentum};
394 const bool selfIsUpstream = lateralFace.directionSign() ==
sign(transportingVelocity);
396 std::array<Scalar, 3> momenta;
400 momenta[1] = faceVars.velocitySelf() * insideVolVars.density();
402 if (ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
403 momenta[0] = faceVars.velocityParallel(localSubFaceIdx, 0) * insideVolVars.density();
405 momenta[0] = getParallelVelocityFromBoundary_(problem, element, fvGeometry, ownScvf,
406 faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
407 localSubFaceIdx) * insideVolVars.density();
410 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
413 if (ownScvf.hasParallelNeighbor(oppositeSubFaceIdx, 0))
414 momenta[2] = faceVars.velocityParallel(oppositeSubFaceIdx, 0) * insideVolVars.density();
416 momenta[2] = getParallelVelocityFromOppositeBoundary_(problem, element, fvGeometry, ownScvf,
417 faceVars, currentScvfBoundaryTypes,
418 oppositeSubFaceIdx) * insideVolVars.density();
422 momenta[0] = faceVars.velocitySelf() * outsideVolVars.density();
423 momenta[1] = faceVars.velocityParallel(localSubFaceIdx, 0) * outsideVolVars.density();
426 if (ownScvf.hasParallelNeighbor(localSubFaceIdx, 1))
427 momenta[2] = faceVars.velocityParallel(localSubFaceIdx, 1) * outsideVolVars.density();
430 const Element& elementParallel = fvGeometry.gridGeometry().element(lateralFace.outsideScvIdx());
431 const SubControlVolumeFace& firstParallelScvf = fvGeometry.scvf(lateralFace.outsideScvIdx(), ownScvf.localFaceIdx());
433 momenta[2] = getParallelVelocityFromOppositeBoundary_(problem, elementParallel, fvGeometry, firstParallelScvf,
434 faceVars, problem.boundaryTypes(elementParallel, firstParallelScvf),
435 localSubFaceIdx) * outsideVolVars.density();
447 static std::array<Scalar, 2> getLateralUpwindingMomenta_(
const Problem& problem,
448 const FVElementGeometry& fvGeometry,
449 const Element& element,
450 const SubControlVolumeFace& ownScvf,
451 const ElementVolumeVariables& elemVolVars,
452 const FaceVariables& faceVars,
453 const Scalar transportingVelocity,
454 const int localSubFaceIdx,
455 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
456 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
460 const SubControlVolumeFace& lateralFace = fvGeometry.scvf(ownScvf.insideScvIdx(), ownScvf.pairData(localSubFaceIdx).localLateralFaceIdx);
463 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
464 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
466 const Scalar momentumParallel = ownScvf.hasParallelNeighbor(localSubFaceIdx, 0)
467 ? faceVars.velocityParallel(localSubFaceIdx, 0) * outsideVolVars.density()
468 : (getParallelVelocityFromBoundary_(problem, element, fvGeometry, ownScvf,
469 faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
471 * insideVolVars.density());
475 if (!ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
476 return std::array<Scalar, 2>{momentumParallel, momentumParallel};
478 const bool selfIsUpstream = lateralFace.directionSign() ==
sign(transportingVelocity);
479 const Scalar momentumSelf = faceVars.velocitySelf() * insideVolVars.density();
481 return selfIsUpstream ? std::array<Scalar, 2>{momentumParallel, momentumSelf}
482 : std::array<Scalar, 2>{momentumSelf, momentumParallel};
490 static Scalar doLateralMomentumUpwinding_(
const FVElementGeometry& fvGeometry,
491 const SubControlVolumeFace& scvf,
492 const std::array<Scalar, 2>& momenta,
493 const Scalar transportingVelocity,
494 const int localSubFaceIdx,
495 const GridFluxVariablesCache& gridFluxVarsCache)
497 return doFrontalMomentumUpwinding_(scvf, momenta, transportingVelocity, gridFluxVarsCache);
509 static Scalar doLateralMomentumUpwinding_(
const FVElementGeometry& fvGeometry,
510 const SubControlVolumeFace& scvf,
511 const std::array<Scalar, 3>& momenta,
512 const Scalar transportingVelocity,
513 const int localSubFaceIdx,
514 const GridFluxVariablesCache& gridFluxVarsCache)
516 const auto eIdx = scvf.insideScvIdx();
517 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
519 const bool selfIsUpstream = ( lateralFace.directionSign() ==
sign(transportingVelocity) );
520 const auto& upwindScheme = gridFluxVarsCache.staggeredUpwindMethods();
521 std::array<Scalar,3> distances = getLateralDistances_(scvf, localSubFaceIdx, selfIsUpstream);
522 return upwindScheme.tvd(momenta, distances, selfIsUpstream, upwindScheme.tvdApproach());
532 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
533 static std::array<Scalar, 3> getLateralDistances_(
const SubControlVolumeFace& ownScvf,
534 const int localSubFaceIdx,
535 const bool selfIsUpstream)
538 std::array<Scalar, 3> distances;
543 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
545 distances[0] = ownScvf.parallelDofsDistance(localSubFaceIdx, 0);
546 distances[1] = ownScvf.parallelDofsDistance(oppositeSubFaceIdx, 0);
547 if (ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
548 distances[2] = ownScvf.pairData(localSubFaceIdx).parallelCellWidths[0];
550 distances[2] = ownScvf.area() / 2.0;
554 distances[0] = ownScvf.parallelDofsDistance(localSubFaceIdx, 0);
555 distances[1] = ownScvf.parallelDofsDistance(localSubFaceIdx, 1);
556 distances[2] = ownScvf.pairData(localSubFaceIdx).parallelCellWidths[0];
574 static Scalar getParallelVelocityFromBoundary_(
const Problem& problem,
575 const Element& element,
576 const FVElementGeometry& fvGeometry,
577 const SubControlVolumeFace& scvf,
578 const FaceVariables& faceVars,
579 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
580 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
581 const int localSubFaceIdx)
584 const bool useZeroGradient = lateralFaceBoundaryTypes && (lateralFaceBoundaryTypes->isSymmetry()
585 || lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
586 const bool lateralFaceHasBJS = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex()));
587 const bool lateralFaceHasDirichletVelocity = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex()));
588 const Scalar velocitySelf = faceVars.velocitySelf();
595 if (lateralFaceHasBJS)
597 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
599 else if(lateralFaceHasDirichletVelocity)
610 const auto ghostFace = makeParallelGhostFace_(scvf, localSubFaceIdx);
611 return problem.dirichlet(element, ghostFace)[Indices::velocity(scvf.directionIndex())];
616 DUNE_THROW(Dune::InvalidStateException,
"Something went wrong with the boundary conditions for the momentum equations at global position " << fvGeometry.scvf(scvf.insideScvIdx(), scvf.pairData(localSubFaceIdx).localLateralFaceIdx).center());
629 static Scalar getParallelVelocityFromOppositeBoundary_(
const Problem& problem,
630 const Element& element,
631 const FVElementGeometry& fvGeometry,
632 const SubControlVolumeFace& scvf,
633 const FaceVariables& faceVars,
634 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
635 const int localOppositeSubFaceIdx)
638 const SubControlVolumeFace& lateralOppositeScvf = fvGeometry.scvf(scvf.insideScvIdx(), scvf.pairData(localOppositeSubFaceIdx).localLateralFaceIdx);
639 GlobalPosition center = scvf.pairData(localOppositeSubFaceIdx).lateralStaggeredFaceCenter + lateralOppositeScvf.center();
654 const auto lateralOppositeFaceBoundaryTypes = problem.boundaryTypes(element, lateralOppositeScvf.makeBoundaryFace(center));
656 return getParallelVelocityFromBoundary_(problem, element, fvGeometry, scvf,
657 faceVars, currentScvfBoundaryTypes, lateralOppositeFaceBoundaryTypes,
658 localOppositeSubFaceIdx);
662 static SubControlVolumeFace makeParallelGhostFace_(
const SubControlVolumeFace& ownScvf,
const int localSubFaceIdx)
664 return ownScvf.makeBoundaryFace(ownScvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter);
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.
This file contains different higher order methods for approximating the velocity.
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
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 density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
The upwinding variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: staggeredupwindfluxvariables.hh:47
static FacePrimaryVariables computeUpwindedLateralMomentum(const Problem &problem, const FVElementGeometry &fvGeometry, const Element &element, const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const FaceVariables &faceVars, const GridFluxVariablesCache &gridFluxVarsCache, const int localSubFaceIdx, const std::optional< BoundaryTypes > ¤tScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes)
Returns the momentum in the lateral directon.
Definition: staggeredupwindfluxvariables.hh:119
static FacePrimaryVariables computeUpwindedFrontalMomentum(const SubControlVolumeFace &scvf, const ElementFaceVariables &elemFaceVars, const ElementVolumeVariables &elemVolVars, const GridFluxVariablesCache &gridFluxVarsCache, const Scalar transportingVelocity)
Returns the momentum in the frontal directon.
Definition: staggeredupwindfluxvariables.hh:88
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: velocitygradients.hh:39
static Scalar beaversJosephVelocityAtLateralScvf(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > ¤tScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes, const std::size_t localSubFaceIdx)
Returns the Beavers-Jospeh slip velocity for a lateral scvf which lies on the boundary.
Definition: velocitygradients.hh:301
Declares all properties used in Dumux.