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 | * \copydoc Dumux::NavierStokesStaggeredProblem | ||
11 | */ | ||
12 | #ifndef DUMUX_NAVIERSTOKES_STAGGERED_PROBLEM_HH | ||
13 | #define DUMUX_NAVIERSTOKES_STAGGERED_PROBLEM_HH | ||
14 | |||
15 | #include <dune/common/exceptions.hh> | ||
16 | #include <dune/common/typetraits.hh> | ||
17 | #include <dumux/common/properties.hh> | ||
18 | #include <dumux/common/staggeredfvproblem.hh> | ||
19 | #include <dumux/discretization/method.hh> | ||
20 | #include <dumux/freeflow/navierstokes/momentum/boundarytypes.hh> | ||
21 | |||
22 | namespace Dumux { | ||
23 | |||
24 | /*! | ||
25 | * \ingroup NavierStokesModel | ||
26 | * \brief Navier-Stokes staggered problem base class. | ||
27 | * | ||
28 | * This implements gravity (if desired) and a function returning the temperature. | ||
29 | * Includes a specialized method used only by the staggered grid discretization. | ||
30 | */ | ||
31 | template<class TypeTag> | ||
32 | 4 | class NavierStokesStaggeredProblem : public StaggeredFVProblem<TypeTag> | |
33 | { | ||
34 | using ParentType = StaggeredFVProblem<TypeTag>; | ||
35 | using Implementation = GetPropType<TypeTag, Properties::Problem>; | ||
36 | |||
37 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
38 | using GridView = typename GridGeometry::GridView; | ||
39 | using Element = typename GridView::template Codim<0>::Entity; | ||
40 | |||
41 | using GridVariables = GetPropType<TypeTag, Properties::GridVariables>; | ||
42 | using GridFaceVariables = typename GridVariables::GridFaceVariables; | ||
43 | using ElementFaceVariables = typename GridFaceVariables::LocalView; | ||
44 | using FaceVariables = typename GridFaceVariables::FaceVariables; | ||
45 | using GridVolumeVariables = typename GridVariables::GridVolumeVariables; | ||
46 | using ElementVolumeVariables = typename GridVolumeVariables::LocalView; | ||
47 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
48 | |||
49 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
50 | using SubControlVolume = typename FVElementGeometry::SubControlVolume; | ||
51 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
52 | using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>; | ||
53 | using Indices = typename GetPropType<TypeTag, Properties::ModelTraits>::Indices; | ||
54 | |||
55 | enum { | ||
56 | dim = GridView::dimension, | ||
57 | dimWorld = GridView::dimensionworld | ||
58 | }; | ||
59 | |||
60 | using GlobalPosition = typename SubControlVolumeFace::GlobalPosition; | ||
61 | using VelocityVector = Dune::FieldVector<Scalar, dimWorld>; | ||
62 | using GravityVector = Dune::FieldVector<Scalar, dimWorld>; | ||
63 | |||
64 | public: | ||
65 | /*! | ||
66 | * \brief The constructor | ||
67 | * \param gridGeometry The finite volume grid geometry | ||
68 | * \param paramGroup The parameter group in which to look for runtime parameters first (default is "") | ||
69 | */ | ||
70 | 50 | NavierStokesStaggeredProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "") | |
71 | : ParentType(gridGeometry, paramGroup) | ||
72 |
2/6✓ Branch 2 taken 50 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 50 times.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
100 | , gravity_(0.0) |
73 | { | ||
74 |
3/4✓ Branch 1 taken 50 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2 times.
✓ Branch 4 taken 48 times.
|
50 | if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity")) |
75 | 4 | gravity_[dim-1] = -9.81; | |
76 | |||
77 |
1/2✓ Branch 1 taken 50 times.
✗ Branch 2 not taken.
|
50 | enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms"); |
78 | 50 | } | |
79 | |||
80 | |||
81 | /*! | ||
82 | * \brief Returns the acceleration due to gravity. | ||
83 | * | ||
84 | * If the <tt>Problem.EnableGravity</tt> parameter is true, this means | ||
85 | * \f$\boldsymbol{g} = ( 0,\dots,\ -9.81)^T \f$, else \f$\boldsymbol{g} = ( 0,\dots, 0)^T \f$ | ||
86 | */ | ||
87 | const GravityVector& gravity() const | ||
88 | 150537765 | { return gravity_; } | |
89 | |||
90 | /*! | ||
91 | * \brief Returns whether interia terms should be considered. | ||
92 | */ | ||
93 | ✗ | bool enableInertiaTerms() const | |
94 | ✗ | { return enableInertiaTerms_; } | |
95 | |||
96 | //! Applies the initial face solution (velocities on the faces). Specialization for staggered grid discretization. | ||
97 | template <class SolutionVector, class G = GridGeometry> | ||
98 | typename std::enable_if<G::discMethod == DiscretizationMethods::staggered, void>::type | ||
99 | ✗ | applyInitialFaceSolution(SolutionVector& sol, | |
100 | const SubControlVolumeFace& scvf, | ||
101 | const PrimaryVariables& initSol) const | ||
102 | { | ||
103 | 287840 | sol[GridGeometry::faceIdx()][scvf.dofIndex()][0] = initSol[Indices::velocity(scvf.directionIndex())]; | |
104 | ✗ | } | |
105 | |||
106 | /*! | ||
107 | * \brief An additional drag term can be included as source term for the momentum balance | ||
108 | * to mimic 3D flow behavior in 2D: | ||
109 | * \f[ | ||
110 | * f_{drag} = -(8 \mu / h^2)v | ||
111 | * \f] | ||
112 | * Here, \f$h\f$ corresponds to the extruded height that is | ||
113 | * bounded by the imaginary walls. See Flekkoy et al. (1995) \cite flekkoy1995a<BR> | ||
114 | * A value of 8.0 is used as a default factor, corresponding | ||
115 | * to the velocity profile at the center plane | ||
116 | * of the virtual height (maximum velocity). Setting this value to 12.0 corresponds | ||
117 | * to an depth-averaged velocity (Venturoli and Boek, 2006) \cite venturoli2006a. | ||
118 | */ | ||
119 | Scalar pseudo3DWallFriction(const Scalar velocity, | ||
120 | const Scalar viscosity, | ||
121 | const Scalar height, | ||
122 | const Scalar factor = 8.0) const | ||
123 | { | ||
124 | static_assert(dim == 2, "Pseudo 3D wall friction may only be used in 2D"); | ||
125 | return -factor * velocity * viscosity / (height*height); | ||
126 | } | ||
127 | |||
128 | //! Convenience function for staggered grid implementation. | ||
129 | template <class ElementVolumeVariables, class ElementFaceVariables, class G = GridGeometry> | ||
130 | typename std::enable_if<G::discMethod == DiscretizationMethods::staggered, Scalar>::type | ||
131 | pseudo3DWallFriction(const SubControlVolumeFace& scvf, | ||
132 | const ElementVolumeVariables& elemVolVars, | ||
133 | const ElementFaceVariables& elemFaceVars, | ||
134 | const Scalar height, | ||
135 | const Scalar factor = 8.0) const | ||
136 | { | ||
137 | const Scalar velocity = elemFaceVars[scvf].velocitySelf(); | ||
138 | const Scalar viscosity = elemVolVars[scvf.insideScvIdx()].effectiveViscosity(); | ||
139 | return pseudo3DWallFriction(velocity, viscosity, height, factor); | ||
140 | } | ||
141 | |||
142 | /*! | ||
143 | * \brief Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boundary condition | ||
144 | * | ||
145 | * This member function must be overloaded in the problem implementation, if the BJS boundary condition is used. | ||
146 | */ | ||
147 | ✗ | Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const | |
148 | { | ||
149 | ✗ | DUNE_THROW(Dune::NotImplemented, | |
150 | "When using the Beavers-Joseph-Saffman boundary condition, " | ||
151 | "the permeability must be returned in the actual problem" | ||
152 | ); | ||
153 | } | ||
154 | |||
155 | /*! | ||
156 | * \brief Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition | ||
157 | * | ||
158 | * This member function must be overloaded in the problem implementation, if the BJS boundary condition is used. | ||
159 | */ | ||
160 | ✗ | Scalar alphaBJ(const SubControlVolumeFace& scvf) const | |
161 | { | ||
162 | ✗ | DUNE_THROW(Dune::NotImplemented, | |
163 | "When using the Beavers-Joseph-Saffman boundary condition, " | ||
164 | "the alpha value must be returned in the actual problem" | ||
165 | ); | ||
166 | } | ||
167 | |||
168 | /*! | ||
169 | * \brief Returns the beta value which is the alpha value divided by the square root of the (scalar-valued) interface permeability. | ||
170 | */ | ||
171 | 146832 | Scalar betaBJ(const Element& element, const SubControlVolumeFace& scvf, const GlobalPosition& tangentialVector) const | |
172 | { | ||
173 | 146832 | const Scalar interfacePermeability = interfacePermeability_(element, scvf, tangentialVector); | |
174 | using std::sqrt; | ||
175 | 146832 | return asImp_().alphaBJ(scvf) / sqrt(interfacePermeability); | |
176 | } | ||
177 | |||
178 | /*! | ||
179 | * \brief Returns the velocity in the porous medium (which is 0 by default according to Saffmann). | ||
180 | */ | ||
181 | ✗ | VelocityVector porousMediumVelocity(const Element& element, const SubControlVolumeFace& scvf) const | |
182 | { | ||
183 | ✗ | return VelocityVector(0.0); | |
184 | } | ||
185 | |||
186 | /*! | ||
187 | * \brief Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition. | ||
188 | */ | ||
189 | 863288 | const Scalar beaversJosephVelocity(const Element& element, | |
190 | const SubControlVolume& scv, | ||
191 | const SubControlVolumeFace& ownScvf, | ||
192 | const SubControlVolumeFace& faceOnPorousBoundary, | ||
193 | const Scalar velocitySelf, | ||
194 | const Scalar tangentialVelocityGradient) const | ||
195 | { | ||
196 | // create a unit normal vector oriented in positive coordinate direction | ||
197 | 863288 | GlobalPosition orientation = ownScvf.unitOuterNormal(); | |
198 | 863288 | orientation[ownScvf.directionIndex()] = 1.0; | |
199 | |||
200 | // du/dy + dv/dx = alpha/sqrt(K) * (u_boundary-uPM) | ||
201 | // beta = alpha/sqrt(K) | ||
202 | 863288 | const Scalar betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary, orientation); | |
203 | 3453152 | const Scalar distanceNormalToBoundary = (faceOnPorousBoundary.center() - scv.center()).two_norm(); | |
204 | |||
205 | 863288 | return (tangentialVelocityGradient*distanceNormalToBoundary | |
206 | 1726576 | + asImp_().porousMediumVelocity(element, faceOnPorousBoundary) * orientation * betaBJ * distanceNormalToBoundary | |
207 | 1726576 | + velocitySelf) / (betaBJ*distanceNormalToBoundary + 1.0); | |
208 | } | ||
209 | |||
210 | private: | ||
211 | //! Returns a scalar permeability value at the coupling interface | ||
212 | 146832 | Scalar interfacePermeability_(const Element& element, const SubControlVolumeFace& scvf, const GlobalPosition& tangentialVector) const | |
213 | { | ||
214 | 146832 | const auto& K = asImp_().permeability(element, scvf); | |
215 | |||
216 | // use t*K*t for permeability tensors | ||
217 | if constexpr (Dune::IsNumber<std::decay_t<decltype(K)>>::value) | ||
218 | ✗ | return K; | |
219 | else | ||
220 | 146832 | return vtmv(tangentialVector, K, tangentialVector); | |
221 | } | ||
222 | |||
223 | //! Returns the implementation of the problem (i.e. static polymorphism) | ||
224 | Implementation &asImp_() | ||
225 | { return *static_cast<Implementation *>(this); } | ||
226 | |||
227 | //! \copydoc asImp_() | ||
228 | const Implementation &asImp_() const | ||
229 | { return *static_cast<const Implementation *>(this); } | ||
230 | |||
231 | GravityVector gravity_; | ||
232 | bool enableInertiaTerms_; | ||
233 | }; | ||
234 | |||
235 | } // end namespace Dumux | ||
236 | |||
237 | #endif | ||
238 |