Loading [MathJax]/extensions/tex2jax.js
3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
test_1pproblem.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_1P_PROBLEM_HH
25#define DUMUX_TEST_1P_PROBLEM_HH
26
27#include <dune/grid/yaspgrid.hh>
28
31
35
37
38namespace Dumux
39{
40
41template<class TypeTag>
42class TestProblemOneP;
43
45// Specify the properties
47namespace Properties
48{
49NEW_TYPE_TAG(TestOneP, INHERITS_FROM(FVPressureOneP));
50
51// Set the grid type
52SET_TYPE_PROP(TestOneP, Grid, Dune::YaspGrid<2>);
53
54// the fluid system
55template<class TypeTag>
56struct FluidSystem<TypeTag, TTag::TestOneP>
57{
58 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
60};
61
62// Set the spatial parameters
64
65//Set the problem
67}
68
74template<class TypeTag>
75class TestProblemOneP: public DiffusionProblem1P<TypeTag >
76{
78 using TimeManager = typename GET_PROP_TYPE(TypeTag, TimeManager);
79
80 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
81 using Grid = typename GridView::Grid;
82
83 using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
84
85 using PrimaryVariables = typename GET_PROP_TYPE(TypeTag, PrimaryVariables);
86 using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
87
88 enum
89 {
90 dim = GridView::dimension, dimWorld = GridView::dimensionworld
91 };
92
93 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
94
95 using Element = typename GridView::Traits::template Codim<0>::Entity;
96 using Intersection = typename GridView::Intersection;
97 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
98 using LocalPosition = Dune::FieldVector<Scalar, dim>;
99
100
101public:
102 TestProblemOneP(TimeManager& timeManager, Grid& grid) :
103 ParentType(grid), velocity_(*this)
104 {
105 delta_ = getParam<Scalar>("Problem.Delta", 1e-3);
106
107 this->spatialParams().initialize(delta_);
108 }
109
113 // \{
114
120 std::string name() const
121 {
122 return "test_1p";
123 }
124
126 { return false; }
127
129 {
130 velocity_.calculateVelocity();
131 velocity_.addOutputVtkFields(this->resultWriter());
132 }
133
139 Scalar temperatureAtPos(const GlobalPosition& globalPos) const
140 {
141 return 273.15 + 10; // -> 10°C
142 }
143
144 // \}
145
147 Scalar referencePressureAtPos(const GlobalPosition& globalPos) const
148 {
149 return 1e5; // -> 10°C
150 }
151
153 void source(PrimaryVariables &values, const Element& element) const
154 {
155 values = 0;
156 values = integratedSource_(element, 4);
157 }
158
164 void boundaryTypes(BoundaryTypes &bcType,
165 const Intersection& intersection) const
166 {
167 bcType.setAllDirichlet();
168 }
169
171 void dirichletAtPos(PrimaryVariables &values,
172 const GlobalPosition &globalPos) const
173 {
174 values = exact(globalPos);
175 }
176
177
179 void neumann(PrimaryVariables &values, const Intersection& intersection) const
180 {
181 values = 0;
182 }
183
184private:
185 Scalar exact (const GlobalPosition& globalPos) const
186 {
187 double pi = 4.0*atan(1.0);
188 using std::sin;
189 return (sin(pi*globalPos[0])*sin(pi*globalPos[1]));
190 }
191
192 Dune::FieldVector<Scalar,dim> exactGrad (const GlobalPosition& globalPos) const
193 {
194 Dune::FieldVector<Scalar,dim> grad(0);
195 using std::sin;
196 using std::cos;
197 using std::atan;
198 double pi = 4.0*atan(1.0);
199 grad[0] = pi*cos(pi*globalPos[0])*sin(pi*globalPos[1]);
200 grad[1] = pi*cos(pi*globalPos[1])*sin(pi*globalPos[0]);
201
202 return grad;
203 }
204
205 Scalar integratedSource_(const Element& element, int integrationPoints) const
206 {
207 Scalar source = 0.;
208 LocalPosition localPos(0.0);
209 GlobalPosition globalPos(0.0);
210 Scalar halfInterval = 1.0/double(integrationPoints)/2.;
211 for (int i = 1; i <= integrationPoints; i++)
212 {
213 for (int j = 1; j <= integrationPoints; j++)
214 {
215 localPos[0] = double(i)/double(integrationPoints) - halfInterval;
216 localPos[1] = double(j)/double(integrationPoints) - halfInterval;
217 globalPos = element.geometry().global(localPos);
218 source += 1./(integrationPoints*integrationPoints) * evaluateSource_(globalPos);
219 }
220 }
221
222 return source;
223 }
224
225 Scalar evaluateSource_(const GlobalPosition& globalPos) const
226 {
227 Scalar temp = temperatureAtPos(globalPos);
228 Scalar referencePress = referencePressureAtPos(globalPos);
229
230 using std::sin;
231 using std::cos;
232 using std::atan;
233
234 Scalar pi = 4.0 * atan(1.0);
235 Scalar x = globalPos[0];
236 Scalar y = globalPos[1];
237
238 Scalar dpdx = pi * cos(pi * x) * sin(pi * y);
239 Scalar dpdy = pi * sin(pi * x) * cos(pi * y);
240 Scalar dppdxx = -pi * pi * sin(pi * x) * sin(pi * y);
241 Scalar dppdxy = pi * pi * cos(pi * x) * cos(pi * y);
242 Scalar dppdyx = dppdxy;
243 Scalar dppdyy = dppdxx;
244 Scalar kxx = (delta_* x*x + y*y)/(x*x + y*y);
245 Scalar kxy = -(1.0 - delta_) * x * y / (x*x + y*y);
246 Scalar kyy = (x*x + delta_*y*y)/(x*x + y*y);
247 Scalar dkxxdx = 2 * x * y*y * (delta_ - 1.0)/((x*x + y*y) * (x*x + y*y));
248 Scalar dkyydy = 2 * x*x * y * (delta_ - 1.0)/((x*x + y*y) * (x*x + y*y));
249 Scalar dkxydx = (1.0 - delta_) * y * (x*x - y*y) /((x*x + y*y) * (x*x + y*y));
250 Scalar dkxydy = (1.0 - delta_) * x * (y*y - x*x) /((x*x + y*y) * (x*x + y*y));
251
252 Scalar fx = dkxxdx * dpdx + kxx * dppdxx + dkxydx * dpdy + kxy * dppdyx;
253 Scalar fy = dkxydy * dpdx + kxy * dppdxy + dkyydy * dpdy + kyy * dppdyy;
254
255 return -(fx + fy) / FluidSystem::viscosity(temp, referencePress) * FluidSystem::density(temp, referencePress);
256 }
257
258 double delta_;
259 FVVelocity<TypeTag, typename GET_PROP_TYPE(TypeTag, Velocity) > velocity_;
260};
261} //end namespace
262
263#endif
#define GET_PROP_TYPE(TypeTag, PropTagName)
Definition: propertysystemmacros.hh:283
#define NEW_TYPE_TAG(...)
Definition: propertysystemmacros.hh:130
A liquid phase consisting of a single component.
Setting constant fluid properties via the input file.
spatial parameters for the test problem for diffusion models.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
SET_TYPE_PROP(FVPressureOneP, Velocity, FVVelocity1P< TypeTag >)
Set velocity reconstruction implementation standard cell centered finite volume schemes as default.
Property tag Velocity
The type velocity reconstruction.
Definition: porousmediumflow/sequential/properties.hh:67
Type tag FVPressureOneP INHERITS_FROM(PressureOneP))
The type tag for the one-phase problems using a standard finite volume model.
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
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 type of the spatial parameters object.
Definition: common/properties.hh:221
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
Base class for all single phase diffusion problems.
Definition: dumux/porousmediumflow/1p/sequential/diffusion/problem.hh:45
SpatialParams & spatialParams()
Returns the spatial parameters object.
Definition: dumux/porousmediumflow/1p/sequential/diffusion/problem.hh:214
void addOutputVtkFields(MultiWriter &writer)
Adds velocity output to the output file.
Definition: sequential/cellcentered/velocity.hh:71
void calculateVelocity()
Function which reconstructs a global velocity field.
Definition: sequential/cellcentered/velocity.hh:96
Base class for definition of an sequential diffusion (pressure) or transport problem.
Definition: onemodelproblem.hh:46
TimeManager & timeManager()
Returns TimeManager object used by the simulation.
Definition: onemodelproblem.hh:551
VtkMultiWriter & resultWriter()
Definition: onemodelproblem.hh:639
test problem for the sequential one-phase model.
Definition: test_1pproblem.hh:76
Scalar referencePressureAtPos(const GlobalPosition &globalPos) const
Returns the reference pressure for evaluation of constitutive relations.
Definition: test_1pproblem.hh:147
void boundaryTypes(BoundaryTypes &bcType, const Intersection &intersection) const
Returns the type of boundary condition.
Definition: test_1pproblem.hh:164
bool shouldWriteRestartFile() const
Definition: test_1pproblem.hh:125
void neumann(PrimaryVariables &values, const Intersection &intersection) const
return neumann condition (flux, [kg/(m^2 s)])
Definition: test_1pproblem.hh:179
std::string name() const
The problem name.
Definition: test_1pproblem.hh:120
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
return dirichlet condition (pressure, [Pa])
Definition: test_1pproblem.hh:171
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature within the domain.
Definition: test_1pproblem.hh:139
void addOutputVtkFields()
Definition: test_1pproblem.hh:128
void source(PrimaryVariables &values, const Element &element) const
source term [kg/(m^3 s)]
Definition: test_1pproblem.hh:153
TestProblemOneP(TimeManager &timeManager, Grid &grid)
Definition: test_1pproblem.hh:102
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_1pproblem.hh:58
spatial parameters for the test problem for 1-p diffusion models.
Definition: test_1pspatialparams.hh:38
Base class for all single phase diffusion problems.
Defines the properties required for finite volume pressure models.
Finite volume velocity reconstruction.