3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test_transportproblem.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 *****************************************************************************/
24#ifndef DUMUX_TEST_TRANSPORT_PROBLEM_HH
25#define DUMUX_TEST_TRANSPORT_PROBLEM_HH
26
27#include <dune/grid/yaspgrid.hh>
28#include <dune/grid/io/file/dgfparser/dgfyasp.hh>
29
32
35
37
38namespace Dumux
39{
44template<class TypeTag>
45class TestTransportProblem;
46
48// Specify the properties
50namespace Properties
51{
52NEW_TYPE_TAG(TransportTest, INHERITS_FROM(FVTransportTwoP, TestTransportSpatialParams));
53
54// Set the grid type
55SET_TYPE_PROP(TransportTest, Grid, Dune::YaspGrid<2>);
56
57// Set the problem property
59
60// Set the fluid system
61template<class TypeTag>
62struct FluidSystem<TypeTag, TTag::TransportTest>
63{
64 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
68};
69
71}
72
89template<class TypeTag>
91{
93 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
94 using Grid = typename GridView::Grid;
95
96 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
97
98 using TimeManager = typename GET_PROP_TYPE(TypeTag, TimeManager);
99 using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
100 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
101 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
102
103 using CellData = typename GET_PROP_TYPE(TypeTag, CellData);
104
105 enum
106 {
107 dimWorld = GridView::dimensionworld
108 };
109
110 enum
111 {
112 wPhaseIdx = Indices::wPhaseIdx, nPhaseIdx = Indices::nPhaseIdx
113 };
114
115 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
116
117 using Element = typename GridView::template Codim<0>::Entity;
118 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
119 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
120public:
121 TestTransportProblem(TimeManager& timeManager, Grid& grid) :
123 {}
124
125
127 void init()
128 {
130
131 VelocityVector vel(0);
132 vel[0] = 1e-5;
133
134 // compute update vector
135 for (const auto& element : elements(this->gridView()))
136 {
137 // cell index
138 int eIdxGlobal = this->elementMapper().index(element);
139
140 CellData& cellData = this->variables().cellData(eIdxGlobal);
141
142 // run through all intersections with neighbors and boundary
143 for (const auto& intersection : intersections(this->gridView(), element))
144 {
145 // local number of facet
146 int indexInInside = intersection.indexInInside();
147
148 cellData.fluxData().setVelocity(wPhaseIdx, indexInInside, vel);
149
150 const GlobalPosition& unitOuterNormal = intersection.centerUnitOuterNormal();
151
152 Scalar pot = vel * unitOuterNormal;
153
154 cellData.fluxData().setUpwindPotential(wPhaseIdx, indexInInside, pot);
155 cellData.fluxData().setUpwindPotential(nPhaseIdx, indexInInside, pot);
156 }
157 }
158 }
159
163 // \{
164
170 std::string name() const
171 {
172 return "test_transport";
173 }
174
176 {
177 return false;
178 }
179
185 Scalar temperatureAtPos(const GlobalPosition& globalPos) const
186 {
187 return 273.15 + 10; // -> 10°C
188 }
189
190 // \}
191
192
194 Scalar referencePressureAtPos(const GlobalPosition& globalPos) const
195 {
196 return 1e5; // -> 10°C
197 }
198
199 void sourceAtPos(PrimaryVariables &values,const GlobalPosition& globalPos) const
200 {
201 values = 0;
202 }
203
210 void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition& globalPos) const
211 {
212 if (globalPos[0] < eps_)
213 {
214 bcTypes.setAllDirichlet();
215 }
216 else if (globalPos[0] > this->bBoxMax()[0] - eps_)
217 {
218 bcTypes.setAllOutflow();
219 }
220 // all other boundaries
221 else
222 {
223 bcTypes.setAllNeumann();
224 }
225 }
226
228 void dirichletAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
229 {
230 values = 0;
231 if (globalPos[0] < eps_)
232 {
233 values = 1.0;
234 }
235 }
236
238 void neumannAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
239 {
240 values = 0;
241 }
242
244 void initialAtPos(PrimaryVariables &values,
245 const GlobalPosition &globalPos) const
246 {
247 values = 0;
248 }
249
250private:
251 static constexpr Scalar eps_ = 1e-6;
252};
253} //end namespace
254
255#endif
#define GET_PROP(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:281
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
#define NEW_TYPE_TAG(...)
Definition: propertysystemmacros.hh:130
Setting constant fluid properties via the input file.
A liquid phase consisting of a single component.
spatial parameters for the explicit transport test
Dune::IteratorRange< typename MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper >::Intersections::const_iterator > intersections(const MultiDomainGlue< DomainGridView, TargetGridView, DomainMapper, TargetMapper > &glue)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: glue.hh:62
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
SET_INT_PROP(SequentialOneP, NumEq, 1)
Set number of equations to 1 for isothermal one-phase models.
SET_TYPE_PROP(FVPressureOneP, Velocity, FVVelocity1P< TypeTag >)
Set velocity reconstruction implementation standard cell centered finite volume schemes as default.
Type tag FVPressureOneP INHERITS_FROM(PressureOneP))
The type tag for the one-phase problems using a standard finite volume model.
Property tag VelocityFormulation
The type of velocity reconstructed for the transport model.
Definition: porousmediumflow/2p/sequential/properties.hh:54
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 type of the fluid system to use.
Definition: common/properties.hh:223
A liquid phase consisting of a single component.
Definition: 1pliquid.hh:46
A fluid system for two-phase models assuming immiscibility and thermodynamic equilibrium.
Definition: 2pimmiscible.hh:59
static const int velocityTotal
Indicates total velocity.
Definition: porousmediumflow/2p/sequential/indices.hh:63
Base class for a sequential two-phase transport problem.
Definition: dumux/porousmediumflow/2p/sequential/transport/problem.hh:51
Base class for definition of an sequential diffusion (pressure) or transport problem.
Definition: onemodelproblem.hh:46
const ElementMapper & elementMapper() const
Returns the mapper for elements to indices.
Definition: onemodelproblem.hh:531
TimeManager & timeManager()
Returns TimeManager object used by the simulation.
Definition: onemodelproblem.hh:551
void init()
Called by the TimeManager in order to initialize the problem.
Definition: onemodelproblem.hh:327
const GlobalPosition & bBoxMax() const
The coordinate of the corner of the GridView's bounding box with the largest values.
Definition: onemodelproblem.hh:545
const GridView & gridView() const
The GridView which used by the problem.
Definition: onemodelproblem.hh:519
Variables & variables()
Returns variables object.
Definition: onemodelproblem.hh:563
test problem for the explicit transport model
Definition: test_transportproblem.hh:91
void neumannAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set neumann condition for phases (flux, [kg/(m^2 s)])
Definition: test_transportproblem.hh:238
bool shouldWriteRestartFile() const
Definition: test_transportproblem.hh:175
void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition &globalPos) const
Returns the type of boundary condition.
Definition: test_transportproblem.hh:210
std::string name() const
The problem name.
Definition: test_transportproblem.hh:170
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature within the domain.
Definition: test_transportproblem.hh:185
Scalar referencePressureAtPos(const GlobalPosition &globalPos) const
Returns the reference pressure for evaluation of constitutive relations.
Definition: test_transportproblem.hh:194
void sourceAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Definition: test_transportproblem.hh:199
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set dirichlet condition (saturation [-])
Definition: test_transportproblem.hh:228
TestTransportProblem(TimeManager &timeManager, Grid &grid)
Definition: test_transportproblem.hh:121
void init()
set initial velocity field -> v_total is assumed to be constant in this test!
Definition: test_transportproblem.hh:127
void initialAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
return initial solution
Definition: test_transportproblem.hh:244
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_transportproblem.hh:64
spatial parameters for the explicit transport test
Definition: test_transportspatialparams.hh:68
Specifies the properties for immiscible 2p transport.
Base class for two-phase transport problems.