3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
1p2c_1p2c/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 StokesOnePTwoC { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
45} // end namespace TTag
46
47// The fluid system
48template<class TypeTag>
49struct FluidSystem<TypeTag, TTag::StokesOnePTwoC>
50{
52 static constexpr auto phaseIdx = H2OAir::liquidPhaseIdx; // simulate the water phase
54};
55
56// Set the grid type
57template<class TypeTag>
58struct Grid<TypeTag, TTag::StokesOnePTwoC> { using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
59
60// Set the problem property
61template<class TypeTag>
62struct Problem<TypeTag, TTag::StokesOnePTwoC> { using type = Dumux::StokesSubProblem<TypeTag> ; };
63
64template<class TypeTag>
65struct EnableGridGeometryCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
66template<class TypeTag>
67struct EnableGridFluxVariablesCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
68template<class TypeTag>
69struct EnableGridVolumeVariablesCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
70
71// Use moles
72template<class TypeTag>
73struct UseMoles<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
74
75// Do not replace one equation with a total mass balance
76template<class TypeTag>
77struct ReplaceCompEqIdx<TypeTag, TTag::StokesOnePTwoC> { static constexpr int value = 3; };
78} // end namespace Properties
79
86template <class TypeTag>
87class StokesSubProblem : public NavierStokesProblem<TypeTag>
88{
89 using ParentType = NavierStokesProblem<TypeTag>;
90 using GridView = GetPropType<TypeTag, Properties::GridView>;
91 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
93 using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
94 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
95 using FVElementGeometry = typename GridGeometry::LocalView;
96 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
97 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
98 using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
99
100 using Element = typename GridView::template Codim<0>::Entity;
101 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
102
103 using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>;
104 using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;
105
106public:
107 StokesSubProblem(std::shared_ptr<const GridGeometry> gridGeometry, std::shared_ptr<CouplingManager> couplingManager)
108 : ParentType(gridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
109 {
110 problemName_ = getParam<std::string>("Vtk.OutputName") + "_" + getParamFromGroup<std::string>(this->paramGroup(), "Problem.Name");
111
112 // determine whether to simulate a vertical or horizontal flow configuration
113 verticalFlow_ = problemName_.find("vertical") != std::string::npos;
114 }
115
119 const std::string& name() const
120 {
121 return problemName_;
122 }
123
127 // \{
128
134 Scalar temperature() const
135 { return 273.15 + 10; } // 10°C
136
142 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
143 { return NumEqVector(0.0); }
144 // \}
145
149 // \{
150
158 BoundaryTypes boundaryTypes(const Element& element,
159 const SubControlVolumeFace& scvf) const
160 {
161 BoundaryTypes values;
162
163 const auto& globalPos = scvf.dofPosition();
164
165 if (verticalFlow_)
166 {
167 // inflow
168 if(onUpperBoundary_(globalPos))
169 {
170 values.setDirichlet(Indices::velocityXIdx);
171 values.setDirichlet(Indices::velocityYIdx);
172 values.setDirichlet(Indices::conti0EqIdx + 1);
173 }
174
175 // left/right wall
176 if (onRightBoundary_(globalPos) || (onLeftBoundary_(globalPos)))
177 {
178 values.setDirichlet(Indices::velocityXIdx);
179 values.setDirichlet(Indices::velocityYIdx);
180 values.setNeumann(Indices::conti0EqIdx);
181 values.setNeumann(Indices::conti0EqIdx + 1);
182 }
183
184 }
185 else // horizontal flow
186 {
187 if (onLeftBoundary_(globalPos))
188 {
189 values.setDirichlet(Indices::conti0EqIdx + 1);
190 values.setDirichlet(Indices::velocityXIdx);
191 values.setDirichlet(Indices::velocityYIdx);
192 }
193 else if (onRightBoundary_(globalPos))
194 {
195 values.setDirichlet(Indices::pressureIdx);
196 values.setOutflow(Indices::conti0EqIdx + 1);
197 }
198 else
199 {
200 values.setDirichlet(Indices::velocityXIdx);
201 values.setDirichlet(Indices::velocityYIdx);
202 values.setNeumann(Indices::conti0EqIdx);
203 values.setNeumann(Indices::conti0EqIdx + 1);
204 }
205 }
206
207 if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
208 {
209 values.setCouplingNeumann(Indices::conti0EqIdx);
210 values.setCouplingNeumann(Indices::conti0EqIdx + 1);
211 values.setCouplingNeumann(Indices::momentumYBalanceIdx);
212 values.setBJS(Indices::momentumXBalanceIdx);
213 }
214
215 return values;
216 }
217
221 PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
222 {
223 PrimaryVariables values(0.0);
224 values = initialAtPos(globalPos);
225
226 if (verticalFlow_)
227 {
228 // Check if this a pure diffusion problem.
229 static const bool isDiffusionProblem = problemName_.find("diffusion") != std::string::npos;
230
231 Scalar topMoleFraction = 1e-3;
232
233 if (isDiffusionProblem)
234 {
235 // For the diffusion problem, change the top mole fraction after some time
236 // in order to revert the concentration gradient.
237 if (time() >= 1e10)
238 topMoleFraction = 0.0;
239 }
240 else // advection problem
241 {
242 // reverse the flow direction after some time for the advection problem
243 if (time() >= 3e5)
244 values[Indices::velocityYIdx] *= -1.0;
245 }
246
247 if(globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_)
248 values[Indices::conti0EqIdx + 1] = topMoleFraction;
249 }
250 else // horizontal flow
251 {
252 static const Scalar inletMoleFraction = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.InletMoleFraction");
253 if(globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_)
254 values[Indices::conti0EqIdx + 1] = inletMoleFraction;
255 }
256
257
258 return values;
259 }
260
270 template<class ElementVolumeVariables, class ElementFaceVariables>
271 NumEqVector neumann(const Element& element,
272 const FVElementGeometry& fvGeometry,
273 const ElementVolumeVariables& elemVolVars,
274 const ElementFaceVariables& elemFaceVars,
275 const SubControlVolumeFace& scvf) const
276 {
277 NumEqVector values(0.0);
278
279 if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
280 {
281 values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
282
283 const auto tmp = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
284 values[Indices::conti0EqIdx] = tmp[0];
285 values[Indices::conti0EqIdx + 1] = tmp[1];
286 }
287 return values;
288 }
289
290 // \}
291
293 const CouplingManager& couplingManager() const
294 { return *couplingManager_; }
295
299 // \{
300
306 PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const
307 {
308 PrimaryVariables values(0.0);
309 values[Indices::pressureIdx] = 1e5;
310
311 static const Scalar vMax = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.Velocity", 0.0);
312
313 auto parabolicProfile = [&](const GlobalPosition& globalPos, int coord)
314 {
315 return vMax * (globalPos[coord] - this->gridGeometry().bBoxMin()[coord])
316 * (this->gridGeometry().bBoxMax()[coord] - globalPos[coord])
317 / (0.25 * (this->gridGeometry().bBoxMax()[coord] - this->gridGeometry().bBoxMin()[coord])
318 * (this->gridGeometry().bBoxMax()[coord] - this->gridGeometry().bBoxMin()[coord]));
319 };
320
321 if (verticalFlow_)
322 values[Indices::velocityYIdx] = parabolicProfile(globalPos, 0);
323 else // horizontal flow
324 values[Indices::velocityXIdx] = parabolicProfile(globalPos, 1);
325
326 return values;
327 }
328
333 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
334 {
335 return couplingManager().couplingData().darcyPermeability(element, scvf);
336 }
337
342 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
343 {
344 return 1.0;
345 }
346
350 void setTimeLoop(TimeLoopPtr timeLoop)
351 { timeLoop_ = timeLoop; }
352
356 Scalar time() const
357 { return timeLoop_->time(); }
358
359 // \}
360
361private:
362 bool onLeftBoundary_(const GlobalPosition &globalPos) const
363 { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
364
365 bool onRightBoundary_(const GlobalPosition &globalPos) const
366 { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
367
368 bool onLowerBoundary_(const GlobalPosition &globalPos) const
369 { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
370
371 bool onUpperBoundary_(const GlobalPosition &globalPos) const
372 { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
373
374 Scalar eps_;
375 bool verticalFlow_;
376 std::string problemName_;
377 std::shared_ptr<CouplingManager> couplingManager_;
378 TimeLoopPtr timeLoop_;
379};
380} // end namespace Dumux
381
382#endif // DUMUX_STOKES_SUBPROBLEM_HH
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
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
The type of the fluid system to use.
Definition: common/properties.hh:223
An adapter for multi-phase fluid systems to be used with (compositional) one-phase models.
Definition: 1padapter.hh:46
A compositional two-phase fluid system with water and air as components in both, the liquid and the g...
Definition: h2oair.hh:75
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: 1p2c_1p2c/problem_stokes.hh:134
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: 1p2c_1p2c/problem_stokes.hh:271
void setTimeLoop(TimeLoopPtr timeLoop)
Sets the time loop pointer.
Definition: 1p2c_1p2c/problem_stokes.hh:350
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: 1p2c_1p2c/problem_stokes.hh:342
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: 1p2c_1p2c/problem_stokes.hh:158
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: 1p2c_1p2c/problem_stokes.hh:107
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: 1p2c_1p2c/problem_stokes.hh:142
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: 1p2c_1p2c/problem_stokes.hh:221
Scalar time() const
Returns the time.
Definition: 1p2c_1p2c/problem_stokes.hh:356
const std::string & name() const
The problem name.
Definition: 1p2c_1p2c/problem_stokes.hh:119
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: 1p2c_1p2c/problem_stokes.hh:333
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:267
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:57
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:57
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: 1p2c_1p2c/diffusionlawcomparison/problem_stokes.hh:71
Defines a type tag and some properties for ree-flow models using the staggered scheme....