Line | Branch | Exec | Source |
---|---|---|---|
1 | // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*- | ||
2 | // vi: set et ts=4 sw=4 sts=4: | ||
3 | // | ||
4 | // SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup NavierStokesModel | ||
10 | * \brief Navier Stokes slip condition | ||
11 | */ | ||
12 | #ifndef DUMUX_FREEFLOW_NAVIERSTOKES_SLIPCONDITION_HH | ||
13 | #define DUMUX_FREEFLOW_NAVIERSTOKES_SLIPCONDITION_HH | ||
14 | |||
15 | #include <dumux/common/tag.hh> | ||
16 | #include <dumux/discretization/method.hh> | ||
17 | |||
18 | namespace Dumux::NavierStokes::SlipConditions { | ||
19 | |||
20 | /*! | ||
21 | * \ingroup NavierStokesModel | ||
22 | * \brief Tag for the Beavers-Joseph slip condition | ||
23 | */ | ||
24 | struct BJ : public Utility::Tag<BJ> { | ||
25 | static std::string name() { return "Beavers-Joseph"; } | ||
26 | }; | ||
27 | |||
28 | /*! | ||
29 | * \ingroup NavierStokesModel | ||
30 | * \brief Tag for the Beavers-Joseph-Saffman slip condition | ||
31 | */ | ||
32 | struct BJS : public Utility::Tag<BJS> { | ||
33 | static std::string name() { return "Beavers-Joseph-Saffman"; } | ||
34 | }; | ||
35 | |||
36 | /*! | ||
37 | * \ingroup NavierStokesModel | ||
38 | * \brief Tag for the Beavers-Joseph slip condition | ||
39 | */ | ||
40 | inline constexpr BJ bj{}; | ||
41 | |||
42 | /*! | ||
43 | * \ingroup NavierStokesModel | ||
44 | * \brief Tag for the Beavers-Joseph-Saffman slip condition | ||
45 | */ | ||
46 | inline constexpr BJS bjs{}; | ||
47 | |||
48 | } // end namespace Dumux::NavierStokes::SlipConditions | ||
49 | |||
50 | namespace Dumux { | ||
51 | |||
52 | /*! | ||
53 | * \ingroup NavierStokesModel | ||
54 | * \brief Navier Stokes slip velocity policy | ||
55 | */ | ||
56 | template<class DiscretizationMethod, class SlipCondition> | ||
57 | class NavierStokesSlipVelocity; | ||
58 | |||
59 | /*! | ||
60 | * \ingroup NavierStokesModel | ||
61 | * \brief Navier Stokes slip velocity helper for fcstaggered discretization | ||
62 | * | ||
63 | * For now, this class implements the Beavers-Joseph or the Beavers-Joseph-Saffman condition | ||
64 | * which models the slip velocity at a porous boundary. The condition is chosen by passing | ||
65 | * the corresponding SlipCondition tag to the class. | ||
66 | */ | ||
67 | template<class SlipCondition> | ||
68 | struct NavierStokesSlipVelocity<DiscretizationMethods::FCStaggered, SlipCondition> | ||
69 | { | ||
70 | static constexpr SlipCondition slipCondition{}; | ||
71 | |||
72 | /*! | ||
73 | * \brief Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition. | ||
74 | * \note This only returns a vector filled with one component of the slip velocity (corresponding to the dof axis of the scv the svf belongs to) | ||
75 | * \param problem The problem | ||
76 | * \param fvGeometry The finite-volume geometry | ||
77 | * \param scvf The sub control volume face | ||
78 | * \param elemVolVars The volume variables for the element | ||
79 | * \param tangentialVelocityDeriv Pre-calculated velocity derivative | ||
80 | */ | ||
81 | template<class Problem, class FVElementGeometry, class ElementVolumeVariables, class Scalar> | ||
82 | 147310 | static auto velocity(const Problem& problem, | |
83 | const FVElementGeometry& fvGeometry, | ||
84 | const typename FVElementGeometry::SubControlVolumeFace& scvf, | ||
85 | const ElementVolumeVariables& elemVolVars, | ||
86 | Scalar tangentialVelocityDeriv) | ||
87 | { | ||
88 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 147310 times.
|
147310 | assert(scvf.isLateral()); |
89 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 147310 times.
|
147310 | assert(scvf.boundary()); |
90 | |||
91 | static constexpr int dimWorld = FVElementGeometry::GridGeometry::GridView::dimensionworld; | ||
92 | using Vector = Dune::FieldVector<Scalar, dimWorld>; | ||
93 | |||
94 | 147310 | Vector porousMediumVelocity(0.0); | |
95 | |||
96 | if constexpr (slipCondition == NavierStokes::SlipConditions::bj) | ||
97 | 3606 | porousMediumVelocity = problem.porousMediumVelocity(fvGeometry, scvf); | |
98 | else if (!(slipCondition == NavierStokes::SlipConditions::bjs)) | ||
99 | DUNE_THROW(Dune::NotImplemented, "Fcstaggered currently only implements BJ or BJS slip conditions"); | ||
100 | |||
101 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 147310 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 147310 times.
|
294620 | const auto& scv = fvGeometry.scv(scvf.insideScvIdx()); |
102 | |||
103 | // create a unit normal vector oriented in positive coordinate direction | ||
104 | 147310 | Vector tangent(0.0); | |
105 | 147310 | tangent[scv.dofAxis()] = 1.0; | |
106 | |||
107 | // du/dy + dv/dx = beta * (u_boundary-uPM) | ||
108 | // beta = alpha/sqrt(K) | ||
109 | 147310 | const Scalar betaBJ = problem.betaBJ(fvGeometry, scvf, tangent); | |
110 | 589240 | const Scalar distanceNormalToBoundary = (scvf.ipGlobal() - scv.dofPosition()).two_norm(); | |
111 | |||
112 |
5/8✓ Branch 0 taken 16 times.
✓ Branch 1 taken 147294 times.
✓ Branch 3 taken 16 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 16 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 16 times.
✗ Branch 10 not taken.
|
147310 | static const bool onlyNormalGradient = getParamFromGroup<bool>(problem.paramGroup(), "FreeFlow.EnableUnsymmetrizedVelocityGradientForBeaversJoseph", false); |
113 |
2/2✓ Branch 0 taken 63680 times.
✓ Branch 1 taken 83630 times.
|
147310 | if (onlyNormalGradient) |
114 | 63680 | tangentialVelocityDeriv = 0.0; | |
115 | |||
116 | 294620 | const Scalar scalarSlipVelocity = (tangentialVelocityDeriv*distanceNormalToBoundary | |
117 | 147310 | + porousMediumVelocity * tangent * betaBJ * distanceNormalToBoundary | |
118 | 147310 | + elemVolVars[scv].velocity()) / (betaBJ*distanceNormalToBoundary + 1.0); | |
119 | |||
120 | 294620 | return scalarSlipVelocity*tangent; | |
121 | } | ||
122 | }; | ||
123 | |||
124 | } // end namespace Dumux | ||
125 | |||
126 | #endif | ||
127 |