3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
1p3c_1p3c/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_ONEPTHREEC_HH
26#define DUMUX_STOKES_SUBPROBLEM_ONEPTHREEC_HH
27
28#include <dune/grid/yaspgrid.hh>
29
30#include "h2n2co2fluidsystem.hh"
32
36
37// for StokesDarcyCouplingOptions
39
40namespace Dumux {
41template <class TypeTag>
42class StokesSubProblem;
43
44namespace Properties {
45// Create new type tags
46namespace TTag {
47struct StokesOnePThreeC { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
48} // end namespace TTag
49
50// Set the fluid system
51template<class TypeTag>
53
54// Set the grid type
55template<class TypeTag>
56struct Grid<TypeTag, TTag::StokesOnePThreeC> { using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
57
58// Set the problem property
59template<class TypeTag>
60struct Problem<TypeTag, TTag::StokesOnePThreeC> { using type = Dumux::StokesSubProblem<TypeTag> ; };
61
62template<class TypeTag>
63struct EnableGridGeometryCache<TypeTag, TTag::StokesOnePThreeC> { static constexpr bool value = true; };
64template<class TypeTag>
65struct EnableGridFluxVariablesCache<TypeTag, TTag::StokesOnePThreeC> { static constexpr bool value = true; };
66template<class TypeTag>
67struct EnableGridVolumeVariablesCache<TypeTag, TTag::StokesOnePThreeC> { static constexpr bool value = true; };
68
69// Use moles
70template<class TypeTag>
71struct UseMoles<TypeTag, TTag::StokesOnePThreeC> { static constexpr bool value = true; };
72
73// Set the grid type
74template<class TypeTag>
75struct MolecularDiffusionType<TypeTag, TTag::StokesOnePThreeC> { using type = MaxwellStefansLaw<TypeTag>; };
76
77// Do not replace one equation with a total mass balance
78template<class TypeTag>
79struct ReplaceCompEqIdx<TypeTag, TTag::StokesOnePThreeC> { static constexpr int value = 3; };
80} // end namespace Properties
81
88template <class TypeTag>
89class StokesSubProblem : public NavierStokesProblem<TypeTag>
90{
91 using ParentType = NavierStokesProblem<TypeTag>;
98 using FVElementGeometry = typename GridGeometry::LocalView;
99 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
102 using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;
103
104 using Element = typename GridView::template Codim<0>::Entity;
105 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
106
108
109public:
110 StokesSubProblem(std::shared_ptr<const GridGeometry> gridGeometry, std::shared_ptr<CouplingManager> couplingManager)
111 : ParentType(gridGeometry, "Stokes"), eps_(1e-6),
112 couplingManager_(couplingManager)
113 {
114 inletVelocity_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Velocity");
115 pressure_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Pressure");
116 }
117
121 // \{
122
124 { return false; }
125
131 Scalar temperature() const
132 { return 273.15 + 10; } // 10°C
133
139 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
140 { return NumEqVector(0.0); }
141 // \}
142
146 // \{
147
155 BoundaryTypes boundaryTypes(const Element& element,
156 const SubControlVolumeFace& scvf) const
157 {
158 BoundaryTypes values;
159
160 const auto& globalPos = scvf.dofPosition();
161
162 if(onLeftBoundary_(globalPos))
163 {
164 values.setDirichlet(Indices::conti0EqIdx + 2);
165 values.setDirichlet(Indices::conti0EqIdx + 1);
166 values.setDirichlet(Indices::velocityXIdx);
167 values.setDirichlet(Indices::velocityYIdx);
168 }
169
170 else if(onRightBoundary_(globalPos))
171 {
172 values.setDirichlet(Indices::pressureIdx);
173 values.setOutflow(Indices::conti0EqIdx + 1);
174 values.setOutflow(Indices::conti0EqIdx + 2);
175
176 }
177 else
178 {
179 values.setDirichlet(Indices::velocityXIdx);
180 values.setDirichlet(Indices::velocityYIdx);
181 values.setNeumann(Indices::conti0EqIdx);
182 values.setNeumann(Indices::conti0EqIdx + 1);
183 values.setNeumann(Indices::conti0EqIdx + 2);
184 }
185
186 if(couplingManager().isCoupledEntity(couplingManager().stokesIdx, scvf))
187 {
188 values.setNeumann(Indices::conti0EqIdx);
189 values.setNeumann(Indices::conti0EqIdx+1);
190 values.setNeumann(Indices::conti0EqIdx+2);
191 values.setNeumann(Indices::momentumYBalanceIdx);
192 values.setBJS(Indices::velocityXIdx);
193 }
194
195 return values;
196 }
197
201 PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
202 {
203 PrimaryVariables values(0.0);
204 values = initialAtPos(globalPos);
205 return values;
206 }
207
217 template<class ElementVolumeVariables, class ElementFaceVariables>
218 NumEqVector neumann(const Element& element,
219 const FVElementGeometry& fvGeometry,
220 const ElementVolumeVariables& elemVolVars,
221 const ElementFaceVariables& elemFaceVars,
222 const SubControlVolumeFace& scvf) const
223 {
224 NumEqVector values(0.0);
225
226 if(couplingManager().isCoupledEntity(couplingManager().stokesIdx, scvf))
227 {
228 values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
229
230 const auto tmp = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, DiffusionCoefficientAveragingType::harmonic);
231 values[Indices::conti0EqIdx] = tmp[0];
232 values[Indices::conti0EqIdx + 1] = tmp[1];
233 values[Indices::conti0EqIdx + 2] = tmp[2];
234 }
235 return values;
236 }
237
238 // \}
239
241 const CouplingManager& couplingManager() const
242 { return *couplingManager_; }
243
247 // \{
248
254 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
255 {
256 PrimaryVariables values(0.0);
257 values[Indices::pressureIdx] = pressure_;
258 values[Indices::velocityXIdx] = inletVelocity_ * (globalPos[1] - this->gridGeometry().bBoxMin()[1])
259 * (this->gridGeometry().bBoxMax()[1] - globalPos[1])
260 / (0.25 * (this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1])
261 * (this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]));
262 return values;
263 }
264
269 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
270 {
271 return couplingManager().couplingData().darcyPermeability(element, scvf);
272 }
273
278 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
279 {
280 return 1.0;
281 }
282
283 // \}
284
285private:
286 bool onLeftBoundary_(const GlobalPosition &globalPos) const
287 { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
288
289 bool onRightBoundary_(const GlobalPosition &globalPos) const
290 { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
291
292 bool onLowerBoundary_(const GlobalPosition &globalPos) const
293 { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
294
295 bool onUpperBoundary_(const GlobalPosition &globalPos) const
296 { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
297
298 Scalar eps_;
299 Scalar inletVelocity_;
300 Scalar pressure_;
301 std::shared_ptr<CouplingManager> couplingManager_;
302};
303} // end namespace Dumux
304
305#endif // DUMUX_STOKES_SUBPROBLEM_HH
This file contains the data which is required to calculate diffusive mass fluxes due to molecular dif...
A fluid system for one phase with the components h2, n2 and co2.
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
Property whether to use moles or kg as amount unit for balance equations.
Definition: common/properties.hh:102
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:104
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 for the calculation of the molecular diffusion fluxes.
Definition: common/properties.hh:212
The type of the fluid system to use.
Definition: common/properties.hh:223
Definition: maxwellstefanslaw.hh:37
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficiencients (moleculat diffusion or thermal conduc...
Definition: couplingdata.hh:52
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: 1p3c_1p3c/problem_stokes.hh:131
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: 1p3c_1p3c/problem_stokes.hh:218
bool shouldWriteRestartFile() const
Definition: 1p3c_1p3c/problem_stokes.hh:123
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: 1p3c_1p3c/problem_stokes.hh:278
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: 1p3c_1p3c/problem_stokes.hh:155
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: 1p3c_1p3c/problem_stokes.hh:110
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: 1p3c_1p3c/problem_stokes.hh:139
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: 1p3c_1p3c/problem_stokes.hh:201
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: 1p3c_1p3c/problem_stokes.hh:269
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:267
A simple fluid system with one Maxwell-Stefan component.
Definition: h2n2co2fluidsystem.hh:39
Definition: 1p3c_1p3c/problem_stokes.hh:47
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: 1p3c_1p3c/problem_stokes.hh:47
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: 1p3c_1p3c/problem_stokes.hh:56
Defines a type tag and some properties for ree-flow models using the staggered scheme....