version 3.8
staggeredupwindhelper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
13#define DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
14
15#include <array>
16#include <optional>
17
18#include <dumux/common/math.hh>
23
26#include "velocitygradients.hh"
27
28namespace Dumux {
29
34template<class TypeTag, int upwindSchemeOrder>
36{
38
39 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
40 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
41 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
42
43 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
44 using UpwindScheme = typename GridFluxVariablesCache::UpwindScheme;
45 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
46
47 using GridFaceVariables = typename GridVariables::GridFaceVariables;
48 using ElementFaceVariables = typename GridFaceVariables::LocalView;
49 using FaceVariables = typename GridFaceVariables::FaceVariables;
50
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;
61 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
63 using BoundaryTypes = typename ProblemTraits<Problem>::BoundaryTypes;
65
66 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
67 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
68
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!");
78
79public:
80 StaggeredUpwindHelper(const Element& element,
81 const FVElementGeometry& fvGeometry,
82 const SubControlVolumeFace& scvf,
83 const ElementFaceVariables& elemFaceVars,
84 const ElementVolumeVariables& elemVolVars,
85 const UpwindScheme& upwindScheme)
86 : element_(element)
87 , fvGeometry_(fvGeometry)
88 , scvf_(scvf)
89 , elemFaceVars_(elemFaceVars)
90 , faceVars_(elemFaceVars[scvf])
91 , elemVolVars_(elemVolVars)
92 , upwindScheme_(upwindScheme)
93 {}
94
103 FacePrimaryVariables computeUpwindFrontalMomentum(const bool selfIsUpstream) const
104 {
105 const auto density = elemVolVars_[scvf_.insideScvIdx()].density();
106
107 // for higher order schemes do higher order upwind reconstruction
108 if constexpr (useHigherOrder)
109 {
110 // only second order is implemented so far
111 if (canDoFrontalSecondOrder_(selfIsUpstream))
112 {
113 const auto distances = getFrontalDistances_(selfIsUpstream);
114 const auto upwindMomenta = getFrontalSecondOrderUpwindMomenta_(density, selfIsUpstream);
115 return upwindScheme_.tvd(upwindMomenta, distances, selfIsUpstream, upwindScheme_.tvdApproach());
116 }
117 }
118
119 // otherwise apply first order upwind scheme
120 const auto upwindMomenta = getFrontalFirstOrderUpwindMomenta_(density, selfIsUpstream);
121 return upwindScheme_.upwind(upwindMomenta[0], upwindMomenta[1]);
122 }
123
133 FacePrimaryVariables computeUpwindLateralMomentum(const bool selfIsUpstream,
134 const SubControlVolumeFace& lateralFace,
135 const int localSubFaceIdx,
136 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
137 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes) const
138 {
139 // Check whether the own or the neighboring element is upstream.
140 // Get the transporting velocity, located at the scvf perpendicular to the current scvf where the dof
141 // of interest is located.
142 const Scalar insideDensity = elemVolVars_[lateralFace.insideScvIdx()].density();
143 const Scalar outsideDensity = elemVolVars_[lateralFace.outsideScvIdx()].density();
144
145 // for higher order schemes do higher order upwind reconstruction
146 if constexpr (useHigherOrder)
147 {
148 // only second order is implemented so far
149 if (canDoLateralSecondOrder_(selfIsUpstream, localSubFaceIdx))
150 {
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());
154 }
155 }
156
157 // otherwise apply first order upwind scheme
158 const auto upwindMomenta = getLateralFirstOrderUpwindMomenta_(insideDensity, outsideDensity, selfIsUpstream, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes);
159 return upwindScheme_.upwind(upwindMomenta[0], upwindMomenta[1]);
160 }
161
162private:
171 bool canDoFrontalSecondOrder_(bool selfIsUpstream) const
172 {
173 // Depending on selfIsUpstream I have to check if I have a forward or a backward neighbor to retrieve
174 return selfIsUpstream ? scvf_.hasForwardNeighbor(0) : scvf_.hasBackwardNeighbor(0);
175 }
176
180 std::array<Scalar, 3> getFrontalSecondOrderUpwindMomenta_(const Scalar density, bool selfIsUpstream) const
181 {
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;
185 if (selfIsUpstream)
186 return {momentumOpposite, momentumSelf, faceVars_.velocityForward(0)*density};
187 else
188 return {momentumSelf, momentumOpposite, faceVars_.velocityBackward(0)*density};
189 }
190
194 std::array<Scalar, 2> getFrontalFirstOrderUpwindMomenta_(const Scalar density, bool selfIsUpstream) const
195 {
196 const Scalar momentumSelf = faceVars_.velocitySelf() * density;
197 const Scalar momentumOpposite = faceVars_.velocityOpposite() * density;
198 if (selfIsUpstream)
199 return {momentumOpposite, momentumSelf};
200 else
201 return {momentumSelf, momentumOpposite};
202 }
203
209 std::array<Scalar, 3> getFrontalDistances_(const bool selfIsUpstream) const
210 {
211 static_assert(useHigherOrder, "Should only be reached if higher order methods are enabled");
212 if (selfIsUpstream)
213 {
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]);
218 return distances;
219 }
220 else
221 {
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]);
226 return distances;
227 }
228 }
229
236 bool canDoLateralSecondOrder_(const bool selfIsUpstream, const int localSubFaceIdx) const
237 {
238 // If the self velocity is upstream, the downstream velocity can be assigned or retrieved
239 // from the boundary, even if there is no parallel neighbor.
240 // If the self velocity is downstream and there is no parallel neighbor I cannot use a second order approximation.
241 return selfIsUpstream || scvf_.hasParallelNeighbor(localSubFaceIdx, 0);
242 }
243
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
290 {
291 static_assert(useHigherOrder, "Should only be reached if higher order methods are enabled");
292
293 // If the lateral face lies on a boundary, we assume that the parallel velocity on the boundary is actually known,
294 // thus we always use this value for the computation of the transported momentum.
295 if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0) || scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
296 {
297 if ((scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx)) && dirichletParallelNeighbor_(localSubFaceIdx))
298 {
299 const Scalar boundaryVelocity = getParallelVelocityFromCorner_(localSubFaceIdx);
300 const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
301
302 return {boundaryMomentum, boundaryMomentum, boundaryMomentum};
303 }
304 else if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0))
305 {
306 const Scalar boundaryVelocity = getParallelVelocityFromBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
307 const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
308
309 return {boundaryMomentum, boundaryMomentum, boundaryMomentum};
310 }
311 }
312
313 if (selfIsUpstream)
314 {
315 std::array<Scalar, 3> momenta;
316 momenta[1] = faceVars_.velocitySelf()*insideDensity;
317 momenta[0] = faceVars_.velocityParallel(localSubFaceIdx, 0)*insideDensity;
318
319 // The local index of the faces that is opposite to localSubFaceIdx
320 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
321
322 // The "upstream-upstream" velocity is retrieved from the other parallel neighbor or from the boundary.
323 if (scvf_.hasParallelNeighbor(oppositeSubFaceIdx, 0))
324 momenta[2] = faceVars_.velocityParallel(oppositeSubFaceIdx, 0)*insideDensity;
325 else
326 momenta[2] = getParallelVelocityFromOppositeBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, oppositeSubFaceIdx)*insideDensity;
327 return momenta;
328 }
329 else
330 {
331 std::array<Scalar, 3> momenta;
332 momenta[0] = faceVars_.velocitySelf()*outsideDensity;
333 momenta[1] = faceVars_.velocityParallel(localSubFaceIdx, 0)*outsideDensity;
334
335 // If there is another parallel neighbor I can assign the "upstream-upstream" velocity, otherwise I retrieve it from the boundary.
336 if (scvf_.hasParallelNeighbor(localSubFaceIdx, 1))
337 momenta[2] = faceVars_.velocityParallel(localSubFaceIdx, 1)*outsideDensity;
338 else
339 {
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;
346 }
347 return momenta;
348 }
349 }
350
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
360 {
361 // If the lateral face lies on a boundary, we assume that the parallel velocity on the boundary is actually known,
362 // thus we always use this value for the computation of the transported momentum.
363 if ((scvf_.hasHalfParallelNeighbor(localSubFaceIdx) || scvf_.hasCornerParallelNeighbor(localSubFaceIdx)) && dirichletParallelNeighbor_(localSubFaceIdx))
364 {
365 const Scalar boundaryVelocity = getParallelVelocityFromCorner_(localSubFaceIdx);
366 const Scalar boundaryMomentum = boundaryVelocity*outsideDensity;
367 return {boundaryMomentum, boundaryMomentum};
368 }
369 else if (!scvf_.hasParallelNeighbor(localSubFaceIdx, 0))
370 {
371 const Scalar boundaryVelocity = getParallelVelocityFromBoundary_(element_, scvf_, faceVars_, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
372 const Scalar boundaryMomentum = boundaryVelocity*insideDensity;
373 return {boundaryMomentum, boundaryMomentum};
374 }
375
376 const Scalar momentumParallel = faceVars_.velocityParallel(localSubFaceIdx, 0)*outsideDensity;
377 const Scalar momentumSelf = faceVars_.velocitySelf()*insideDensity;
378 if (selfIsUpstream)
379 return {momentumParallel, momentumSelf};
380 else
381 return {momentumSelf, momentumParallel};
382 }
383
388 std::array<Scalar, 3> getLateralDistances_(const int localSubFaceIdx, const bool selfIsUpstream) const
389 {
390 static_assert(useHigherOrder, "Should only be reached if higher order methods are enabled");
391 if (selfIsUpstream)
392 {
393 // The local index of the faces that is opposite to localSubFaceIdx
394 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
395
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];
401 else
402 distances[2] = scvf_.area() / 2.0;
403 return distances;
404 }
405 else
406 {
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];
411 return distances;
412 }
413 }
414
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
425 {
426 // If there is a Dirichlet condition for the pressure we assume zero gradient for the velocity,
427 // so the velocity at the boundary equal to that on the scvf.
428 const bool useZeroGradient = lateralFaceBoundaryTypes && (lateralFaceBoundaryTypes->isSymmetry() || lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
429 if (useZeroGradient)
430 return faceVars.velocitySelf();
431
432 const bool lateralFaceHasBJS = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex()));
433 if (lateralFaceHasBJS)
434 return VelocityGradients::beaversJosephVelocityAtLateralScvf(elemVolVars_.gridVolVars().problem(), element, fvGeometry_, scvf, faceVars,
435 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
436
437 const bool lateralFaceHasDirichletVelocity = lateralFaceBoundaryTypes && lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex()));
438 if (lateralFaceHasDirichletVelocity)
439 {
440 // ________________
441 // ---------------o || frontal face of staggered half-control-volume
442 // | || # current scvf # ghostFace of interest, lies on boundary
443 // | || # x dof position
444 // | || x~~~~> vel.Self -- element boundaries
445 // | || # __ domain boundary
446 // | || # o position at which the boundary conditions will be evaluated
447 // ---------------#
448
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())];
453 }
454
455 // Neumann conditions are not well implemented
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());
458 }
459
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
468 {
469 // A ghost subface at the boundary is created, featuring the location of the sub face's center
470 const auto& lateralOppositeScvf = fvGeometry_.scvf(scvf.insideScvIdx(), scvf.pairData(localOppositeSubFaceIdx).localLateralFaceIdx);
471 GlobalPosition center = scvf.pairData(localOppositeSubFaceIdx).lateralStaggeredFaceCenter + lateralOppositeScvf.center();
472 center *= 0.5;
473
474 // lateral face # lateralFace currently being assembled
475 // --------######## || frontal face of staggered half-control-volume
476 // | || | current scvf % lateralOppositeBoundaryFace of interest, lies on boundary
477 // | || | x dof position
478 // | || x~~~~> vel.Self -- element boundaries
479 // | || | __ domain boundary
480 // | || | o position at which the boundary conditions will be evaluated
481 // --------%%%c%%%o
482 // ________________ c center of opposite boundary face
483
484
485 // Get the boundary types of the lateral opposite boundary face
486 const auto& problem = elemVolVars_.gridVolVars().problem();
487 const auto lateralOppositeBoundaryFace = makeStaggeredBoundaryFace(lateralOppositeScvf, center);
488 const auto lateralOppositeFaceBoundaryTypes = problem.boundaryTypes(element, lateralOppositeBoundaryFace);
489 return getParallelVelocityFromBoundary_(element, scvf, faceVars, currentScvfBoundaryTypes, lateralOppositeFaceBoundaryTypes, localOppositeSubFaceIdx);
490 }
491
511 const SubControlVolumeFace& boundaryScvf_(const int localSubFaceIdx) const
512 {
513 if (scvf_.hasHalfParallelNeighbor(localSubFaceIdx))
514 {
515 return fvGeometry_.scvf(scvf_.outsideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
516 }
517 else if (scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
518 {
519 /*
520 * ------------
521 * | xxxxxxxx o
522 * | xxxxxxxx o
523 * | xxxxxxxx o
524 * lllllllllll-bbbbbbbbbbb
525 * | yyyyyyyy p |
526 * | yyyyyyyy p |
527 * | yyyyyyyy p |
528 * -----------------------
529 *
530 * o: scvf_, l: lateralFace, p: parallelFace, b: returned scvf, x: scvf_ inside scv, y: lateralFace
531 * outside scv
532 */
533 const SubControlVolumeFace& lateralFace = fvGeometry_.scvf(scvf_.insideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
534 const SubControlVolumeFace& parallelFace = fvGeometry_.scvf(lateralFace.outsideScvIdx(), scvf_.localFaceIdx());
535
536 const auto& localLateralIdx = scvf_.pairData(localSubFaceIdx).localLateralFaceIdx;
537 const auto& localLateralOppositeIdx = (localLateralIdx % 2) ? (localLateralIdx - 1) : (localLateralIdx + 1);
538
539 return fvGeometry_.scvf(parallelFace.outsideScvIdx(), localLateralOppositeIdx);
540 }
541 else
542 {
543 DUNE_THROW(Dune::InvalidStateException, "The function boundaryScvf_ should only be called when hasHalfParallelNeighbor or hasCornerParallelNeighbor.");
544 }
545 }
546
566 Element boundaryElement_(const int localSubFaceIdx) const
567 {
568 if (scvf_.hasHalfParallelNeighbor(localSubFaceIdx))
569 {
570 return fvGeometry_.gridGeometry().element(scvf_.outsideScvIdx());
571 }
572 else if (scvf_.hasCornerParallelNeighbor(localSubFaceIdx))
573 {
574 const SubControlVolumeFace& lateralFace = fvGeometry_.scvf(scvf_.insideScvIdx(), scvf_.pairData(localSubFaceIdx).localLateralFaceIdx);
575 const SubControlVolumeFace& parallelFace = fvGeometry_.scvf(lateralFace.outsideScvIdx(), scvf_.localFaceIdx());
576
577 return fvGeometry_.gridGeometry().element(parallelFace.outsideScvIdx());
578 }
579 else
580 {
581 DUNE_THROW(Dune::InvalidStateException, "When entering boundaryElement_ scvf_ should have either hasHalfParallelNeighbor or hasCornerParallelNeighbor true. Not the case here.");
582 }
583 }
584
604 bool dirichletParallelNeighbor_(const int localSubFaceIdx) const
605 {
606 const auto& problem = elemVolVars_.gridVolVars().problem();
607 const Element& boundaryElement = boundaryElement_(localSubFaceIdx);
608 const SubControlVolumeFace& boundaryScvf = boundaryScvf_(localSubFaceIdx);
609
610 return problem.boundaryTypes(boundaryElement, boundaryScvf).isDirichlet(Indices::velocity(scvf_.directionIndex()));
611 }
612
631 Scalar getParallelVelocityFromCorner_(const int localSubFaceIdx) const
632 {
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);
637
638 return problem.dirichlet(boundaryElement, ghostFace)[Indices::velocity(scvf_.directionIndex())];
639 }
640
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_;
648};
649
650} // end namespace Dumux
651
652#endif
The upwinding variables class for the Navier-Stokes model using the staggered grid discretization.
Definition: staggeredupwindhelper.hh:36
FacePrimaryVariables computeUpwindFrontalMomentum(const bool selfIsUpstream) const
Returns the momentum in the frontal direction.
Definition: staggeredupwindhelper.hh:103
FacePrimaryVariables computeUpwindLateralMomentum(const bool selfIsUpstream, const SubControlVolumeFace &lateralFace, const int localSubFaceIdx, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, const std::optional< BoundaryTypes > &lateralFaceBoundaryTypes) const
Returns the momentum in the lateral direction.
Definition: staggeredupwindhelper.hh:133
StaggeredUpwindHelper(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const ElementFaceVariables &elemFaceVars, const ElementVolumeVariables &elemVolVars, const UpwindScheme &upwindScheme)
Definition: staggeredupwindhelper.hh:80
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:46
static Scalar beaversJosephVelocityAtLateralScvf(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf, const FaceVariables &faceVars, const std::optional< BoundaryTypes > &currentScvfBoundaryTypes, 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: staggered/velocitygradients.hh:311
Defines all properties used in Dumux.
Type traits for problem classes.
Some exceptions thrown in DuMux
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:30
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
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:62
Define some often used mathematical functions.
The available discretization methods in Dumux.
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:53
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
This file contains different higher order methods for approximating the velocity.
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::BoundaryTypes BoundaryTypes
Definition: common/typetraits/problem.hh:34