3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/navierstokes/problem.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_PROBLEM_HH
25#define DUMUX_NAVIERSTOKES_PROBLEM_HH
26
27#include <dune/common/exceptions.hh>
28#include <dune/common/typetraits.hh>
32
33namespace Dumux {
34
36template<class TypeTag, DiscretizationMethod discMethod> struct NavierStokesParentProblemImpl;
37
38template<class TypeTag>
40{
42};
43
45template<class TypeTag>
47 typename NavierStokesParentProblemImpl<TypeTag,
49
58template<class TypeTag>
60{
61 using ParentType = NavierStokesParentProblem<TypeTag>;
62 using Implementation = GetPropType<TypeTag, Properties::Problem>;
63
65 using GridView = typename GridGeometry::GridView;
66 using Element = typename GridView::template Codim<0>::Entity;
67
69 using GridFaceVariables = typename GridVariables::GridFaceVariables;
70 using ElementFaceVariables = typename GridFaceVariables::LocalView;
71 using FaceVariables = typename GridFaceVariables::FaceVariables;
72 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
73 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
75
76 using FVElementGeometry = typename GridGeometry::LocalView;
77 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
78 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
81
82 enum {
83 dim = GridView::dimension,
84 dimWorld = GridView::dimensionworld
85 };
86
87 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
88 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
89 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
90
91public:
97 NavierStokesProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
98 : ParentType(gridGeometry, paramGroup)
99 , gravity_(0.0)
100 {
101 if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity"))
102 gravity_[dim-1] = -9.81;
103
104 enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms");
105 }
106
115 Scalar temperatureAtPos(const GlobalPosition &globalPos) const
116 { return asImp_().temperature(); }
117
123 Scalar temperature() const
124 { DUNE_THROW(Dune::NotImplemented, "temperature() method not implemented by the actual problem"); }
125
132 const GravityVector& gravity() const
133 { return gravity_; }
134
139 { return enableInertiaTerms_; }
140
142 template <class SolutionVector, class G = GridGeometry>
143 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, void>::type
144 applyInitialFaceSolution(SolutionVector& sol,
145 const SubControlVolumeFace& scvf,
146 const PrimaryVariables& initSol) const
147 {
148 sol[GridGeometry::faceIdx()][scvf.dofIndex()][0] = initSol[Indices::velocity(scvf.directionIndex())];
149 }
150
164 Scalar pseudo3DWallFriction(const Scalar velocity,
165 const Scalar viscosity,
166 const Scalar height,
167 const Scalar factor = 8.0) const
168 {
169 static_assert(dim == 2, "Pseudo 3D wall friction may only be used in 2D");
170 return -factor * velocity * viscosity / (height*height);
171 }
172
174 template <class ElementVolumeVariables, class ElementFaceVariables, class G = GridGeometry>
175 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, Scalar>::type
176 pseudo3DWallFriction(const SubControlVolumeFace& scvf,
177 const ElementVolumeVariables& elemVolVars,
178 const ElementFaceVariables& elemFaceVars,
179 const Scalar height,
180 const Scalar factor = 8.0) const
181 {
182 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
183 const Scalar viscosity = elemVolVars[scvf.insideScvIdx()].effectiveViscosity();
184 return pseudo3DWallFriction(velocity, viscosity, height, factor);
185 }
186
192 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
193 {
194 DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the permeability must be returned in the acutal problem");
195 }
196
202 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
203 {
204 DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the alpha value must be returned in the acutal problem");
205 }
206
210 Scalar betaBJ(const Element& element, const SubControlVolumeFace& scvf, const GlobalPosition& tangentialVector) const
211 {
212 const Scalar interfacePermeability = interfacePermeability_(element, scvf, tangentialVector);
213 using std::sqrt;
214 return asImp_().alphaBJ(scvf) / sqrt(interfacePermeability);
215 }
216
220 VelocityVector porousMediumVelocity(const Element& element, const SubControlVolumeFace& scvf) const
221 {
222 return VelocityVector(0.0);
223 }
224
228 const Scalar beaversJosephVelocity(const Element& element,
229 const SubControlVolume& scv,
230 const SubControlVolumeFace& ownScvf,
231 const SubControlVolumeFace& faceOnPorousBoundary,
232 const Scalar velocitySelf,
233 const Scalar tangentialVelocityGradient) const
234 {
235 // create a unit normal vector oriented in positive coordinate direction
236 GlobalPosition orientation = ownScvf.unitOuterNormal();
237 orientation[ownScvf.directionIndex()] = 1.0;
238
239 // du/dy + dv/dx = alpha/sqrt(K) * (u_boundary-uPM)
240 // beta = alpha/sqrt(K)
241 const Scalar betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary, orientation);
242 const Scalar distanceNormalToBoundary = (faceOnPorousBoundary.center() - scv.center()).two_norm();
243
244 return (tangentialVelocityGradient*distanceNormalToBoundary
245 + asImp_().porousMediumVelocity(element, faceOnPorousBoundary) * orientation * betaBJ * distanceNormalToBoundary
246 + velocitySelf) / (betaBJ*distanceNormalToBoundary + 1.0);
247 }
248
249private:
251 Scalar interfacePermeability_(const Element& element, const SubControlVolumeFace& scvf, const GlobalPosition& tangentialVector) const
252 {
253 const auto& K = asImp_().permeability(element, scvf);
254
255 // use t*K*t for permeability tensors
256 if constexpr (Dune::IsNumber<std::decay_t<decltype(K)>>::value)
257 return K;
258 else
259 return vtmv(tangentialVector, K, tangentialVector);
260 }
261
263 Implementation &asImp_()
264 { return *static_cast<Implementation *>(this); }
265
267 const Implementation &asImp_() const
268 { return *static_cast<const Implementation *>(this); }
269
270 GravityVector gravity_;
271 bool enableInertiaTerms_;
272};
273
274} // end namespace Dumux
275
276#endif
Base class for all staggered fv problems.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:863
Definition: adapt.hh:29
typename NavierStokesParentProblemImpl< TypeTag, GetPropType< TypeTag, Properties::GridGeometry >::discMethod >::type NavierStokesParentProblem
The actual NavierStokesParentProblem.
Definition: freeflow/navierstokes/problem.hh:48
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
Base class for all staggered finite-volume problems.
Definition: staggeredfvproblem.hh:48
The implementation is specialized for the different discretizations.
Definition: freeflow/navierstokes/problem.hh:36
Navier-Stokes problem base class.
Definition: freeflow/navierstokes/problem.hh:60
Scalar temperature() const
Returns the temperature within the domain.
Definition: freeflow/navierstokes/problem.hh:123
Scalar pseudo3DWallFriction(const Scalar velocity, const Scalar viscosity, const Scalar height, const Scalar factor=8.0) const
An additional drag term can be included as source term for the momentum balance to mimic 3D flow beha...
Definition: freeflow/navierstokes/problem.hh:164
Scalar permeability(const Element &element, const SubControlVolumeFace &scvf) const
Returns the intrinsic permeability of required as input parameter for the Beavers-Joseph-Saffman boun...
Definition: freeflow/navierstokes/problem.hh:192
VelocityVector porousMediumVelocity(const Element &element, const SubControlVolumeFace &scvf) const
Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
Definition: freeflow/navierstokes/problem.hh:220
std::enable_if< G::discMethod==DiscretizationMethod::staggered, void >::type applyInitialFaceSolution(SolutionVector &sol, const SubControlVolumeFace &scvf, const PrimaryVariables &initSol) const
Applys the initial face solution (velocities on the faces). Specialization for staggered grid discret...
Definition: freeflow/navierstokes/problem.hh:144
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: freeflow/navierstokes/problem.hh:202
const GravityVector & gravity() const
Returns the acceleration due to gravity.
Definition: freeflow/navierstokes/problem.hh:132
std::enable_if< G::discMethod==DiscretizationMethod::staggered, Scalar >::type pseudo3DWallFriction(const SubControlVolumeFace &scvf, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const Scalar height, const Scalar factor=8.0) const
Convenience function for staggered grid implementation.
Definition: freeflow/navierstokes/problem.hh:176
Scalar betaBJ(const Element &element, const SubControlVolumeFace &scvf, const GlobalPosition &tangentialVector) const
Returns the beta value which is the alpha value divided by the square root of the (scalar-valued) int...
Definition: freeflow/navierstokes/problem.hh:210
const Scalar beaversJosephVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &ownScvf, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf, const Scalar tangentialVelocityGradient) const
Returns the slip velocity at a porous boundary based on the Beavers-Joseph(-Saffman) condition.
Definition: freeflow/navierstokes/problem.hh:228
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature at a given global position.
Definition: freeflow/navierstokes/problem.hh:115
bool enableInertiaTerms() const
Returns whether interia terms should be considered.
Definition: freeflow/navierstokes/problem.hh:138
NavierStokesProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor.
Definition: freeflow/navierstokes/problem.hh:97
Declares all properties used in Dumux.