39    using GridVolumeVariables = 
typename GridVariables::GridVolumeVariables;
 
   40    using ElementVolumeVariables = 
typename GridVolumeVariables::LocalView;
 
   41    using VolumeVariables = 
typename GridVolumeVariables::VolumeVariables;
 
   43    using GridFluxVariablesCache = 
typename GridVariables::GridFluxVariablesCache;
 
   44    using UpwindScheme = 
typename GridFluxVariablesCache::UpwindScheme;
 
   45    using FluxVariablesCache = 
typename GridFluxVariablesCache::FluxVariablesCache;
 
   47    using GridFaceVariables = 
typename GridVariables::GridFaceVariables;
 
   48    using ElementFaceVariables = 
typename GridFaceVariables::LocalView;
 
   49    using FaceVariables = 
typename GridFaceVariables::FaceVariables;
 
   53    using FVElementGeometry = 
typename GridGeometry::LocalView;
 
   54    using GridView = 
typename GridGeometry::GridView;
 
   55    using Element = 
typename GridView::template Codim<0>::Entity;
 
   58    using Indices = 
typename ModelTraits::Indices;
 
   59    using SubControlVolume = 
typename FVElementGeometry::SubControlVolume;
 
   60    using SubControlVolumeFace = 
typename FVElementGeometry::SubControlVolumeFace;
 
   66    using GlobalPosition = 
typename Element::Geometry::GlobalCoordinate;
 
   67    static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
 
   69    static_assert(upwindSchemeOrder <= 2, 
"Not implemented: Order higher than 2!");
 
   70    static_assert(upwindSchemeOrder <= 1 || GridFluxVariablesCache::cachingEnabled,
 
   71        "Higher order upwind method requires caching enabled for the GridFluxVariablesCache!");
 
   72    static_assert(upwindSchemeOrder <= 1 || GridGeometry::cachingEnabled,
 
   73        "Higher order upwind method requires caching enabled for the GridGeometry!");
 
   74    static_assert(upwindSchemeOrder <= 1 || GridFaceVariables::cachingEnabled,
 
   75        "Higher order upwind method requires caching enabled for the GridFaceVariables!");
 
   76    static_assert(upwindSchemeOrder <= 1 || GridVolumeVariables::cachingEnabled,
 
   77        "Higher order upwind method requires caching enabled for the GridGeometry!");
 
   81                          const FVElementGeometry& fvGeometry,
 
   82                          const SubControlVolumeFace& scvf,
 
   83                          const ElementFaceVariables& elemFaceVars,
 
   84                          const ElementVolumeVariables& elemVolVars,
 
   85                          const UpwindScheme& upwindScheme)
 
   87    , fvGeometry_(fvGeometry)
 
   89    , elemFaceVars_(elemFaceVars)
 
   90    , faceVars_(elemFaceVars[scvf])
 
   91    , elemVolVars_(elemVolVars)
 
   92    , upwindScheme_(upwindScheme)
 
 
  105        const auto density = elemVolVars_[scvf_.insideScvIdx()].density();
 
  108        if constexpr (useHigherOrder)
 
  111            if (canDoFrontalSecondOrder_(selfIsUpstream))
 
  113                const auto distances = getFrontalDistances_(selfIsUpstream);
 
  114                const auto upwindMomenta = getFrontalSecondOrderUpwindMomenta_(density, selfIsUpstream);
 
  115                return upwindScheme_.tvd(upwindMomenta, distances, selfIsUpstream, upwindScheme_.tvdApproach());
 
  120        const auto upwindMomenta = getFrontalFirstOrderUpwindMomenta_(density, selfIsUpstream);
 
  121        return upwindScheme_.upwind(upwindMomenta[0], upwindMomenta[1]);
 
 
  134                                                      const SubControlVolumeFace& lateralFace,
 
  135                                                      const int localSubFaceIdx,
 
  136                                                      const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
 
  137                                                      const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes)
 const 
  142        const Scalar insideDensity = elemVolVars_[lateralFace.insideScvIdx()].density();
 
  143        const Scalar outsideDensity = elemVolVars_[lateralFace.outsideScvIdx()].density();
 
  146        if constexpr (useHigherOrder)
 
  149            if (canDoLateralSecondOrder_(selfIsUpstream, localSubFaceIdx))
 
  151                const auto upwindMomenta = getLateralSecondOrderUpwindMomenta_(insideDensity, outsideDensity, selfIsUpstream, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes);
 
  152                const auto distances = getLateralDistances_(localSubFaceIdx, selfIsUpstream);
 
  153                return upwindScheme_.tvd(upwindMomenta, distances, selfIsUpstream, upwindScheme_.tvdApproach());
 
  158        const auto upwindMomenta = getLateralFirstOrderUpwindMomenta_(insideDensity, outsideDensity, selfIsUpstream, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes);
 
  159        return upwindScheme_.upwind(upwindMomenta[0], upwindMomenta[1]);
 
 
  171    bool canDoFrontalSecondOrder_(
bool selfIsUpstream)
 const 
  174        return selfIsUpstream ? scvf_.hasForwardNeighbor(0) : scvf_.hasBackwardNeighbor(0);
 
  180    std::array<Scalar, 3> getFrontalSecondOrderUpwindMomenta_(
const Scalar density, 
bool selfIsUpstream)
 const 
  182        static_assert(useHigherOrder, 
"Should only be reached if higher order methods are enabled");
 
  183        const Scalar momentumSelf = faceVars_.velocitySelf() * density;
 
  184        const Scalar momentumOpposite = faceVars_.velocityOpposite() * density;
 
  186            return {momentumOpposite, momentumSelf, faceVars_.velocityForward(0)*density};
 
  188            return {momentumSelf, momentumOpposite, faceVars_.velocityBackward(0)*
density};
 
  194    std::array<Scalar, 2> getFrontalFirstOrderUpwindMomenta_(
const Scalar density, 
bool selfIsUpstream)
 const 
  196        const Scalar momentumSelf = faceVars_.velocitySelf() * 
