3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
dumux/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/deprecated.hh>
28
29#include <dune/common/exceptions.hh>
33
34#include "model.hh"
35
36namespace Dumux {
37
39template<class TypeTag, DiscretizationMethod discMethod> struct NavierStokesParentProblemImpl;
40
41template<class TypeTag>
43{
45};
46
48template<class TypeTag>
50 typename NavierStokesParentProblemImpl<TypeTag,
52
61template<class TypeTag>
63{
64 using ParentType = NavierStokesParentProblem<TypeTag>;
65 using Implementation = GetPropType<TypeTag, Properties::Problem>;
66
68 using GridView = typename GridGeometry::GridView;
69 using Element = typename GridView::template Codim<0>::Entity;
70
72 using GridFaceVariables = typename GridVariables::GridFaceVariables;
73 using ElementFaceVariables = typename GridFaceVariables::LocalView;
74 using FaceVariables = typename GridFaceVariables::FaceVariables;
75 using GridVolumeVariables = typename GridVariables::GridVolumeVariables;
76 using ElementVolumeVariables = typename GridVolumeVariables::LocalView;
78
79 using FVElementGeometry = typename GridGeometry::LocalView;
80 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
81 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
84
85 enum {
86 dim = GridView::dimension,
87 dimWorld = GridView::dimensionworld
88 };
89
90 using GlobalPosition = typename SubControlVolumeFace::GlobalPosition;
91 using GravityVector = Dune::FieldVector<Scalar, dimWorld>;
92
93public:
99 NavierStokesProblem(std::shared_ptr<const GridGeometry> gridGeometry, const std::string& paramGroup = "")
100 : ParentType(gridGeometry, paramGroup)
101 , gravity_(0.0)
102 {
103 if (getParamFromGroup<bool>(paramGroup, "Problem.EnableGravity"))
104 gravity_[dim-1] = -9.81;
105
106 enableInertiaTerms_ = getParamFromGroup<bool>(paramGroup, "Problem.EnableInertiaTerms");
107 }
108
117 Scalar temperatureAtPos(const GlobalPosition &globalPos) const
118 { return asImp_().temperature(); }
119
125 Scalar temperature() const
126 { DUNE_THROW(Dune::NotImplemented, "temperature() method not implemented by the actual problem"); }
127
134 const GravityVector& gravity() const
135 { return gravity_; }
136
141 { return enableInertiaTerms_; }
142
144 template <class SolutionVector, class G = GridGeometry>
145 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, void>::type
146 applyInitialFaceSolution(SolutionVector& sol,
147 const SubControlVolumeFace& scvf,
148 const PrimaryVariables& initSol) const
149 {
150 sol[GridGeometry::faceIdx()][scvf.dofIndex()][0] = initSol[Indices::velocity(scvf.directionIndex())];
151 }
152
153
167 Scalar pseudo3DWallFriction(const Scalar velocity,
168 const Scalar viscosity,
169 const Scalar height,
170 const Scalar factor = 8.0) const
171 {
172 static_assert(dim == 2, "Pseudo 3D wall friction may only be used in 2D");
173 return -factor * velocity * viscosity / (height*height);
174 }
175
177 template <class ElementVolumeVariables, class ElementFaceVariables, class G = GridGeometry>
178 typename std::enable_if<G::discMethod == DiscretizationMethod::staggered, Scalar>::type
179 pseudo3DWallFriction(const SubControlVolumeFace& scvf,
180 const ElementVolumeVariables& elemVolVars,
181 const ElementFaceVariables& elemFaceVars,
182 const Scalar height,
183 const Scalar factor = 8.0) const
184 {
185 const Scalar velocity = elemFaceVars[scvf].velocitySelf();
186 const Scalar viscosity = elemVolVars[scvf.insideScvIdx()].effectiveViscosity();
187 return pseudo3DWallFriction(velocity, viscosity, height, factor);
188 }
189
195 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
196 {
197 DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the permeability must be returned in the acutal problem");
198 }
199
205 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
206 {
207 DUNE_THROW(Dune::NotImplemented, "When using the Beavers-Joseph-Saffman boundary condition, the alpha value must be returned in the acutal problem");
208 }
209
215 Scalar betaBJ(const Element& element, const SubControlVolumeFace& scvf) const
216 {
217 using std::sqrt;
218 return asImp_().alphaBJ(scvf) / sqrt(asImp_().permeability(element, scvf));
219 }
220
224 Scalar velocityPorousMedium(const Element& element, const SubControlVolumeFace& scvf) const
225 { return 0.0; }
226
228 DUNE_DEPRECATED_MSG("Use beaversJosephVelocity(element, scv, faceOnPorousBoundary, velocitySelf, tangentialVelocityGradient) instead")
229 const Scalar bjsVelocity(const Element& element,
230 const SubControlVolume& scv,
231 const SubControlVolumeFace& faceOnPorousBoundary,
232 const Scalar velocitySelf) const
233 {
234 // assume tangential velocity gradient of zero and
235 return beaversJosephVelocity(element, scv, faceOnPorousBoundary, velocitySelf, 0.0);
236 }
237
239 const Scalar beaversJosephVelocity(const Element& element,
240 const SubControlVolume& scv,
241 const SubControlVolumeFace& faceOnPorousBoundary,
242 const Scalar velocitySelf,
243 const Scalar tangentialVelocityGradient) const
244 {
245 // du/dy + dv/dx = alpha/sqrt(K) * (u_boundary-uPM)
246 // beta = alpha/sqrt(K)
247 const Scalar betaBJ = asImp_().betaBJ(element, faceOnPorousBoundary);
248 const Scalar distance = (faceOnPorousBoundary.center() - scv.center()).two_norm();
249 return (tangentialVelocityGradient*distance + asImp_().velocityPorousMedium(element,faceOnPorousBoundary)*betaBJ*distance + velocitySelf) / (betaBJ*distance + 1.0);
250 }
251
252private:
253
255 Implementation &asImp_()
256 { return *static_cast<Implementation *>(this); }
257
259 const Implementation &asImp_() const
260 { return *static_cast<const Implementation *>(this); }
261
262 GravityVector gravity_;
263 bool enableInertiaTerms_;
264};
265
266} // end namespace Dumux
267
268#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
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename NavierStokesParentProblemImpl< TypeTag, GetPropType< TypeTag, Properties::GridGeometry >::discMethod >::type NavierStokesParentProblem
The actual NavierStokesParentProblem.
Definition: dumux/freeflow/navierstokes/problem.hh:51
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
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:45
The implementation is specialized for the different discretizations.
Definition: dumux/freeflow/navierstokes/problem.hh:39
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
Scalar temperature() const
Returns the temperature within the domain.
Definition: dumux/freeflow/navierstokes/problem.hh:125
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: dumux/freeflow/navierstokes/problem.hh:167
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: dumux/freeflow/navierstokes/problem.hh:195
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: dumux/freeflow/navierstokes/problem.hh:146
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: dumux/freeflow/navierstokes/problem.hh:205
const GravityVector & gravity() const
Returns the acceleration due to gravity.
Definition: dumux/freeflow/navierstokes/problem.hh:134
Scalar velocityPorousMedium(const Element &element, const SubControlVolumeFace &scvf) const
Returns the velocity in the porous medium (which is 0 by default according to Saffmann).
Definition: dumux/freeflow/navierstokes/problem.hh:224
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: dumux/freeflow/navierstokes/problem.hh:179
Scalar betaBJ(const Element &element, const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: dumux/freeflow/navierstokes/problem.hh:215
const Scalar beaversJosephVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf, const Scalar tangentialVelocityGradient) const
helper function to evaluate the slip velocity on the boundary when the Beavers-Joseph condition is us...
Definition: dumux/freeflow/navierstokes/problem.hh:239
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature at a given global position.
Definition: dumux/freeflow/navierstokes/problem.hh:117
bool enableInertiaTerms() const
Returns whether interia terms should be considered.
Definition: dumux/freeflow/navierstokes/problem.hh:140
const Scalar bjsVelocity(const Element &element, const SubControlVolume &scv, const SubControlVolumeFace &faceOnPorousBoundary, const Scalar velocitySelf) const
helper function to evaluate the slip velocity on the boundary when the Beavers-Joseph-Saffman conditi...
Definition: dumux/freeflow/navierstokes/problem.hh:229
NavierStokesProblem(std::shared_ptr< const GridGeometry > gridGeometry, const std::string &paramGroup="")
The constructor.
Definition: dumux/freeflow/navierstokes/problem.hh:99
Declares all properties used in Dumux.
Adaption of the fully implicit scheme to the tracer transport model.