3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
1p_1p/problem_stokes.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 *****************************************************************************/
25#ifndef DUMUX_STOKES_SUBPROBLEM_HH
26#define DUMUX_STOKES_SUBPROBLEM_HH
27
28#include <dune/grid/yaspgrid.hh>
29
32
36
37namespace Dumux {
38template <class TypeTag>
39class StokesSubProblem;
40
41namespace Properties {
42// Create new type tags
43namespace TTag {
44struct StokesOneP { using InheritsFrom = std::tuple<NavierStokes, StaggeredFreeFlowModel>; };
45} // end namespace TTag
46
47// the fluid system
48template<class TypeTag>
49struct FluidSystem<TypeTag, TTag::StokesOneP>
50{
53};
54
55// Set the grid type
56template<class TypeTag>
57struct Grid<TypeTag, TTag::StokesOneP> { using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
58
59// Set the problem property
60template<class TypeTag>
61struct Problem<TypeTag, TTag::StokesOneP> { using type = Dumux::StokesSubProblem<TypeTag> ; };
62
63template<class TypeTag>
64struct EnableGridGeometryCache<TypeTag, TTag::StokesOneP> { static constexpr bool value = true; };
65template<class TypeTag>
66struct EnableGridFluxVariablesCache<TypeTag, TTag::StokesOneP> { static constexpr bool value = true; };
67template<class TypeTag>
68struct EnableGridVolumeVariablesCache<TypeTag, TTag::StokesOneP> { static constexpr bool value = true; };
69} // end namespace Properties
70
77template <class TypeTag>
78class StokesSubProblem : public NavierStokesProblem<TypeTag>
79{
80 using ParentType = NavierStokesProblem<TypeTag>;
81
84
86
88
90 using FVElementGeometry = typename GridGeometry::LocalView;
91 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
92 using Element = typename GridView::template Codim<0>::Entity;
93
94 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
95
98
100
101public:
102 StokesSubProblem(std::shared_ptr<const GridGeometry> gridGeometry, std::shared_ptr<CouplingManager> couplingManager)
103 : ParentType(gridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
104 {
105 problemName_ = getParam<std::string>("Vtk.OutputName") + "_" + getParamFromGroup<std::string>(this->paramGroup(), "Problem.Name");
106
107 // determine whether to simulate a vertical or horizontal flow configuration
108 verticalFlow_ = problemName_.find("vertical") != std::string::npos;
109 }
110
114 const std::string& name() const
115 {
116 return problemName_;
117 }
118
122 // \{
123
129 Scalar temperature() const
130 { return 273.15 + 10; } // 10°C
131
137 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
138 { return NumEqVector(0.0); }
139 // \}
140
144 // \{
145
153 BoundaryTypes boundaryTypes(const Element& element,
154 const SubControlVolumeFace& scvf) const
155 {
156 BoundaryTypes values;
157
158 const auto& globalPos = scvf.dofPosition();
159
160 if (verticalFlow_)
161 {
162 // inflow
163 if(onUpperBoundary_(globalPos))
164 {
165 values.setDirichlet(Indices::velocityXIdx);
166 values.setDirichlet(Indices::velocityYIdx);
167 }
168
169 // left/right wall
170 if (onRightBoundary_(globalPos) || (onLeftBoundary_(globalPos)))
171 {
172 values.setDirichlet(Indices::velocityXIdx);
173 values.setDirichlet(Indices::velocityYIdx);
174 }
175 }
176 else // horizontal flow
177 {
178 // inlet / outlet
179 if(onLeftBoundary_(globalPos) || onRightBoundary_(globalPos))
180 values.setDirichlet(Indices::pressureIdx);
181 else // wall
182 {
183 values.setDirichlet(Indices::velocityXIdx);
184 values.setDirichlet(Indices::velocityYIdx);
185 }
186
187 }
188
189 if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
190 {
191 values.setCouplingNeumann(Indices::conti0EqIdx);
192 values.setCouplingNeumann(Indices::momentumYBalanceIdx);
193 values.setBJS(Indices::momentumXBalanceIdx);
194 }
195
196 return values;
197 }
198
202 PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
203 {
204 return initialAtPos(globalPos);
205 }
206
216 template<class ElementVolumeVariables, class ElementFaceVariables>
217 NumEqVector neumann(const Element& element,
218 const FVElementGeometry& fvGeometry,
219 const ElementVolumeVariables& elemVolVars,
220 const ElementFaceVariables& elemFaceVars,
221 const SubControlVolumeFace& scvf) const
222 {
223 NumEqVector values(0.0);
224
225 if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
226 {
227 values[Indices::conti0EqIdx] = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
228 values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
229 }
230 return values;
231 }
232
233 // \}
234
236 const CouplingManager& couplingManager() const
237 { return *couplingManager_; }
238
242 // \{
243
249 PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const
250 {
251 PrimaryVariables values(0.0);
252
253 if (verticalFlow_)
254 {
255 static const Scalar inletVelocity = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Velocity");
256 values[Indices::velocityYIdx] = inletVelocity * globalPos[0] * (this->gridGeometry().bBoxMax()[0] - globalPos[0]);
257 }
258 else // horizontal flow
259 {
260 static const Scalar deltaP = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.PressureDifference");
261 if(onLeftBoundary_(globalPos))
262 values[Indices::pressureIdx] = deltaP;
263 }
264
265 return values;
266 }
267
272 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
273 {
274 return couplingManager().couplingData().darcyPermeability(element, scvf);
275 }
276
281 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
282 {
283 return couplingManager().problem(CouplingManager::darcyIdx).spatialParams().beaversJosephCoeffAtPos(scvf.center());
284 }
285
286 // \}
287
288private:
289 bool onLeftBoundary_(const GlobalPosition &globalPos) const
290 { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
291
292 bool onRightBoundary_(const GlobalPosition &globalPos) const
293 { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
294
295 bool onLowerBoundary_(const GlobalPosition &globalPos) const
296 { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
297
298 bool onUpperBoundary_(const GlobalPosition &globalPos) const
299 { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
300
301 Scalar eps_;
302 std::string problemName_;
303 bool verticalFlow_;
304
305 std::shared_ptr<CouplingManager> couplingManager_;
306};
307} // end namespace Dumux
308
309#endif // DUMUX_STOKES_SUBPROBLEM_HH
A much simpler (and thus potentially less buggy) version of pure water.
A liquid phase consisting of a single component.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
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
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
The DUNE grid type.
Definition: common/properties.hh:57
Property to specify the type of a problem which has to be solved.
Definition: common/properties.hh:69
Definition: common/properties.hh:169
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
The type of the fluid system to use.
Definition: common/properties.hh:223
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
A liquid phase consisting of a single component.
Definition: 1pliquid.hh:46
Definition: multidomain/couplingmanager.hh:46
Test problem for the 1pnc (Navier-) Stokes problem.
Definition: 1p_2p/problem_stokes.hh:68
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: 1p_1p/problem_stokes.hh:129
NumEqVector neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Evaluates the boundary conditions for a Neumann control volume.
Definition: 1p_1p/problem_stokes.hh:217
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: 1p_1p/problem_stokes.hh:281
BoundaryTypes boundaryTypes(const Element &element, const SubControlVolumeFace &scvf) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: 1p_1p/problem_stokes.hh:153
const CouplingManager & couplingManager() const
Get the coupling manager.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:254
StokesSubProblem(std::shared_ptr< const GridGeometry > gridGeometry, std::shared_ptr< CouplingManager > couplingManager)
Definition: 1p_1p/problem_stokes.hh:102
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: 1p_1p/problem_stokes.hh:137
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: 1p_1p/problem_stokes.hh:202
const std::string & name() const
The problem name.
Definition: 1p_1p/problem_stokes.hh:114
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: 1p_1p/problem_stokes.hh:272
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:267
Definition: 1p_1p/problem_stokes.hh:44
std::tuple< NavierStokes, StaggeredFreeFlowModel > InheritsFrom
Definition: 1p_1p/problem_stokes.hh:44
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: 1p_1p/problem_stokes.hh:51
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: 1p_1p/problem_stokes.hh:57
Defines a type tag and some properties for ree-flow models using the staggered scheme....
A single-phase, isothermal Navier-Stokes model.