density;
 
  197        const Scalar momentumOpposite = faceVars_.velocityOpposite() * 
density;
 
  199            return {momentumOpposite, momentumSelf};
 
  201            return {momentumSelf, momentumOpposite};
 
  209    std::array<Scalar, 3> getFrontalDistances_(
const bool selfIsUpstream)
 const 
  211        static_assert(useHigherOrder, 
"Should only be reached if higher order methods are enabled");
 
  214            std::array<Scalar, 3> distances;
 
  215            distances[0] = scvf_.selfToOppositeDistance();
 
  216            distances[1] = scvf_.axisData().inAxisForwardDistances[0];
 
  217            distances[2] = 0.5 * (scvf_.axisData().selfToOppositeDistance + scvf_.axisData().inAxisBackwardDistances[0]);
 
  222            std::array<Scalar, 3> distances;
 
  223            distances[0] = scvf_.selfToOppositeDistance();
 
  224            distances[1] = scvf_.axisData().inAxisBackwardDistances[0];
 
  225            distances[2] = 0.5 * (scvf_.axisData().selfToOppositeDistance + scvf_.axisData().inAxisForwardDistances[0]);
 
  236    bool canDoLateralSecondOrder_(
const bool selfIsUpstream, 
const int localSubFaceIdx)
 const 
  241        return selfIsUpstream || scvf_.hasParallelNeighbor(localSubFaceIdx, 0);
 
  284    std::array<Scalar, 3> getLateralSecondOrderUpwindMomenta_(
const Scalar insideDensity,
 
  285                                                              const Scalar outsideDensity,
 
  286                                                              const bool selfIsUpstream,
 
  287                                                              const int localSubFaceIdx,
 
  288                                                              const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
 
  289                                                              const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes)
 const 
  291        static_assert(useHigherOrder, 
"Should only be reached if higher order methods are enabled");
 
  295        if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0) || scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
 
  297            if ((scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx)) && dirichletParallelNeighbor_(localSubFaceIdx))
 
  299                const Scalar boundaryVelocity = getParallelVelocityFromCorner_(localSubFaceIdx);
 
  300                const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
 
  302                return {boundaryMomentum, boundaryMomentum, boundaryMomentum};
 
  304            else if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0))
 
  306                const Scalar boundaryVelocity = getParallelVelocityFromBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
 
  307                const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
 
  309                return {boundaryMomentum, boundaryMomentum, boundaryMomentum};
 
  315            std::array<Scalar, 3> momenta;
 
  316            momenta[1] = faceVars_.velocitySelf()*insideDensity;
 
  317            momenta[0] = faceVars_.velocityParallel(localSubFaceIdx, 0)*insideDensity;
 
  320            const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
 
  323            if (scvf_.hasParallelNeighbor(oppositeSubFaceIdx, 0))
 
  324                momenta[2] = faceVars_.velocityParallel(oppositeSubFaceIdx, 0)*insideDensity;
 
  326                momenta[2] = getParallelVelocityFromOppositeBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, oppositeSubFaceIdx)*insideDensity;
 
  331            std::array<Scalar, 3> momenta;
 
  332            momenta[0] = faceVars_.velocitySelf()*outsideDensity;
 
  333            momenta[1] = faceVars_.velocityParallel(localSubFaceIdx, 0)*outsideDensity;
 
  336            if (scvf_.hasParallelNeighbor(localSubFaceIdx, 1))
 
  337                momenta[2] = faceVars_.velocityParallel(localSubFaceIdx, 1)*outsideDensity;
 
  340                const auto& lateralFace = fvGeometry_.scvf(scvf_.insideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
 
  341                const auto& elementParallel = fvGeometry_.gridGeometry().element(lateralFace.outsideScvIdx());
 
  342                const auto& firstParallelScvf = fvGeometry_.scvf(lateralFace.outsideScvIdx(), scvf_.localFaceIdx());
 
  343                const auto& problem = elemVolVars_.gridVolVars().problem();
 
  344                const auto& boundaryTypes = problem.boundaryTypes(elementParallel, firstParallelScvf);
 
  345                momenta[2] = getParallelVelocityFromOppositeBoundary_(elementParallel, firstParallelScvf, elemFaceVars_[firstParallelScvf], boundaryTypes, localSubFaceIdx)*outsideDensity;
 
  354    std::array<Scalar, 2> getLateralFirstOrderUpwindMomenta_(
const Scalar insideDensity,
 
  355                                                             const Scalar outsideDensity,
 
  356                                                             const bool selfIsUpstream,
 
  357                                                             const int localSubFaceIdx,
 
  358                                                             const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
 
  359                                                             const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes)
 const 
  363        if ((scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx)) && dirichletParallelNeighbor_(localSubFaceIdx))
 
  365            const Scalar boundaryVelocity =  getParallelVelocityFromCorner_(localSubFaceIdx);
 
  366            const Scalar boundaryMomentum = boundaryVelocity*outsideDensity;
 
  367            return {boundaryMomentum, boundaryMomentum};
 
  369        else if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0))
 
  371            const Scalar boundaryVelocity = getParallelVelocityFromBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
 
  372            const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
 
  373            return {boundaryMomentum, boundaryMomentum};
 
  376        const Scalar momentumParallel = faceVars_.velocityParallel(localSubFaceIdx, 0)*outsideDensity;
 
  377        const Scalar momentumSelf = faceVars_.velocitySelf()*insideDensity;
 
  379            return {momentumParallel, momentumSelf};
 
  381            return {momentumSelf, momentumParallel};
 
  388    std::array<Scalar, 3> getLateralDistances_(
const int localSubFaceIdx, 
const bool selfIsUpstream)
 const 
  390        static_assert(useHigherOrder, 
"Should only be reached if higher order methods are enabled");
 
  394            const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
 
  396            std::array<Scalar, 3> distances;
 
  397            distances[0] = scvf_.parallelDofsDistance(localSubFaceIdx, 0);
 
  398            distances[1] = scvf_.parallelDofsDistance(oppositeSubFaceIdx, 0);
 
  399            if (scvf_.hasParallelNeighbor(localSubFaceIdx, 0))
 
  400                distances[2] = scvf_.pairData(localSubFaceIdx).parallelCellWidths[0];
 
  402                distances[2] = scvf_.area() / 2.0;
 
  407            std::array<Scalar, 3> distances;
 
  408            distances[0] = scvf_.parallelDofsDistance(localSubFaceIdx, 0);
 
  409            distances[1] = scvf_.parallelDofsDistance(localSubFaceIdx, 1);
 
  410            distances[2] = scvf_.pairData(localSubFaceIdx).parallelCellWidths[0];
 
  419    Scalar getParallelVelocityFromBoundary_(
const Element& element,
 
  420                                            const SubControlVolumeFace& scvf,
 
  421                                            const FaceVariables& faceVars,
 
  422                                            const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
 
  423                                            const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
 
  424                                            const int localSubFaceIdx)
 const 
  428        const bool useZeroGradient = lateralFaceBoundaryTypes && (lateralFaceBoundaryTypes->isSymmetry() || lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
 
  430            return faceVars.velocitySelf();
 
  432        const bool lateralFaceHasBJS = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex()));
 
  433        if (lateralFaceHasBJS)
 
  435                                                                         currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
 
  437        const bool lateralFaceHasDirichletVelocity = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex()));
 
  438        if (lateralFaceHasDirichletVelocity)
 
  449            const auto& lateralFace = fvGeometry_.scvf(scvf.insideScvIdx(), scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
 
  450            const auto ghostFace = 
makeStaggeredBoundaryFace(lateralFace, scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter);
 
  451            const auto& problem = elemVolVars_.gridVolVars().problem();
 
  452            return problem.dirichlet(element, ghostFace)[Indices::velocity(scvf.directionIndex())];
 
  456        DUNE_THROW(Dune::InvalidStateException, 
"Something went wrong with the boundary conditions for the momentum equations at global position " 
  457                    << fvGeometry_.scvf(scvf.insideScvIdx(), scvf.pairData(localSubFaceIdx).localLateralFaceIdx).center());
 
  463    Scalar getParallelVelocityFromOppositeBoundary_(
const Element& element,
 
  464                                                    const SubControlVolumeFace& scvf,
 
  465                                                    const FaceVariables& faceVars,
 
  466                                                    const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
 
  467                                                    const int localOppositeSubFaceIdx)
 const 
  470        const auto& lateralOppositeScvf = fvGeometry_.scvf(scvf.insideScvIdx(), scvf.pairData(localOppositeSubFaceIdx).localLateralFaceIdx);
 
  471        GlobalPosition 
center = scvf.pairData(localOppositeSubFaceIdx).lateralStaggeredFaceCenter + lateralOppositeScvf.center();
 
  486        const auto& problem = elemVolVars_.gridVolVars().problem();
 
  488        const auto lateralOppositeFaceBoundaryTypes = problem.boundaryTypes(element, lateralOppositeBoundaryFace);
 
  489        return getParallelVelocityFromBoundary_(element, scvf, faceVars, currentScvfBoundaryTypes, lateralOppositeFaceBoundaryTypes, localOppositeSubFaceIdx);
 
  511    const SubControlVolumeFace& boundaryScvf_(
const int localSubFaceIdx)
 const 
  513        if (scvf_.hasHalfParallelNeighbor(localSubFaceIdx))
 
  515            return fvGeometry_.scvf(scvf_.outsideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
 
  517        else if (scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
 
  533            const SubControlVolumeFace& lateralFace = fvGeometry_.scvf(scvf_.insideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
 
  534            const SubControlVolumeFace& parallelFace = fvGeometry_.scvf(lateralFace.outsideScvIdx(), scvf_.localFaceIdx());
 
  536            const auto& localLateralIdx = scvf_.pairData(localSubFaceIdx).localLateralFaceIdx;
 
  537            const auto& localLateralOppositeIdx =  (localLateralIdx % 2) ? (localLateralIdx - 1) : (localLateralIdx + 1);
 
  539            return fvGeometry_.scvf(parallelFace.outsideScvIdx(), localLateralOppositeIdx);
 
  543            DUNE_THROW(Dune::InvalidStateException, 
"The function boundaryScvf_ should only be called when hasHalfParallelNeighbor or hasCornerParallelNeighbor.");
 
  566    Element boundaryElement_(
const int localSubFaceIdx)
 const 
  568        if (scvf_.hasHalfParallelNeighbor(localSubFaceIdx))
 
  570            return fvGeometry_.gridGeometry().element(scvf_.outsideScvIdx());
 
  572        else if (scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
 
  574            const SubControlVolumeFace& lateralFace = fvGeometry_.scvf(scvf_.insideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
 
  575            const SubControlVolumeFace& parallelFace = fvGeometry_.scvf(lateralFace.outsideScvIdx(), scvf_.localFaceIdx());
 
  577            return fvGeometry_.gridGeometry().element(parallelFace.outsideScvIdx());
 
  581            DUNE_THROW(Dune::InvalidStateException, 
"When entering boundaryElement_ scvf_ should have either hasHalfParallelNeighbor or hasCornerParallelNeighbor true. Not the case here.");
 
  604    bool dirichletParallelNeighbor_(
const int localSubFaceIdx)
 const 
  606        const auto& problem = elemVolVars_.gridVolVars().problem();
 
  607        const Element& boundaryElement = boundaryElement_(localSubFaceIdx);
 
  608        const SubControlVolumeFace& boundaryScvf = boundaryScvf_(localSubFaceIdx);
 
  610        return problem.boundaryTypes(boundaryElement, boundaryScvf).isDirichlet(Indices::velocity(scvf_.directionIndex()));
 
  631    Scalar getParallelVelocityFromCorner_(
const int localSubFaceIdx)
 const 
  633        const auto& problem = elemVolVars_.gridVolVars().problem();
 
  634        const Element& boundaryElement = boundaryElement_(localSubFaceIdx);
 
  635        const SubControlVolumeFace& boundaryScvf = boundaryScvf_(localSubFaceIdx);
 
  636        const auto ghostFace = 
makeStaggeredBoundaryFace(boundaryScvf, scvf_.pairData(localSubFaceIdx).lateralStaggeredFaceCenter);
 
  638        return problem.dirichlet(boundaryElement, ghostFace)[Indices::velocity(scvf_.directionIndex())];
 
  641    const Element& element_;
 
  642    const FVElementGeometry& fvGeometry_;
 
  643    const SubControlVolumeFace& scvf_;
 
  644    const ElementFaceVariables& elemFaceVars_;
 
  645    const FaceVariables& faceVars_;
 
  646    const ElementVolumeVariables& elemVolVars_;
 
  647    const UpwindScheme& upwindScheme_;