3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokesnc/densitydrivenflow/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 *****************************************************************************/
25#ifndef DUMUX_DENSITY_FLOW_NC_TEST_PROBLEM_HH
26#define DUMUX_DENSITY_FLOW_NC_TEST_PROBLEM_HH
27
28#include <dune/grid/yaspgrid.hh>
29
33
37
38
39namespace Dumux {
40
41template <class TypeTag>
42class DensityDrivenFlowProblem;
43
44namespace Properties {
45
46// Create new type tags
47namespace TTag {
48struct DensityDrivenFlow { using InheritsFrom = std::tuple<NavierStokesNC, StaggeredFreeFlowModel>; };
49} // end namespace TTag
50
51// Select the fluid system
52template<class TypeTag>
53struct FluidSystem<TypeTag, TTag::DensityDrivenFlow>
54{
56 static constexpr int phaseIdx = H2OAir::liquidPhaseIdx;
58};
59
60template<class TypeTag>
61struct ReplaceCompEqIdx<TypeTag, TTag::DensityDrivenFlow> { static constexpr int value = 0; };
62
63// Set the grid type
64template<class TypeTag>
65struct Grid<TypeTag, TTag::DensityDrivenFlow> { using type = Dune::YaspGrid<2>; };
66
67// Set the problem property
68template<class TypeTag>
69struct Problem<TypeTag, TTag::DensityDrivenFlow> { using type = Dumux::DensityDrivenFlowProblem<TypeTag> ; };
70
71template<class TypeTag>
72struct EnableGridGeometryCache<TypeTag, TTag::DensityDrivenFlow> { static constexpr bool value = true; };
73
74template<class TypeTag>
75struct EnableGridFluxVariablesCache<TypeTag, TTag::DensityDrivenFlow> { static constexpr bool value = true; };
76template<class TypeTag>
77struct EnableGridVolumeVariablesCache<TypeTag, TTag::DensityDrivenFlow> { static constexpr bool value = true; };
78
79template<class TypeTag>
80struct UseMoles<TypeTag, TTag::DensityDrivenFlow> { static constexpr bool value = true; };
81} // end namespace Properties
82
94template <class TypeTag>
96{
97 using ParentType = NavierStokesProblem<TypeTag>;
98
106
107 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
108 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
109
110 using TimeLoopPtr = std::shared_ptr<CheckPointTimeLoop<Scalar>>;
111
112 static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
113
114 static constexpr auto transportCompIdx = Indices::conti0EqIdx + 1;
115 static constexpr auto transportEqIdx = Indices::conti0EqIdx + 1;
116
117public:
118 DensityDrivenFlowProblem(std::shared_ptr<const GridGeometry> gridGeometry)
119 : ParentType(gridGeometry), eps_(1e-6)
120 {
121 useWholeLength_ = getParam<bool>("Problem.UseWholeLength");
122 FluidSystem::init();
123 deltaRho_.resize(this->gridGeometry().numCellCenterDofs());
124 }
125
129 // \{
130
132 {
133 return false;
134 }
135
141 Scalar temperature() const
142 { return 273.15 + 10; } // 10C
143
149 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
150 {
151 return NumEqVector(0.0);
152 }
153
154 // \}
158 // \{
159
166 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
167 {
168 BoundaryTypes values;
169
170 // set Dirichlet values for the velocity everywhere
171 values.setDirichlet(Indices::velocityXIdx);
172 values.setDirichlet(Indices::velocityYIdx);
173 values.setNeumann(Indices::conti0EqIdx);
174 values.setNeumann(transportEqIdx);
175
176 if(globalPos[1] > this->gridGeometry().bBoxMax()[1] - eps_)
177 {
178 if(useWholeLength_)
179 values.setDirichlet(transportCompIdx);
180 else
181 if(globalPos[0] > 0.4 && globalPos[0] < 0.6)
182 values.setDirichlet(transportCompIdx);
183 }
184
185 return values;
186 }
187
196 template<class FVElementGeometry, class SubControlVolume>
197 bool isDirichletCell(const Element& element,
198 const FVElementGeometry& fvGeometry,
199 const SubControlVolume& scv,
200 int pvIdx) const
201 {
202 // set a fixed pressure in one cell
203 return (isLowerLeftCell_(scv) && pvIdx == Indices::pressureIdx);
204 }
205
211 PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
212 {
213 PrimaryVariables values;
214
215 values[Indices::pressureIdx] = 1.1e+5;
216 values[transportCompIdx] = 1e-3;
217 values[Indices::velocityXIdx] = 0.0;
218 values[Indices::velocityYIdx] = 0.0;
219
220 return values;
221 }
222
223 // \}
224
228 // \{
229
235 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
236 {
237 PrimaryVariables values;
238 values[Indices::pressureIdx] = 1.1e+5;
239 values[transportCompIdx] = 0.0;
240 values[Indices::velocityXIdx] = 0.0;
241 values[Indices::velocityYIdx] = 0.0;
242
243 return values;
244 }
245
254 template<class GridVariables, class SolutionVector>
255 void calculateDeltaRho(const GridVariables& gridVariables, const SolutionVector& sol)
256 {
257 for (const auto& element : elements(this->gridGeometry().gridView()))
258 {
259 auto fvGeometry = localView(this->gridGeometry());
260 fvGeometry.bindElement(element);
261 for (auto&& scv : scvs(fvGeometry))
262 {
263 auto ccDofIdx = scv.dofIndex();
264
265 auto elemVolVars = localView(gridVariables.curGridVolVars());
266 elemVolVars.bind(element, fvGeometry, sol);
267
268 deltaRho_[ccDofIdx] = elemVolVars[scv].density() - 999.694;
269 }
270 }
271 }
272
273 const std::vector<Scalar>& getDeltaRho() const
274 { return deltaRho_; }
275
276
277
278 // \}
279
280private:
281
282 template<class SubControlVolume>
283 bool isLowerLeftCell_(const SubControlVolume& scv) const
284 { return scv.dofIndex() == 0; }
285
286 const Scalar eps_;
287 bool useWholeLength_;
288 std::vector<Scalar> deltaRho_;
289};
290} // end namespace Dumux
291
292#endif
A much simpler (and thus potentially less buggy) version of pure water.
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...
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
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 of the fluid system to use.
Definition: common/properties.hh:223
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
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 one-phase model.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:96
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:211
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:141
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:149
bool shouldWriteRestartFile() const
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:131
DensityDrivenFlowProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:118
bool isDirichletCell(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, int pvIdx) const
Returns whether a fixed Dirichlet value shall be used at a given cell.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:197
void calculateDeltaRho(const GridVariables &gridVariables, const SolutionVector &sol)
Adds additional VTK output data to the VTKWriter.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:255
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:235
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary cont...
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:166
const std::vector< Scalar > & getDeltaRho() const
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:273
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:48
std::tuple< NavierStokesNC, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:48
Dune::YaspGrid< 2 > type
Definition: test/freeflow/navierstokesnc/densitydrivenflow/problem.hh:65
Defines a type tag and some properties for ree-flow models using the staggered scheme....