3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
1p2c_2p2c/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_STOKES1P2C_SUBPROBLEM_HH
26#define DUMUX_STOKES1P2C_SUBPROBLEM_HH
27
28#include <dune/grid/yaspgrid.hh>
29
32
37
38namespace Dumux {
39template <class TypeTag>
40class StokesSubProblem;
41
42namespace Properties {
43// Create new type tags
44namespace TTag {
45#if !NONISOTHERMAL
46struct StokesOnePTwoC { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
47#else
48struct StokesOnePTwoC { using InheritsFrom = std::tuple<NavierStokesNCNI, StaggeredFreeFlowModel>; };
49#endif
50} // end namespace TTag
51
52
53// Set the grid type
54template<class TypeTag>
55struct Grid<TypeTag, TTag::StokesOnePTwoC> { using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
56
57// The fluid system
58template<class TypeTag>
59struct FluidSystem<TypeTag, TTag::StokesOnePTwoC>
60{
62 static constexpr auto phaseIdx = H2OAir::gasPhaseIdx; // simulate the water phase
64};
65
66template<class TypeTag>
67struct ReplaceCompEqIdx<TypeTag, TTag::StokesOnePTwoC> { static constexpr int value = 3; };
68
69// Use formulation based on mass fractions
70template<class TypeTag>
71struct UseMoles<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
72
73// Set the problem property
74template<class TypeTag>
75struct Problem<TypeTag, TTag::StokesOnePTwoC> { using type = Dumux::StokesSubProblem<TypeTag> ; };
76
77template<class TypeTag>
78struct EnableGridGeometryCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
79template<class TypeTag>
80struct EnableGridFluxVariablesCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
81template<class TypeTag>
82struct EnableGridVolumeVariablesCache<TypeTag, TTag::StokesOnePTwoC> { static constexpr bool value = true; };
83} // end namespace Properties
84
91template <class TypeTag>
92class StokesSubProblem : public NavierStokesProblem<TypeTag>
93{
94 using ParentType = NavierStokesProblem<TypeTag>;
95
96 using GridView = GetPropType<TypeTag, Properties::GridView>;
97 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
98 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
100 using BoundaryTypes = GetPropType<TypeTag, Properties::BoundaryTypes>;
101
102 using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>;
103 using FVElementGeometry = typename GridGeometry::LocalView;
104 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
105 using Element = typename GridView::template Codim<0>::Entity;
106 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
107 using ElementFaceVariables = typename GetPropType<TypeTag, Properties::GridFaceVariables>::LocalView;
108 using FluidState = GetPropType<TypeTag, Properties::FluidState>;
109
110 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
111
112 using PrimaryVariables = GetPropType<TypeTag, Properties::PrimaryVariables>;
113 using NumEqVector = GetPropType<TypeTag, Properties::NumEqVector>;
114
115 using CouplingManager = GetPropType<TypeTag, Properties::CouplingManager>;
116 using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;
117
118 using DiffusionCoefficientAveragingType = typename StokesDarcyCouplingOptions::DiffusionCoefficientAveragingType;
119
120 static constexpr bool useMoles = GetPropType<TypeTag, Properties::ModelTraits>::useMoles();
121
122public:
123 StokesSubProblem(std::shared_ptr<const GridGeometry> gridGeometry, std::shared_ptr<CouplingManager> couplingManager)
124 : ParentType(gridGeometry, "Stokes"), eps_(1e-6), couplingManager_(couplingManager)
125 {
126 refVelocity_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefVelocity");
127 refPressure_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefPressure");
128 refMoleFrac_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.refMoleFrac");
129 refTemperature_ = getParamFromGroup<Scalar>(this->paramGroup(), "Problem.RefTemperature");
130
131 diffCoeffAvgType_ = StokesDarcyCouplingOptions::stringToEnum(DiffusionCoefficientAveragingType{},
132 getParamFromGroup<std::string>(this->paramGroup(), "Problem.InterfaceDiffusionCoefficientAvg"));
133 problemName_ = getParam<std::string>("Vtk.OutputName") + "_" + getParamFromGroup<std::string>(this->paramGroup(), "Problem.Name");
134 }
135
139 const std::string& name() const
140 {
141 return problemName_;
142 }
143
147 // \{
148
152 Scalar temperature() const
153 { return refTemperature_; }
154
160 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
161 { return NumEqVector(0.0); }
162
163 // \}
167 // \{
168
176 BoundaryTypes boundaryTypes(const Element& element,
177 const SubControlVolumeFace& scvf) const
178 {
179 BoundaryTypes values;
180
181 const auto& globalPos = scvf.center();
182
183#if NONISOTHERMAL
184 values.setNeumann(Indices::energyEqIdx);
185#endif
186
187 if (onUpperBoundary_(globalPos) || onLeftBoundary_(globalPos))
188 {
189 values.setDirichlet(Indices::velocityXIdx);
190 values.setDirichlet(Indices::velocityYIdx);
191 values.setNeumann(Indices::conti0EqIdx);
192 values.setNeumann(Indices::conti0EqIdx + 1);
193 }
194
195 if (onRightBoundary_(globalPos))
196 {
197 values.setDirichlet(Indices::pressureIdx);
198 values.setOutflow(Indices::conti0EqIdx + 1);
199
200#if NONISOTHERMAL
201 values.setOutflow(Indices::energyEqIdx);
202#endif
203 }
204
205 if (couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
206 {
207 values.setCouplingNeumann(Indices::conti0EqIdx);
208 values.setCouplingNeumann(Indices::conti0EqIdx + 1);
209 values.setCouplingNeumann(Indices::momentumYBalanceIdx);
210 values.setBJS(Indices::momentumXBalanceIdx);
211 }
212 return values;
213 }
214
218 PrimaryVariables dirichletAtPos(const GlobalPosition& pos) const
219 {
220 PrimaryVariables values(0.0);
221 values = initialAtPos(pos);
222
223 return values;
224 }
225
235 NumEqVector neumann(const Element& element,
236 const FVElementGeometry& fvGeometry,
237 const ElementVolumeVariables& elemVolVars,
238 const ElementFaceVariables& elemFaceVars,
239 const SubControlVolumeFace& scvf) const
240 {
241 PrimaryVariables values(0.0);
242 const auto& globalPos = scvf.dofPosition();
243 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
244
245 FluidState fluidState;
246 updateFluidStateForBC_(fluidState, elemVolVars[scv].pressure());
247
248 const Scalar density = useMoles ? fluidState.molarDensity(0) : fluidState.density(0);
249 const Scalar xVelocity = xVelocity_(globalPos);
250
251 if (onLeftBoundary_(globalPos))
252 {
253 // rho*v*X at inflow
254 values[Indices::conti0EqIdx + 1] = -xVelocity * density * refMoleFrac();
255 values[Indices::conti0EqIdx] = -xVelocity * density * (1.0 - refMoleFrac());
256
257#if NONISOTHERMAL
258 values[Indices::energyEqIdx] = -xVelocity * fluidState.density(0) * fluidState.enthalpy(0);
259#endif
260 }
261
262 if(couplingManager().isCoupledEntity(CouplingManager::stokesIdx, scvf))
263 {
264 values[Indices::momentumYBalanceIdx] = couplingManager().couplingData().momentumCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf);
265
266 const auto massFlux = couplingManager().couplingData().massCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
267 values[Indices::conti0EqIdx] = massFlux[0];
268 values[Indices::conti0EqIdx + 1] = massFlux[1];
269
270#if NONISOTHERMAL
271 values[Indices::energyEqIdx] = couplingManager().couplingData().energyCouplingCondition(element, fvGeometry, elemVolVars, elemFaceVars, scvf, diffCoeffAvgType_);
272#endif
273
274 }
275 return values;
276 }
277
278 // \}
279
281 const CouplingManager& couplingManager() const
282 { return *couplingManager_; }
283
287 // \{
288
294 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
295 {
296 FluidState fluidState;
297 updateFluidStateForBC_(fluidState, refPressure());
298
299 const Scalar density = FluidSystem::density(fluidState, 0);
300
301 PrimaryVariables values(0.0);
302 values[Indices::pressureIdx] = refPressure() + density*this->gravity()[1]*(globalPos[1] - this->gridGeometry().bBoxMin()[1]);
303 values[Indices::conti0EqIdx + 1] = refMoleFrac();
304 values[Indices::velocityXIdx] = xVelocity_(globalPos);
305
306#if NONISOTHERMAL
307 values[Indices::temperatureIdx] = refTemperature();
308#endif
309
310 return values;
311 }
312
314 const Scalar refVelocity() const
315 { return refVelocity_ ;}
316
318 const Scalar refPressure() const
319 { return refPressure_; }
320
322 const Scalar refMoleFrac() const
323 { return refMoleFrac_; }
324
326 const Scalar refTemperature() const
327 { return refTemperature_; }
328
329
330 void setTimeLoop(TimeLoopPtr timeLoop)
331 { timeLoop_ = timeLoop; }
332
333 Scalar time() const
334 { return timeLoop_->time(); }
335
340 Scalar permeability(const Element& element, const SubControlVolumeFace& scvf) const
341 {
342 return couplingManager().couplingData().darcyPermeability(element, scvf);
343 }
344
349 Scalar alphaBJ(const SubControlVolumeFace& scvf) const
350 {
351 return couplingManager().problem(CouplingManager::darcyIdx).spatialParams().beaversJosephCoeffAtPos(scvf.center());
352 }
353
354 // \}
355
356private:
357 bool onLeftBoundary_(const GlobalPosition &globalPos) const
358 { return globalPos[0] < this->gridGeometry().bBoxMin()[0] + eps_; }
359
360 bool onRightBoundary_(const GlobalPosition &globalPos) const
361 { return globalPos[0] > this->gridGeometry().bBoxMax()[0] - eps_; }
362
363 bool onLowerBoundary_(const GlobalPosition &globalPos) const
364 { return globalPos[1] < this->gridGeometry().bBoxMin()[1] + eps_; }
365
366 bool onUpperBoundary_(const GlobalPosition &globalPos) const
367 { return globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_; }
368
370 void updateFluidStateForBC_(FluidState& fluidState, const Scalar pressure) const
371 {
372 fluidState.setTemperature(refTemperature());
373 fluidState.setPressure(0, pressure);
374 fluidState.setSaturation(0, 1.0);
375 fluidState.setMoleFraction(0, 1, refMoleFrac());
376 fluidState.setMoleFraction(0, 0, 1.0 - refMoleFrac());
377
378 typename FluidSystem::ParameterCache paramCache;
379 paramCache.updatePhase(fluidState, 0);
380
381 const Scalar density = FluidSystem::density(fluidState, paramCache, 0);
382 fluidState.setDensity(0, density);
383
384 const Scalar molarDensity = FluidSystem::molarDensity(fluidState, paramCache, 0);
385 fluidState.setMolarDensity(0, molarDensity);
386
387 const Scalar enthalpy = FluidSystem::enthalpy(fluidState, paramCache, 0);
388 fluidState.setEnthalpy(0, enthalpy);
389 }
390
392 const Scalar xVelocity_(const GlobalPosition &globalPos) const
393 {
394 const Scalar vmax = refVelocity();
395 return 4 * vmax * (globalPos[1] - this->gridGeometry().bBoxMin()[1]) * (this->gridGeometry().bBoxMax()[1] - globalPos[1])
396 / (height_() * height_());
397 }
398
399 // the height of the free-flow domain
400 const Scalar height_() const
401 { return this->gridGeometry().bBoxMax()[1] - this->gridGeometry().bBoxMin()[1]; }
402
403 Scalar eps_;
404
405 Scalar refVelocity_;
406 Scalar refPressure_;
407 Scalar refMoleFrac_;
408 Scalar refTemperature_;
409 std::string problemName_;
410 TimeLoopPtr timeLoop_;
411
412 std::shared_ptr<CouplingManager> couplingManager_;
413
414 DiffusionCoefficientAveragingType diffCoeffAvgType_;
415};
416} // end namespace Dumux
417
418#endif // DUMUX_STOKES1P2C_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
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string pressure(int phaseIdx) noexcept
I/O name of pressure for multiphase systems.
Definition: name.hh:34
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
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
The component balance index that should be replaced by the total mass/mole balance.
Definition: common/properties.hh:104
Definition: common/properties.hh:169
The type of the fluid system to use.
Definition: common/properties.hh:223
const GravityVector & gravity() const
Returns the acceleration due to gravity.
Definition: dumux/freeflow/navierstokes/problem.hh:134
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
DiffusionCoefficientAveragingType
Defines which kind of averanging of diffusion coefficiencients (moleculat diffusion or thermal conduc...
Definition: couplingdata.hh:52
static DiffusionCoefficientAveragingType stringToEnum(DiffusionCoefficientAveragingType, const std::string &diffusionCoefficientAveragingType)
Convenience function to convert user input given as std::string to the corresponding enum class used ...
Definition: couplingdata.hh:60
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_2p2c/problem_stokes.hh:152
const Scalar refPressure() const
Returns the reference pressure.
Definition: 1p2c_2p2c/problem_stokes.hh:318
PrimaryVariables dirichletAtPos(const GlobalPosition &pos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: 1p2c_2p2c/problem_stokes.hh:218
void setTimeLoop(TimeLoopPtr timeLoop)
Definition: 1p2c_2p2c/problem_stokes.hh:330
Scalar alphaBJ(const SubControlVolumeFace &scvf) const
Returns the alpha value required as input parameter for the Beavers-Joseph-Saffman boundary condition...
Definition: 1p2c_2p2c/problem_stokes.hh:349
const Scalar refMoleFrac() const
Returns the reference mass fraction.
Definition: 1p2c_2p2c/problem_stokes.hh:322
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_2p2c/problem_stokes.hh:176
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_2p2c/problem_stokes.hh:123
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: 1p2c_2p2c/problem_stokes.hh:160
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_2p2c/problem_stokes.hh:235
const Scalar refTemperature() const
Returns the reference temperature.
Definition: 1p2c_2p2c/problem_stokes.hh:326
Scalar time() const
Definition: 1p2c_2p2c/problem_stokes.hh:333
const Scalar refVelocity() const
Returns the reference velocity.
Definition: 1p2c_2p2c/problem_stokes.hh:314
const std::string & name() const
The problem name.
Definition: 1p2c_2p2c/problem_stokes.hh:139
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_2p2c/problem_stokes.hh:340
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....