3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
staggeredupwindfluxvariables.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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
26
27#include <array>
28#include <optional>
29
30#include <dumux/common/math.hh>
34
37#include "velocitygradients.hh"
38
39namespace Dumux {
40
45template<class TypeTag, int upwindSchemeOrder>
47{
49
50 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
51 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
52 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
53
54 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
55 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
56
57 using GridFaceVariables = typename GridVariables::GridFaceVariables;
58 using ElementFaceVariables = typename GridFaceVariables::LocalView;
59 using FaceVariables = typename GridFaceVariables::FaceVariables;
60
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;
71 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
75
76 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
77 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
78
79public:
88 static FacePrimaryVariables computeUpwindedFrontalMomentum(const SubControlVolumeFace& scvf,
89 const ElementFaceVariables& elemFaceVars,
90 const ElementVolumeVariables& elemVolVars,
91 const GridFluxVariablesCache& gridFluxVarsCache,
92 const Scalar transportingVelocity)
93 {
94 const bool canHigherOrder = canFrontalSecondOrder_(scvf, transportingVelocity, std::integral_constant<bool, useHigherOrder>{});
95
96 if (canHigherOrder)
97 {
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);
101 }
102 else
103 {
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);
107 }
108 }
109
119 static FacePrimaryVariables computeUpwindedLateralMomentum(const Problem& problem,
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)
129 {
130 const auto eIdx = scvf.insideScvIdx();
131 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
132
133 // Get the transporting velocity, located at the scvf perpendicular to the current scvf where the dof
134 // of interest is located.
135 const Scalar transportingVelocity = faceVars.velocityLateralInside(localSubFaceIdx);
136
137 // Check whether the own or the neighboring element is upstream.
138 const bool selfIsUpstream = ( lateralFace.directionSign() == sign(transportingVelocity) );
139
140 const bool canHigherOrder = canLateralSecondOrder_(scvf, selfIsUpstream, localSubFaceIdx, std::integral_constant<bool, useHigherOrder>{});
141
142 if (canHigherOrder)
143 {
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);
148 }
149 else
150 {
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);
155 }
156 }
157
158private:
159
169 static bool canFrontalSecondOrder_(const SubControlVolumeFace& ownScvf,
170 const Scalar transportingVelocity,
171 std::false_type)
172 { return false; }
173
187 static bool canFrontalSecondOrder_(const SubControlVolumeFace& ownScvf,
188 const Scalar transportingVelocity,
189 std::true_type)
190 {
191 // Depending on selfIsUpstream I have to check if I have a forward or a backward neighbor to retrieve
192 const bool selfIsUpstream = ownScvf.directionSign() != sign(transportingVelocity);
193 return selfIsUpstream ? ownScvf.hasForwardNeighbor(0) : ownScvf.hasBackwardNeighbor(0);
194 }
195
205 static std::array<Scalar, 2> getFrontalUpwindingMomenta_(const SubControlVolumeFace& scvf,
206 const ElementFaceVariables& elemFaceVars,
207 const Scalar density,
208 const Scalar transportingVelocity,
209 std::false_type)
210 {
211 const Scalar momentumSelf = elemFaceVars[scvf].velocitySelf() * density;
212 const Scalar momentumOpposite = elemFaceVars[scvf].velocityOpposite() * density;
213 const bool selfIsUpstream = scvf.directionSign() != sign(transportingVelocity);
214
215 return selfIsUpstream ? std::array<Scalar, 2>{momentumOpposite, momentumSelf}
216 : std::array<Scalar, 2>{momentumSelf, momentumOpposite};
217 }
218
228 static std::array<Scalar, 3> getFrontalUpwindingMomenta_(const SubControlVolumeFace& scvf,
229 const ElementFaceVariables& elemFaceVars,
230 const Scalar density,
231 const Scalar transportingVelocity,
232 std::true_type)
233 {
234 const bool selfIsUpstream = scvf.directionSign() != sign(transportingVelocity);
235 std::array<Scalar, 3> momenta;
236
237 if (selfIsUpstream)
238 {
239 momenta[0] = elemFaceVars[scvf].velocityOpposite() * density;
240 momenta[1] = elemFaceVars[scvf].velocitySelf() * density;
241 momenta[2] = elemFaceVars[scvf].velocityForward(0) * density;
242 }
243 else
244 {
245 momenta[0] = elemFaceVars[scvf].velocitySelf() * density;
246 momenta[1] = elemFaceVars[scvf].velocityOpposite() * density;
247 momenta[2] = elemFaceVars[scvf].velocityBackward(0) * density;
248 }
249
250 return momenta;
251 }
252
261 static Scalar doFrontalMomentumUpwinding_(const SubControlVolumeFace& scvf,
262 const std::array<Scalar, 2>& momenta,
263 const Scalar transportingVelocity,
264 const GridFluxVariablesCache& gridFluxVarsCache)
265 {
266 const auto& upwindScheme = gridFluxVarsCache.staggeredUpwindMethods();
267 return upwindScheme.upwind(momenta[0], momenta[1]);
268 }
269
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)
283 {
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());
288 }
289
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)
299 {
300 // Depending on selfIsUpstream the downstream and the (up)upstream distances are saved.
301 // distances {upstream to downstream distance, up-upstream to upstream distance, downstream staggered cell size}
302 std::array<Scalar, 3> distances;
303
304 if (selfIsUpstream)
305 {
306 distances[0] = ownScvf.selfToOppositeDistance();
307 distances[1] = ownScvf.axisData().inAxisForwardDistances[0];
308 distances[2] = 0.5 * (ownScvf.axisData().selfToOppositeDistance + ownScvf.axisData().inAxisBackwardDistances[0]);
309 }
310 else
311 {
312 distances[0] = ownScvf.selfToOppositeDistance();
313 distances[1] = ownScvf.axisData().inAxisBackwardDistances[0];
314 distances[2] = 0.5 * (ownScvf.axisData().selfToOppositeDistance + ownScvf.axisData().inAxisForwardDistances[0]);
315 }
316 return distances;
317 }
318
330 static bool canLateralSecondOrder_(const SubControlVolumeFace& ownScvf,
331 const bool selfIsUpstream,
332 const int localSubFaceIdx,
333 std::true_type)
334 {
335 if (selfIsUpstream)
336 {
337 // The self velocity is upstream. The downstream velocity can be assigned or retrieved
338 // from the boundary, even if there is no parallel neighbor.
339 return true;
340 }
341 else
342 {
343 // The self velocity is downstream. If there is no parallel neighbor I cannot use a second order approximation.
344 return ownScvf.hasParallelNeighbor(localSubFaceIdx, 0);
345 }
346 }
347
352 static bool canLateralSecondOrder_(const SubControlVolumeFace& ownScvf,
353 const bool selfIsUpstream,
354 const int localSubFaceIdx,
355 std::false_type)
356 { return false; }
357
358
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,
374 std::true_type)
375 {
376 const SubControlVolumeFace& lateralFace = fvGeometry.scvf(ownScvf.insideScvIdx(), ownScvf.pairData(localSubFaceIdx).localLateralFaceIdx);
377
378 // Get the volume variables of the own and the neighboring element
379 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
380 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
381
382 // If the lateral face lies on a boundary, we assume that the parallel velocity on the boundary is actually known,
383 // thus we always use this value for the computation of the transported momentum.
384 if (!ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
385 {
386 const Scalar boundaryMomentum = getParallelVelocityFromBoundary_(problem, element, fvGeometry, ownScvf,
387 faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
388 localSubFaceIdx) * insideVolVars.density();
389
390 return std::array<Scalar, 3>{boundaryMomentum, boundaryMomentum, boundaryMomentum};
391 }
392
393 // Check whether the own or the neighboring element is upstream.
394 const bool selfIsUpstream = lateralFace.directionSign() == sign(transportingVelocity);
395
396 std::array<Scalar, 3> momenta;
397
398 if (selfIsUpstream)
399 {
400 momenta[1] = faceVars.velocitySelf() * insideVolVars.density();
401
402 if (ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
403 momenta[0] = faceVars.velocityParallel(localSubFaceIdx, 0) * insideVolVars.density();
404 else
405 momenta[0] = getParallelVelocityFromBoundary_(problem, element, fvGeometry, ownScvf,
406 faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
407 localSubFaceIdx) * insideVolVars.density();
408
409 // The local index of the faces that is opposite to localSubFaceIdx
410 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
411
412 // The "upstream-upstream" velocity is retrieved from the other parallel neighbor or from the boundary.
413 if (ownScvf.hasParallelNeighbor(oppositeSubFaceIdx, 0))
414 momenta[2] = faceVars.velocityParallel(oppositeSubFaceIdx, 0) * insideVolVars.density();
415 else
416 momenta[2] = getParallelVelocityFromOppositeBoundary_(problem, element, fvGeometry, ownScvf,
417 faceVars, currentScvfBoundaryTypes,
418 oppositeSubFaceIdx) * insideVolVars.density();
419 }
420 else
421 {
422 momenta[0] = faceVars.velocitySelf() * outsideVolVars.density();
423 momenta[1] = faceVars.velocityParallel(localSubFaceIdx, 0) * outsideVolVars.density();
424
425 // If there is another parallel neighbor I can assign the "upstream-upstream" velocity, otherwise I retrieve it from the boundary.
426 if (ownScvf.hasParallelNeighbor(localSubFaceIdx, 1))
427 momenta[2] = faceVars.velocityParallel(localSubFaceIdx, 1) * outsideVolVars.density();
428 else
429 {
430 const Element& elementParallel = fvGeometry.gridGeometry().element(lateralFace.outsideScvIdx());
431 const SubControlVolumeFace& firstParallelScvf = fvGeometry.scvf(lateralFace.outsideScvIdx(), ownScvf.localFaceIdx());
432
433 momenta[2] = getParallelVelocityFromOppositeBoundary_(problem, elementParallel, fvGeometry, firstParallelScvf,
434 faceVars, problem.boundaryTypes(elementParallel, firstParallelScvf),
435 localSubFaceIdx) * outsideVolVars.density();
436 }
437 }
438
439 return momenta;
440 }
441
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,
457 std::false_type)
458 {
459 // Check whether the own or the neighboring element is upstream.
460 const SubControlVolumeFace& lateralFace = fvGeometry.scvf(ownScvf.insideScvIdx(), ownScvf.pairData(localSubFaceIdx).localLateralFaceIdx);
461
462 // Get the volume variables of the own and the neighboring element
463 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
464 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
465
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,
470 localSubFaceIdx)
471 * insideVolVars.density());
472
473 // If the lateral face lies on a boundary, we assume that the parallel velocity on the boundary is actually known,
474 // thus we always use this value for the computation of the transported momentum.
475 if (!ownScvf.hasParallelNeighbor(localSubFaceIdx, 0))
476 return std::array<Scalar, 2>{momentumParallel, momentumParallel};
477
478 const bool selfIsUpstream = lateralFace.directionSign() == sign(transportingVelocity);
479 const Scalar momentumSelf = faceVars.velocitySelf() * insideVolVars.density();
480
481 return selfIsUpstream ? std::array<Scalar, 2>{momentumParallel, momentumSelf}
482 : std::array<Scalar, 2>{momentumSelf, momentumParallel};
483 }
484
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)
496 {
497 return doFrontalMomentumUpwinding_(scvf, momenta, transportingVelocity, gridFluxVarsCache);
498 }
499
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)
515 {
516 const auto eIdx = scvf.insideScvIdx();
517 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
518
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());
523 }
524
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)
536 {
537 // distances {upstream to downstream distance, up-upstream to upstream distance, downstream staggered cell size}
538 std::array<Scalar, 3> distances;
539
540 if (selfIsUpstream)
541 {
542 // The local index of the faces that is opposite to localSubFaceIdx
543 const int oppositeSubFaceIdx = localSubFaceIdx % 2 ? localSubFaceIdx - 1 : localSubFaceIdx + 1;
544
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];
549 else
550 distances[2] = ownScvf.area() / 2.0;
551 }
552 else
553 {
554 distances[0] = ownScvf.parallelDofsDistance(localSubFaceIdx, 0);
555 distances[1] = ownScvf.parallelDofsDistance(localSubFaceIdx, 1);
556 distances[2] = ownScvf.pairData(localSubFaceIdx).parallelCellWidths[0];
557 }
558
559 return distances;
560 }
561
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)
582 {
583 // Find out what boundary type is set on the lateral face
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();
589
590 // If there is a Dirichlet condition for the pressure we assume zero gradient for the velocity,
591 // so the velocity at the boundary equal to that on the scvf.
592 if (useZeroGradient)
593 return velocitySelf;
594
595 if (lateralFaceHasBJS)
596 return VelocityGradients::beaversJosephVelocityAtLateralScvf(problem, element, fvGeometry, scvf, faceVars,
597 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
598
599 else if(lateralFaceHasDirichletVelocity)
600 {
601 // ________________
602 // ---------------o || frontal face of staggered half-control-volume
603 // | || # current scvf # ghostFace of interest, lies on boundary
604 // | || # x dof position
605 // | || x~~~~> vel.Self -- element boundaries
606 // | || # __ domain boundary
607 // | || # o position at which the boundary conditions will be evaluated
608 // ---------------#
609
610 const auto ghostFace = makeParallelGhostFace_(scvf, localSubFaceIdx);
611 return problem.dirichlet(element, ghostFace)[Indices::velocity(scvf.directionIndex())];
612 }
613 else
614 {
615 // Neumann conditions are not well implemented
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());
617 }
618 }
619
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)
636 {
637 // A ghost subface at the boundary is created, featuring the location of the sub face's center
638 const SubControlVolumeFace& lateralOppositeScvf = fvGeometry.scvf(scvf.insideScvIdx(), scvf.pairData(localOppositeSubFaceIdx).localLateralFaceIdx);
639 GlobalPosition center = scvf.pairData(localOppositeSubFaceIdx).lateralStaggeredFaceCenter + lateralOppositeScvf.center();
640 center *= 0.5;
641
642 // lateral face # lateralFace currently being assembled
643 // --------######## || frontal face of staggered half-control-volume
644 // | || | current scvf % lateralOppositeBoundaryFace of interest, lies on boundary
645 // | || | x dof position
646 // | || x~~~~> vel.Self -- element boundaries
647 // | || | __ domain boundary
648 // | || | o position at which the boundary conditions will be evaluated
649 // --------%%%c%%%o
650 // ________________ c center of opposite boundary face
651
652
653 // Get the boundary types of the lateral opposite boundary face
654 const auto lateralOppositeFaceBoundaryTypes = problem.boundaryTypes(element, lateralOppositeScvf.makeBoundaryFace(center));
655
656 return getParallelVelocityFromBoundary_(problem, element, fvGeometry, scvf,
657 faceVars, currentScvfBoundaryTypes, lateralOppositeFaceBoundaryTypes,
658 localOppositeSubFaceIdx);
659 }
660
662 static SubControlVolumeFace makeParallelGhostFace_(const SubControlVolumeFace& ownScvf, const int localSubFaceIdx)
663 {
664 return ownScvf.makeBoundaryFace(ownScvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter);
665 };
666};
667
668} // end namespace Dumux
669
670#endif // DUMUX_NAVIERSTOKES_STAGGERED_UPWINDVARIABLES_HH
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
Definition: adapt.hh:29
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 > &currentScvfBoundaryTypes, 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 > &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: velocitygradients.hh:301
Declares all properties used in Dumux.