3.4
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, DiscretizationMethod discMethod>
53class NavierStokesFluxVariablesImpl;
54
55template<class TypeTag>
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 static Scalar advectiveFluxForCellCenter(const Problem& problem,
116 const ElementVolumeVariables& elemVolVars,
117 const ElementFaceVariables& elemFaceVars,
118 const SubControlVolumeFace &scvf,
119 UpwindTerm upwindTerm)
120 {
121 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
122 const bool insideIsUpstream = scvf.directionSign() == sign(velocity);
123 static const Scalar upwindWeight = getParamFromGroup<Scalar>(problem.paramGroup(), "Flux.UpwindWeight");
124
125 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
126 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
127
128 const auto& upstreamVolVars = insideIsUpstream ? insideVolVars : outsideVolVars;
129 const auto& downstreamVolVars = insideIsUpstream ? outsideVolVars : insideVolVars;
130
131 const Scalar flux = (upwindWeight * upwindTerm(upstreamVolVars) +
132 (1.0 - upwindWeight) * upwindTerm(downstreamVolVars))
133 * velocity * Extrusion::area(scvf) * scvf.directionSign();
134
135 return flux * extrusionFactor_(elemVolVars, scvf);
136 }
137
153 CellCenterPrimaryVariables computeMassFlux(const Problem& problem,
154 const Element& element,
155 const FVElementGeometry& fvGeometry,
156 const ElementVolumeVariables& elemVolVars,
157 const ElementFaceVariables& elemFaceVars,
158 const SubControlVolumeFace& scvf,
159 const FluxVariablesCache& fluxVarsCache)
160 {
161 // The advectively transported quantity (i.e density for a single-phase system).
162 auto upwindTerm = [](const auto& volVars) { return volVars.density(); };
163
164 // Call the generic flux function.
165 CellCenterPrimaryVariables result(0.0);
166 result[Indices::conti0EqIdx - ModelTraits::dim()] = advectiveFluxForCellCenter(problem, elemVolVars, elemFaceVars, scvf, upwindTerm);
167
168 return result;
169 }
170
174 FacePrimaryVariables computeMomentumFlux(const Problem& problem,
175 const Element& element,
176 const SubControlVolumeFace& scvf,
177 const FVElementGeometry& fvGeometry,
178 const ElementVolumeVariables& elemVolVars,
179 const ElementFaceVariables& elemFaceVars,
180 const GridFluxVariablesCache& gridFluxVarsCache)
181 {
182 return computeFrontalMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache) +
183 computeLateralMomentumFlux(problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, gridFluxVarsCache);
184 }
185
203 FacePrimaryVariables computeFrontalMomentumFlux(const Problem& problem,
204 const Element& element,
205 const SubControlVolumeFace& scvf,
206 const FVElementGeometry& fvGeometry,
207 const ElementVolumeVariables& elemVolVars,
208 const ElementFaceVariables& elemFaceVars,
209 const GridFluxVariablesCache& gridFluxVarsCache)
210 {
211 FacePrimaryVariables frontalFlux(0.0);
212
213 // The velocities of the dof at interest and the one of the opposite scvf.
214 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
215 const Scalar velocityOpposite = elemFaceVars[scvf].velocityOpposite();
216 const auto& faceVars = elemFaceVars[scvf];
217
218 // Advective flux.
219 if (problem.enableInertiaTerms())
220 {
221 // Get the average velocity at the center of the element (i.e. the location of the staggered face).
222 const Scalar transportingVelocity = (velocitySelf + velocityOpposite) * 0.5;
223 const bool selfIsUpstream = scvf.directionSign() != sign(transportingVelocity);
224
225 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
226 frontalFlux += upwindHelper.computeUpwindFrontalMomentum(selfIsUpstream)
227 * transportingVelocity * -1.0 * scvf.directionSign();
228 }
229
230 // The volume variables within the current element. We only require those (and none of neighboring elements)
231 // because the fluxes are calculated over the staggered face at the center of the element.
232 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
233
234 // Diffusive flux.
235 const Scalar velocityGrad_ii = VelocityGradients::velocityGradII(scvf, faceVars) * scvf.directionSign();
236
237 static const bool enableUnsymmetrizedVelocityGradient
238 = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
239 const Scalar factor = enableUnsymmetrizedVelocityGradient ? 1.0 : 2.0;
240 frontalFlux -= factor * insideVolVars.effectiveViscosity() * velocityGrad_ii;
241
242 // The pressure term.
243 // If specified, the pressure can be normalized using the initial value on the scfv of interest.
244 // The scvf is used to normalize by the same value from the left and right side.
245 // Can potentially help to improve the condition number of the system matrix.
246 const Scalar pressure = normalizePressure ?
247 insideVolVars.pressure() - problem.initial(scvf)[Indices::pressureIdx]
248 : insideVolVars.pressure();
249
250 // Account for the orientation of the staggered face's normal outer normal vector
251 // (pointing in opposite direction of the scvf's one).
252 frontalFlux += pressure * -1.0 * scvf.directionSign();
253
254 // Account for the staggered face's area. For rectangular elements, this equals the area of the scvf
255 // our velocity dof of interest lives on but with adjusted centroid
256 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
257 FaceFrontalSubControlVolumeFace frontalFace(scv.center(), scvf.area());
258 return frontalFlux * Extrusion::area(frontalFace) * insideVolVars.extrusionFactor();
259 }
260
278 FacePrimaryVariables computeLateralMomentumFlux(const Problem& problem,
279 const Element& element,
280 const SubControlVolumeFace& scvf,
281 const FVElementGeometry& fvGeometry,
282 const ElementVolumeVariables& elemVolVars,
283 const ElementFaceVariables& elemFaceVars,
284 const GridFluxVariablesCache& gridFluxVarsCache)
285 {
286 FacePrimaryVariables lateralFlux(0.0);
287 const auto& faceVars = elemFaceVars[scvf];
288 const std::size_t numSubFaces = scvf.pairData().size();
289
290 // If the current scvf is on a boundary, check if there is a Neumann BC for the stress in tangential direction.
291 // Create a boundaryTypes object (will be empty if not at a boundary).
292 std::optional<BoundaryTypes> currentScvfBoundaryTypes;
293 if (scvf.boundary())
294 currentScvfBoundaryTypes.emplace(problem.boundaryTypes(element, scvf));
295
296 // Account for all sub faces.
297 for (int localSubFaceIdx = 0; localSubFaceIdx < numSubFaces; ++localSubFaceIdx)
298 {
299 const auto eIdx = scvf.insideScvIdx();
300 // Get the face normal to the face the dof lives on. The staggered sub face conincides with half of this lateral face.
301 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
302
303 // Create a boundaryTypes object (will be empty if not at a boundary).
304 std::optional<BoundaryTypes> lateralFaceBoundaryTypes;
305
306 // Check if there is face/element parallel to our face of interest where the dof lives on. If there is no parallel neighbor,
307 // we are on a boundary where we have to check for boundary conditions.
308 if (lateralFace.boundary())
309 {
310 // Retrieve the boundary types that correspond to the center of the lateral scvf. As a convention, we always query
311 // the type of BCs at the center of the element's "actual" lateral scvf (not the face of the staggered control volume).
312 // The value of the BC will be evaluated at the center of the staggered face.
313 // --------T######V || frontal face of staggered half-control-volume
314 // | || | current scvf # lateral staggered face of interest (may lie on a boundary)
315 // | || | x dof position
316 // | || x~~~~> vel.Self -- element boundaries
317 // | || | T position at which the type of boundary conditions will be evaluated
318 // | || | (center of lateral scvf)
319 // ---------------- V position at which the value of the boundary conditions will be evaluated
320 // (center of the staggered lateral face)
321 lateralFaceBoundaryTypes.emplace(problem.boundaryTypes(element, lateralFace));
322 }
323
324 // If the current scvf is on a bounary and if there is a Neumann or Beavers-Joseph-(Saffmann) BC for the stress in tangential direction,
325 // assign this value for the lateral momentum flux. No further calculations are required. We assume that all lateral faces
326 // have the same type of BC (Neumann or Beavers-Joseph-(Saffmann)), but we sample the value at their actual positions.
327 if (currentScvfBoundaryTypes)
328 {
329 // Handle Neumann BCs.
330 if (currentScvfBoundaryTypes->isNeumann(Indices::velocity(lateralFace.directionIndex())))
331 {
332 // Get the location of the lateral staggered face's center.
333 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
334 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
335 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralStaggeredFaceCenter);
336 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())]
337 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf) * lateralFace.directionSign();
338
339 continue;
340 }
341
342 // Treat the edge case when our current scvf has a BJ condition while the lateral face has a Neumann BC.
343 // It is not clear how to evaluate the BJ condition here.
344 // For symmetry reasons, our own scvf should then have the same Neumann flux as the lateral face.
345 // TODO: We should clarify if this is the correct approach.
346 if (currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())) && lateralFaceBoundaryTypes &&
347 lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
348 {
349 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
350 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
351 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter);
352 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
353 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf) * lateralFace.directionSign();
354
355 continue;
356 }
357 }
358
359 // Check if the lateral face (perpendicular to our current scvf) lies on a boundary. If yes, boundary conditions might need to be treated
360 // and further calculations can be skipped.
361 if (lateralFace.boundary())
362 {
363 // Check if we have a symmetry boundary condition. If yes, the tangential part of the momentum flux can be neglected
364 // and we may skip any further calculations for the given sub face.
365 if (lateralFaceBoundaryTypes->isSymmetry())
366 continue;
367
368 // Handle Neumann boundary conditions. No further calculations are then required for the given sub face.
369 if (lateralFaceBoundaryTypes->isNeumann(Indices::velocity(scvf.directionIndex())))
370 {
371 // Construct a temporary scvf which corresponds to the staggered sub face, featuring the location
372 // the staggered faces's center.
373 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
374 const auto& lateralStaggeredFaceCenter = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
375 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter);
376 lateralFlux += problem.neumann(element, fvGeometry, elemVolVars, elemFaceVars, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
377 * elemVolVars[lateralFace.insideScvIdx()].extrusionFactor() * Extrusion::area(lateralScvf);
378
379 continue;
380 }
381
382 // Handle wall-function fluxes (only required for RANS models)
383 if (incorporateWallFunction_(lateralFlux, problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, localSubFaceIdx))
384 continue;
385 }
386
387 // Check the consistency of the boundary conditions, exactly one of the following must be set
388 if (lateralFaceBoundaryTypes)
389 {
390 std::bitset<3> admittableBcTypes;
391 admittableBcTypes.set(0, lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx));
392 admittableBcTypes.set(1, lateralFaceBoundaryTypes->isDirichlet(Indices::velocity(scvf.directionIndex())));
393 admittableBcTypes.set(2, lateralFaceBoundaryTypes->isBeaversJoseph(Indices::velocity(scvf.directionIndex())));
394 if (admittableBcTypes.count() != 1)
395 {
396 DUNE_THROW(Dune::InvalidStateException, "Invalid boundary conditions for lateral scvf "
397 "for the momentum equations at global position " << lateralStaggeredFaceCenter_(scvf, localSubFaceIdx)
398 << ", current scvf global position " << scvf.center());
399 }
400 }
401
402 // If none of the above boundary conditions apply for the given sub face, proceed to calculate the tangential momentum flux.
403 if (problem.enableInertiaTerms())
404 lateralFlux += computeAdvectivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
405 scvf, elemVolVars, elemFaceVars,
406 gridFluxVarsCache,
407 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
408 localSubFaceIdx);
409
410 lateralFlux += computeDiffusivePartOfLateralMomentumFlux_(problem, fvGeometry, element,
411 scvf, elemVolVars, faceVars,
412 currentScvfBoundaryTypes, lateralFaceBoundaryTypes,
413 localSubFaceIdx);
414 }
415 return lateralFlux;
416 }
417
434 FacePrimaryVariables inflowOutflowBoundaryFlux(const Problem& problem,
435 const Element& element,
436 const SubControlVolumeFace& scvf,
437 const ElementVolumeVariables& elemVolVars,
438 const ElementFaceVariables& elemFaceVars) const
439 {
440 FacePrimaryVariables inOrOutflow(0.0);
441 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
442
443 // Advective momentum flux.
444 if (problem.enableInertiaTerms())
445 {
446 const Scalar velocitySelf = elemFaceVars[scvf].velocitySelf();
447 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
448 const auto& upVolVars = (scvf.directionSign() == sign(velocitySelf)) ?
449 insideVolVars : outsideVolVars;
450
451 inOrOutflow += velocitySelf * velocitySelf * upVolVars.density();
452 }
453
454 // Apply a pressure at the boundary.
455 const Scalar boundaryPressure = normalizePressure
456 ? (problem.dirichlet(element, scvf)[Indices::pressureIdx] -
457 problem.initial(scvf)[Indices::pressureIdx])
458 : problem.dirichlet(element, scvf)[Indices::pressureIdx];
459 inOrOutflow += boundaryPressure;
460
461 // Account for the orientation of the face at the boundary,
462 return inOrOutflow * scvf.directionSign() * Extrusion::area(scvf) * insideVolVars.extrusionFactor();
463 }
464
465private:
466
488 FacePrimaryVariables computeAdvectivePartOfLateralMomentumFlux_(const Problem& problem,
489 const FVElementGeometry& fvGeometry,
490 const Element& element,
491 const SubControlVolumeFace& scvf,
492 const ElementVolumeVariables& elemVolVars,
493 const ElementFaceVariables& elemFaceVars,
494 const GridFluxVariablesCache& gridFluxVarsCache,
495 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
496 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
497 const int localSubFaceIdx)
498 {
499 const auto eIdx = scvf.insideScvIdx();
500 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
501
502 // Get the transporting velocity, located at the scvf perpendicular to the current scvf where the dof
503 // of interest is located.
504 const Scalar transportingVelocity = [&]()
505 {
506 const auto& faceVars = elemFaceVars[scvf];
507 if (!scvf.boundary())
508 return faceVars.velocityLateralInside(localSubFaceIdx);
509 else
510 {
511 // Create a boundaryTypes object. Get the boundary conditions. We sample the type of BC at the center of the current scvf.
512 const auto bcTypes = problem.boundaryTypes(element, scvf);
513
514 if (bcTypes.isDirichlet(Indices::velocity(lateralFace.directionIndex())))
515 {
516 // Construct a temporary scvf which corresponds to the staggered sub face, featuring the location
517 // the staggered faces's center.
518 const auto& lateralBoundaryFacePos = lateralStaggeredFaceCenter_(scvf, localSubFaceIdx);
519 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(scvf, lateralBoundaryFacePos);
520 return problem.dirichlet(element, lateralBoundaryFace)[Indices::velocity(lateralFace.directionIndex())];
521 }
522 else if (bcTypes.isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
523 {
524 return VelocityGradients::beaversJosephVelocityAtCurrentScvf(problem, element, fvGeometry, scvf, faceVars,
525 currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
526 }
527 else
528 return faceVars.velocityLateralInside(localSubFaceIdx);
529 }
530 }();
531
532 const bool selfIsUpstream = lateralFace.directionSign() == sign(transportingVelocity);
533 StaggeredUpwindHelper<TypeTag, upwindSchemeOrder> upwindHelper(element, fvGeometry, scvf, elemFaceVars, elemVolVars, gridFluxVarsCache.staggeredUpwindMethods());
534 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
535 return upwindHelper.computeUpwindLateralMomentum(selfIsUpstream, lateralFace, localSubFaceIdx, currentScvfBoundaryTypes, lateralFaceBoundaryTypes)
536 * transportingVelocity * lateralFace.directionSign() * Extrusion::area(lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
537 }
538
561 FacePrimaryVariables computeDiffusivePartOfLateralMomentumFlux_(const Problem& problem,
562 const FVElementGeometry& fvGeometry,
563 const Element& element,
564 const SubControlVolumeFace& scvf,
565 const ElementVolumeVariables& elemVolVars,
566 const FaceVariables& faceVars,
567 const std::optional<BoundaryTypes>& currentScvfBoundaryTypes,
568 const std::optional<BoundaryTypes>& lateralFaceBoundaryTypes,
569 const int localSubFaceIdx)
570 {
571 const auto eIdx = scvf.insideScvIdx();
572 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
573
574 FacePrimaryVariables lateralDiffusiveFlux(0.0);
575
576 static const bool enableUnsymmetrizedVelocityGradient
577 = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradient", false);
578
579 // Get the volume variables of the own and the neighboring element. The neighboring
580 // element is adjacent to the staggered face normal to the current scvf
581 // where the dof of interest is located.
582 const auto& insideVolVars = elemVolVars[lateralFace.insideScvIdx()];
583 const auto& outsideVolVars = elemVolVars[lateralFace.outsideScvIdx()];
584
585 // Get the averaged viscosity at the staggered face normal to the current scvf.
586 const Scalar muAvg = lateralFace.boundary()
587 ? insideVolVars.effectiveViscosity()
588 : (insideVolVars.effectiveViscosity() + outsideVolVars.effectiveViscosity()) * 0.5;
589
590 // Consider the shear stress caused by the gradient of the velocities normal to our face of interest.
591 if (!enableUnsymmetrizedVelocityGradient)
592 {
593 if (!scvf.boundary() ||
594 currentScvfBoundaryTypes->isDirichlet(Indices::velocity(lateralFace.directionIndex())) ||
595 currentScvfBoundaryTypes->isBeaversJoseph(Indices::velocity(lateralFace.directionIndex())))
596 {
597 const Scalar velocityGrad_ji = VelocityGradients::velocityGradJI(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
598 // Account for the orientation of the staggered normal face's outer normal vector.
599 lateralDiffusiveFlux -= muAvg * velocityGrad_ji * lateralFace.directionSign();
600 }
601 }
602
603 // Consider the shear stress caused by the gradient of the velocities parallel to our face of interest.
604 // If we have a Dirichlet condition for the pressure at the lateral face we assume to have a zero velocityGrad_ij velocity gradient
605 // so we can skip the computation.
606 if (!lateralFace.boundary() || !lateralFaceBoundaryTypes->isDirichlet(Indices::pressureIdx))
607 {
608 const Scalar velocityGrad_ij = VelocityGradients::velocityGradIJ(problem, element, fvGeometry, scvf, faceVars, currentScvfBoundaryTypes, lateralFaceBoundaryTypes, localSubFaceIdx);
609 lateralDiffusiveFlux -= muAvg * velocityGrad_ij * lateralFace.directionSign();
610 }
611
612 // Account for the area of the staggered lateral face (0.5 of the coinciding scfv).
613 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
614 return lateralDiffusiveFlux * Extrusion::area(lateralScvf) * extrusionFactor_(elemVolVars, lateralFace);
615 }
616
631 const GlobalPosition& lateralStaggeredFaceCenter_(const SubControlVolumeFace& scvf, const int localSubFaceIdx) const
632 {
633 return scvf.pairData(localSubFaceIdx).lateralStaggeredFaceCenter;
634 };
635
650 GlobalPosition lateralStaggeredSCVFCenter_(const SubControlVolumeFace& lateralFace,
651 const SubControlVolumeFace& currentFace,
652 const int localSubFaceIdx) const
653 {
654 auto pos = lateralStaggeredFaceCenter_(currentFace, localSubFaceIdx) + lateralFace.center();
655 pos *= 0.5;
656 return pos;
657 }
658
660 static Scalar extrusionFactor_(const ElementVolumeVariables& elemVolVars, const SubControlVolumeFace& scvf)
661 {
662 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
663 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
664 return harmonicMean(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor());
665 }
666
668 template<class ...Args, bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<!turbulenceModel, int> = 0>
669 bool incorporateWallFunction_(Args&&... args) const
670 { return false; }
671
673 template<bool turbulenceModel = ModelTraits::usesTurbulenceModel(), std::enable_if_t<turbulenceModel, int> = 0>
674 bool incorporateWallFunction_(FacePrimaryVariables& lateralFlux,
675 const Problem& problem,
676 const Element& element,
677 const FVElementGeometry& fvGeometry,
678 const SubControlVolumeFace& scvf,
679 const ElementVolumeVariables& elemVolVars,
680 const ElementFaceVariables& elemFaceVars,
681 const std::size_t localSubFaceIdx) const
682 {
683 const auto eIdx = scvf.insideScvIdx();
684 const auto& lateralFace = fvGeometry.scvf(eIdx, scvf.pairData(localSubFaceIdx).localLateralFaceIdx);
685
686 if (problem.useWallFunction(element, lateralFace, Indices::velocity(scvf.directionIndex())))
687 {
688 FaceLateralSubControlVolumeFace lateralScvf(lateralStaggeredSCVFCenter_(lateralFace, scvf, localSubFaceIdx), 0.5*lateralFace.area());
689 const auto lateralBoundaryFace = makeStaggeredBoundaryFace(lateralFace, lateralStaggeredFaceCenter_(scvf, localSubFaceIdx));
690 lateralFlux += problem.wallFunction(element, fvGeometry, elemVolVars, elemFaceVars, scvf, lateralBoundaryFace)[Indices::velocity(scvf.directionIndex())]
691 * extrusionFactor_(elemVolVars, lateralFace) * Extrusion::area(lateralScvf);
692 return true;
693 }
694 else
695 return false;
696 }
697};
698
699} // end namespace Dumux
700
701#endif
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.
Helper classes to compute the integration elements.
Base class for the flux variables living on a sub control volume face.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
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
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
typename Detail::template ProblemTraits< Problem, GridGeometry::discMethod >::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
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:278
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:203
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:115
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:174
GetPropType< TypeTag, Properties::HeatConductionType > HeatConductionType
Definition: freeflow/navierstokes/staggered/fluxvariables.hh:104
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:434
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:153
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 directon.
Definition: staggeredupwindhelper.hh:115
Helper class for calculating the velocity gradients for the Navier-Stokes model using the staggered g...
Definition: velocitygradients.hh:39
Declares all properties used in Dumux.
Type traits for problem classes.