3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/navierstokes/staggered/fluxvariables.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_FLUXVARIABLES_HH
25#define DUMUX_NAVIERSTOKES_STAGGERED_FLUXVARIABLES_HH
26
27#include <array>
28#include <optional>
29#include <type_traits>
30
31#include <dumux/common/math.hh>
36
40
42#include "velocitygradients.hh"
43
44namespace Dumux {
45
51// forward declaration
52template<class TypeTag, class DiscretizationMethod>
53class NavierStokesFluxVariablesImpl;
54
55template<class TypeTag>
56class NavierStokesFluxVariablesImpl<TypeTag, DiscretizationMethods::Staggered>
57: public FluxVariablesBase<GetPropType<TypeTag, Properties::Problem>,
58 typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView,
59 typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView,
60 typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView>
61{
63
64 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
65 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
66 using VolumeVariables = typename GridVolumeVariables::VolumeVariables;
67
68 using GridFluxVariablesCache = typename GridVariables::GridFluxVariablesCache;
69 using FluxVariablesCache = typename GridFluxVariablesCache::FluxVariablesCache;
70
71 using GridFaceVariables = typename GridVariables::GridFaceVariables;
72 using ElementFaceVariables = typename GridFaceVariables::LocalView;
73 using FaceVariables = typename GridFaceVariables::FaceVariables;
74
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;
86 using Extrusion = Extrusion_t<GridGeometry>;
87 using CellCenterPrimaryVariables = GetPropType<TypeTag, Properties::CellCenterPrimaryVariables>;
89 using BoundaryTypes = typename ProblemTraits<Problem>::BoundaryTypes;
91
92 static constexpr bool normalizePressure = getPropValue<TypeTag, Properties::NormalizePressure>();
93
94 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
95
96 static constexpr auto cellCenterIdx = GridGeometry::cellCenterIdx();
97 static constexpr auto faceIdx = GridGeometry::faceIdx();
98
99 static constexpr int upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
100 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
101
102public:
103
105
114 template<class UpwindTerm>
115 [[deprecated("Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
116 static Scalar advectiveFluxForCellCenter(const Problem& problem,
117 const ElementVolumeVariables& elemVolVars,
118 const ElementFaceVariables& elemFaceVars,
119 const SubControlVolumeFace &scvf,
120 UpwindTerm upwindTerm)
121 {
122 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
123 const bool insideIsUpstream = scvf.directionSign() == sign(velocity);
124 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
125
126 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
127 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
128
129 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
130 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
131
132 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
133 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
134 * velocity * Extrusion::area(scvf) * scvf.directionSign();
135
136 return flux * extrusionFactor_(elemVolVars, scvf);
137 }
138
148 template<class UpwindTerm>
149 static Scalar advectiveFluxForCellCenter(const Problem& problem,
150 const FVElementGeometry& fvGeometry,
151 const ElementVolumeVariables& elemVolVars,
152 const ElementFaceVariables& elemFaceVars,
153 const SubControlVolumeFace &scvf,
154 UpwindTerm upwindTerm)
155 {
156 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
157 const bool insideIsUpstream = scvf.directionSign() == sign(velocity);
158 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
159
160 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
161 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
162
163 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
164 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
165
166 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
167 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
168 * velocity * Extrusion::area(fvGeometry, scvf) * scvf.directionSign();
169
170 return flux * extrusionFactor_(elemVolVars, scvf);
171 }
172
188 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
189 const Element& element,
190 const FVElementGeometry& fvGeometry,
191 const ElementVolumeVariables& elemVolVars,
192 const ElementFaceVariables& elemFaceVars,
193 const SubControlVolumeFace& scvf,
194 const FluxVariablesCache& fluxVarsCache)
195 {
196 // The advectively transported quantity (i.e density for a single-phase system).
197 auto upwindTerm = [](const auto& volVars) { return volVars.density(); };
198
199 // Call the generic flux function.
200 CellCenterPrimaryVariables result(0.0);
201 result[Indices::conti0EqIdx - ModelTraits::dim()] = advectiveFluxForCellCenter(problem, fvGeometry, elemVolVars, elemFaceVars, scvf, upwindTerm);
202
203 return result;
204 }
205
209 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
210 const Element& element,
211 const SubControlVolumeFace& scvf,
212 const FVElementGeometry& fvGeometry,
213 const ElementVolumeVariables& elemVolVars,
214 const ElementFaceVariables& elemFaceVars,
215 const GridFluxVariablesCache& gridFluxVarsCache)
216 {
217 return computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache) +
218 computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache);
219 }
220
238 FacePrimaryVariables computeFrontalMomentumFlux(const Problem& problem,
239 const Element& element,
240 const SubControlVolumeFace& scvf,
241 const FVElementGeometry& fvGeometry,
242 const ElementVolumeVariables& elemVolVars,
243 const ElementFaceVariables& elemFaceVars,
244 const GridFluxVariablesCache& gridFluxVarsCache)
245 {
246 FacePrimaryVariables frontalFlux(0.0);
247
248 // The velocities of the dof at interest and the one of the opposite scvf.
249 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
250 const Scalar velocityOpposite = elemFaceVars[scvf].velocityOpposite();
251 const auto& faceVars = elemFaceVars[scvf];
252
253 // Advective flux.
254 if (problem.enableInertiaTerms())
255 {
256 // Get the average velocity at the center of the element (i.e. the location of the staggered face).
257 const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5;
258 const bool selfIsUpstream = scvf.directionSign() != sign(transportingVelocity);
259
260 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
261 frontalFlux += upwindHelper.computeUpwindFrontalMomentum(selfIsUpstream)
262 * transportingVelocity * -1.0 * scvf.directionSign();
263 }
264
265 // The volume variables within the current element. We only require those (and none of neighboring elements)
266 // because the fluxes are calculated over the staggered face at the center of the element.
267 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
268
269 // Diffusive flux.
270 const Scalar velocityGrad_ii = VelocityGradients::velocityGradII(scvf, faceVars) * scvf.directionSign();
271
272 static const bool enableUnsymmetrizedVelocityGradient
273 = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
274 const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0;
275 frontalFlux -= factor * insideVolVars.effectiveViscosity() * velocityGrad_ii;
276
277 // The pressure term.
278 // If specified, the pressure can be normalized using the initial value on the scfv of interest.
279 // The scvf is used to normalize by the same value from the left and right side.
280 // Can potentially help to improve the condition number of the system matrix.
281 const Scalar pressure = normalizePressure ?
282 insideVolVars.pressure() - problem.initial(scvf)[Indices::pressureIdx]
283 : insideVolVars.pressure();
284
285 // Account for the orientation of the staggered face's normal outer normal vector
286 // (pointing in opposite direction of the scvf's one).
287 frontalFlux += pressure * -1.0 * scvf.directionSign();
288
289 // Account for the staggered face's area. For rectangular elements, this equals the area of the scvf
290 // our velocity dof of interest lives on but with adjusted centroid
291 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
292 FaceFrontalSubControlVolumeFace frontalFace(scv.center(), scvf.area());
293 return frontalFlux * Extrusion::area(fvGeometry, frontalFace) * insideVolVars.extrusionFactor();
294 }
295
313 FacePrimaryVariables computeLateralMomentumFlux(const Problem& problem,
314 const Element& element,
315 const SubControlVolumeFace& scvf,
316 const FVElementGeometry& fvGeometry,
317 const ElementVolumeVariables& elemVolVars,
318 const ElementFaceVariables& elemFaceVars,
319 const GridFluxVariablesCache& gridFluxVarsCache)
320 {
321 FacePrimaryVariables lateralFlux(0.0);
322 const auto& faceVars = elemFaceVars[scvf];
323 const std::size_t numSubFaces = scvf.pairData().size();
324
325 // If the current scvf is on a boundary, check if there is a Neumann BC for the stress in tangential direction.
326 // Create a boundaryTypes object (will be empty if not at a boundary).
327 std::optional<BoundaryTypes> currentScvfBoundaryTypes;
328 if (scvf.boundary())
329 currentScvfBoundaryTypes.emplace(problem.boundaryTypes(element, scvf));
330
331 // Account for all sub faces.
332 for (int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
333 {
334 const auto eIdx = scvf.insideScvIdx();
335 // Get the face normal to the face the dof lives on. The staggered sub face coincides with half of this lateral face.
336 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
337
338 // Create a boundaryTypes object (will be empty if not at a boundary).
339 std::optional<BoundaryTypes> lateralFaceBoundaryTypes;
340
341 // Check if there is face/element parallel to our face of interest where the dof lives on. If there is no parallel neighbor,
342 // we are on a boundary where we have to check for boundary conditions.
343 if (lateralFace.boundary())
344 {
345 // Retrieve the boundary types that correspond to the center of the lateral scvf. As a convention, we always query
346 // the type of BCs at the center of the element's "actual" lateral scvf (not the face of the staggered control volume).
347 // The value of the BC will be evaluated at the center of the staggered face.
348 // --------T######V || frontal face of staggered half-control-volume
349 // | || | current scvf # lateral staggered face of interest (may lie on a boundary)
350 // | || | x dof position
351 // | || x~~~~> vel.Self -- element boundaries
352 // | || | T position at which the type of boundary conditions will be evaluated
353 // | || | (center of lateral scvf)
354 // ---------------- V position at which the value of the boundary conditions will be evaluated
355 // (center of the staggered lateral face)
356 lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralFace));
357 }
358
359 // If the current scvf is on a boundary and if there is a Neumann or Beavers-Joseph-(Saffmann) BC for the stress in tangential direction,
360 // assign this value for the lateral momentum flux. No further calculations are required. We assume that all lateral faces
361 // have the same type of BC (Neumann or Beavers-Joseph-(Saffmann)), but we sample the value at their actual positions.
362 if (currentScvfBoundaryTypes)
363 {
364 // Handle Neumann BCs.
365 if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralFace.directionIndex())))
366 {
367 // Get the location of the lateral staggered face's center.
368 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
369 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
370 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralStaggeredFaceCenter);
371 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())]
372 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf) * lateralFace.directionSign();
373
374 continue;
375 }
376
377 // Treat the edge case when our current scvf has a BJ condition while the lateral face has a Neumann BC.
378 // It is not clear how to evaluate the BJ condition here.
379 // For symmetry reasons, our own scvf should then have the same Neumann flux as the lateral face.
380 // TODO: We should clarify if this is the correct approach.
381 if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())) && lateralFaceBoundaryTypes &&
382 lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
383 {
384 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
385 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
386 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter);
387 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
388 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf) * lateralFace.directionSign();
389
390 continue;
391 }
392 }
393
394 // Check if the lateral face (perpendicular to our current scvf) lies on a boundary. If yes, boundary conditions might need to be treated
395 // and further calculations can be skipped.
396 if (lateralFace.boundary())
397 {
398 // Check if we have a symmetry boundary condition. If yes, the tangential part of the momentum flux can be neglected
399 // and we may skip any further calculations for the given sub face.
400 if (lateralFaceBoundaryTypes->isSymmetry())
401 continue;
402
403 // Handle Neumann boundary conditions. No further calculations are then required for the given sub face.
404 if (lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
405 {
406 // Construct a temporary scvf which corresponds to the staggered sub face, featuring the location
407 // the staggered faces's center.
408 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
409 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
410 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter);
411 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
412 * elemVolVars[lateralFace.insideScvIdx()].extrusionFactor() * Extrusion::area(fvGeometry, lateralScvf);
413
414 continue;
415 }
416
417 // Handle wall-function fluxes (only required for RANS models)
418 if (incorporateWallFunction_(lateralFlux, problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, localSubFaceIdx))
419 continue;
420 }
421
422 // Check the consistency of the boundary conditions, exactly one of the following must be set
423 if (lateralFaceBoundaryTypes)
424 {
425 std::bitset<3> admittableBcTypes;
426 admittableBcTypes.set(0, lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
427 admittableBcTypes.set(1, lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())));
428 admittableBcTypes.set(2, lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())));
429 if (admittableBcTypes.count() != 1)
430 {
431 DUNE_THROW(Dune::InvalidStateException, "Invalid boundary conditions for lateral scvf "
432 "for the momentum equations at global position " << lateralStaggeredFaceCenter_(scvf, localSubFaceIdx)
433 << ", current scvf global position " << scvf.center());
434 }
435 }
436
437 // If none of the above boundary conditions apply for the given sub face, proceed to calculate the tangential momentum flux.
438 if (problem.enableInertiaTerms())
439 lateralFlux += computeAdvectivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
440 scvf, elemVolVars, elemFaceVars,
441 gridFluxVarsCache,
442 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
443 localSubFaceIdx);
444
445 lateralFlux += computeDiffusivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
446 scvf, elemVolVars, faceVars,
447 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
448 localSubFaceIdx);
449 }
450 return lateralFlux;
451 }
452
469 [[deprecated("Will be removed after release 3.6. Use interface with additional fvGeometry parameter instead.")]]
470 FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem& problem,
471 const Element& element,
472 const SubControlVolumeFace& scvf,
473 const ElementVolumeVariables& elemVolVars,
474 const ElementFaceVariables& elemFaceVars) const
475 {
476 FacePrimaryVariables inOrOutflow(0.0);
477 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
478
479 // Advective momentum flux.
480 if (problem.enableInertiaTerms())
481 {
482 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
483 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
484 const auto& upVolVars = (scvf.directionSign() == sign(velocitySelf)) ?
485 insideVolVars : outsideVolVars;
486
487 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
488 }
489
490 // Apply a pressure at the boundary.
491 const Scalar boundaryPressure = normalizePressure
492 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
493 problem.initial(scvf)[Indices::pressureIdx])
494 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
495 inOrOutflow += boundaryPressure;
496
497 // Account for the orientation of the face at the boundary,
498 return inOrOutflow * scvf.directionSign() * Extrusion::area(scvf) * insideVolVars.extrusionFactor();
499 }
500
517 FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem& problem,
518 const SubControlVolumeFace& scvf,
519 const FVElementGeometry& fvGeometry,
520 const ElementVolumeVariables& elemVolVars,
521 const ElementFaceVariables& elemFaceVars) const
522 {
523 FacePrimaryVariables inOrOutflow(0.0);
524 const auto& element = fvGeometry.element();
525 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
526
527 // Advective momentum flux.
528 if (problem.enableInertiaTerms())
529 {
530 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
531 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
532 const auto& upVolVars = (scvf.directionSign() == sign(velocitySelf)) ?
533 insideVolVars : outsideVolVars;
534
535 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
536 }
537
538 // Apply a pressure at the boundary.
539 const Scalar boundaryPressure = normalizePressure
540 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
541 problem.initial(scvf)[Indices::pressureIdx])
542 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
543 inOrOutflow += boundaryPressure;
544
545 // Account for the orientation of the face at the boundary,
546 return inOrOutflow * scvf.directionSign() * Extrusion::area(fvGeometry, scvf) * insideVolVars.extrusionFactor();
547 }
548
549private:
550
572 FacePrimaryVariables computeAdvectivePartOfLateralMomentumFlux_(const Problem& problem,
573 const FVElementGeometry& fvGeometry,
574 const Element& element,
575 const SubControlVolumeFace& scvf,
576 const ElementVolumeVariables& elemVolVars,
577 const ElementFaceVariables& elemFaceVars,
578 const GridFluxVariablesCache& gridFluxVarsCache,
579 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
580 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
581 const int localSubFaceIdx)
582 {
583 const auto eIdx = scvf.insideScvIdx();
584 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
585
586 // Get the transporting velocity, located at the scvf perpendicular to the current scvf where the dof
587 // of interest is located.
588 const Scalar transportingVelocity = [&]()
589 {
590 const auto& faceVars = elemFaceVars[scvf];
591 if (!scvf.boundary())
592 return faceVars.velocityLateralInside(localSubFaceIdx);
593 else
594 {
595 // Create a boundaryTypes object. Get the boundary conditions. We sample the type of BC at the center of the current scvf.
596 const auto bcTypes = problem.boundaryTypes(element, scvf);
597
598 if (bcTypes.isDirichlet(Indices::velocity(lateralFace.directionIndex())))
599 {
600 // Construct a temporary scvf which corresponds to the staggered sub face, featuring the location
601 // the staggered faces's center.
602 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
603 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralBoundaryFacePos);
604 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())];
605 }
606 else if (bcTypes.isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
607 {
608 return VelocityGradients::beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
609 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
610 }
611 else
612 return faceVars.velocityLateralInside(localSubFaceIdx);
613 }
614 }();
615
616 const bool selfIsUpstream = lateralFace.directionSign() == sign(transportingVelocity);
617 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
618 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
619 return upwindHelper.computeUpwindLateralMomentum(selfIsUpstream, lateralFace, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes)
620 * transportingVelocity * lateralFace.directionSign() * Extrusion::area(fvGeometry, lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
621 }
622
645 FacePrimaryVariables computeDiffusivePartOfLateralMomentumFlux_(const Problem& problem,
646 const FVElementGeometry& fvGeometry,
647 const Element& element,
648 const SubControlVolumeFace& scvf,
649 const ElementVolumeVariables& elemVolVars,
650 const FaceVariables& faceVars,
651 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
652 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
653 const int localSubFaceIdx)
654 {
655 const auto eIdx = scvf.insideScvIdx();
656 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
657
658 FacePrimaryVariables lateralDiffusiveFlux(0.0);
659
660 static const bool enableUnsymmetrizedVelocityGradient
661 = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
662
663 // Get the volume variables of the own and the neighboring element. The neighboring
664 // element is adjacent to the staggered face normal to the current scvf
665 // where the dof of interest is located.
666 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
667 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
668
669 // Get the averaged viscosity at the staggered face normal to the current scvf.
670 const Scalar muAvg = lateralFace.boundary()
671 ? insideVolVars.effectiveViscosity()
672 : (insideVolVars.effectiveViscosity() + outsideVolVars.effectiveViscosity()) * 0.5;
673
674 // Consider the shear stress caused by the gradient of the velocities normal to our face of interest.
675 if (!enableUnsymmetrizedVelocityGradient)
676 {
677 if (!scvf.boundary() ||
678 currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralFace.directionIndex())) ||
679 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
680 {
681 const Scalar velocityGrad_ji = VelocityGradients::velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
682 // Account for the orientation of the staggered normal face's outer normal vector.
683 lateralDiffusiveFlux -= muAvg * velocityGrad_ji * lateralFace.directionSign();
684 }
685 }
686
687 // Consider the shear stress caused by the gradient of the velocities parallel to our face of interest.
688 // If we have a Dirichlet condition for the pressure at the lateral face we assume to have a zero velocityGrad_ij velocity gradient
689 // so we can skip the computation.
690 if (!lateralFace.boundary() || !lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx))
691 {
692 const Scalar velocityGrad_ij = VelocityGradients::velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
693 lateralDiffusiveFlux -= muAvg * velocityGrad_ij * lateralFace.directionSign();
694 }
695
696 // Account for the area of the staggered lateral face (0.5 of the coinciding scfv).
697 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
698 return lateralDiffusiveFlux * Extrusion::area(fvGeometry, lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
699 }
700
715 const GlobalPosition& lateralStaggeredFaceCenter_(const SubControlVolumeFace& scvf, const int localSubFaceIdx) const
716 {
717 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
718 };
719
734 GlobalPosition lateralStaggeredSCVFCenter_(const SubControlVolumeFace& lateralFace,
735 const SubControlVolumeFace& currentFace,
736 const int localSubFaceIdx) const
737 {
738 auto pos = lateralStaggeredFaceCenter_(currentFace, localSubFaceIdx) + lateralFace.center();
739 pos *= 0.5;
740 return pos;
741 }
742
744 static Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf)
745 {
746 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
747 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
748 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
749 }
750
752 template<class ...Args, bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<!turbulenceModel, int> = 0>
753 bool incorporateWallFunction_(Args&&... args) const
754 { return false; }
755
757 template<bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<turbulenceModel, int> = 0>
758 bool incorporateWallFunction_(FacePrimaryVariables& lateralFlux,
759 const Problem& problem,
760 const Element& element,
761 const FVElementGeometry& fvGeometry,
762 const SubControlVolumeFace& scvf,
763 const ElementVolumeVariables& elemVolVars,
764 const ElementFaceVariables& elemFaceVars,
765 const std::size_t localSubFaceIdx) const
766 {
767 const auto eIdx = scvf.insideScvIdx();
768 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
769
770 if (problem.useWallFunction(element, lateralFace, Indices::velocity(scvf.directionIndex())))
771 {
772 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
773 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter_(scvf, localSubFaceIdx));
774 lateralFlux += problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
775 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(fvGeometry, lateralScvf);
776 return true;
777 }
778 else
779 return false;
780 }
781};
782
783} // end namespace Dumux
784
785#endif
Some exceptions thrown in DuMux
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Base class for the flux variables living on a sub control volume face.
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:69
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
typename Detail::template ProblemTraits< Problem, typename GridGeometry::DiscretizationMethod >::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
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: momentum/velocitygradients.hh:58
static Scalar advectiveFluxForCellCenter(const Problem &problem, const FVElementGeometry &fvGeometry, 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:149
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:313
FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars) const
Returns the momentum flux over an inflow or outflow boundary face.
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:517
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:238
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:116
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:209
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:470
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:104
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:188
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 direction.
Definition: staggeredupwindhelper.hh:115
Declares all properties used in Dumux.
Type traits for problem classes.