3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test_diffusionproblem.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_2P_PROBLEM_HH
25#define DUMUX_TEST_2P_PROBLEM_HH
26
27#include <dune/grid/yaspgrid.hh>
28#include <dune/grid/utility/structuredgridfactory.hh>
29
31
37
39
40namespace Dumux {
41
45template<class TypeTag>
46class TestDiffusionProblem;
47
49// Specify the properties
51namespace Properties
52{
54NEW_TYPE_TAG(FVVelocity2PTestTypeTag, INHERITS_FROM(FVPressureTwoP, TestDiffusionSpatialParams));
56
57// Set the grid type
58SET_TYPE_PROP(FVVelocity2PTestTypeTag, Grid, Dune::YaspGrid<2>);
59
60// Set the fluid system
61template<class TypeTag>
62struct FluidSystem<TypeTag, TTag::FVVelocity2PTestTypeTag>
63{
64 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
68};
69
70// set the types for the MPFA-O FV method
71NEW_TYPE_TAG(FVMPFAOVelocity2PTestTypeTag, INHERITS_FROM(FvMpfaO2dPressureTwoP, TestDiffusionSpatialParams));
72//SET_TYPE_PROP(FVMPFAOVelocity2PTestTypeTag, LinearSolver, ILUnBiCGSTABBackend);
73SET_TYPE_PROP(FVMPFAOVelocity2PTestTypeTag, LinearSolver, SSORBiCGSTABBackend);
74SET_TYPE_PROP(FVMPFAOVelocity2PTestTypeTag, Problem, TestDiffusionProblem<TypeTag>);
75// Set the grid type
76SET_TYPE_PROP(FVMPFAOVelocity2PTestTypeTag, Grid, Dune::YaspGrid<2>);
77
78// Set the fluid system
79template<class TypeTag>
80struct FluidSystem<TypeTag, TTag::FVMPFAOVelocity2PTestTypeTag>
81{
82 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
86};
87
88// set the types for the mimetic FD method
89NEW_TYPE_TAG(MimeticPressure2PTestTypeTag, INHERITS_FROM(MimeticPressureTwoP, TestDiffusionSpatialParams));
90SET_TYPE_PROP(MimeticPressure2PTestTypeTag, Problem, TestDiffusionProblem<TypeTag>);
91
92// Set the grid type
93SET_TYPE_PROP(MimeticPressure2PTestTypeTag, Grid, Dune::YaspGrid<2>);
94
95// Set the fluid system
96template<class TypeTag>
97struct FluidSystem<TypeTag, TTag::MimeticPressure2PTestTypeTag>
98{
99 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
103};
104
105} // end namespace Properties
106
115template<class TypeTag>
117{
119 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
120 using Grid = typename GridView::Grid;
121
122 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
123
124 using WettingPhase = typename GET_PROP(TypeTag, FluidSystem)::WettingPhase;
125
126 using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
127
128 enum
129 {
130 dim = GridView::dimension, dimWorld = GridView::dimensionworld
131 };
132
133 enum
134 {
135 wPhaseIdx = Indices::wPhaseIdx,
136 pwIdx = Indices::pwIdx,
137 swIdx = Indices::swIdx
138 };
139
140 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
141
142 using Element = typename GridView::Traits::template Codim<0>::Entity;
143 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
144 using LocalPosition = Dune::FieldVector<Scalar, dim>;
145
146public:
147 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
148 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
149 using ScalarSolution = typename SolutionTypes::ScalarSolution;
150
152 ParentType(grid), velocity_(*this)
153 {
154 delta_ = getParam<Scalar>("Problem.Delta", 1e-3);
155 }
156
158 void init()
159 {
160 this->variables().initialize();
161 this->spatialParams().initialize(delta_);
162 for (int i = 0; i < this->gridView().size(0); i++)
163 {
164 this->variables().cellData(i).setSaturation(wPhaseIdx, 1.0);
165 }
166 this->model().initialize();
167 velocity_.initialize();
168 }
169
173 // \{
174
176 { return false; }
177
178
180 {
181 velocity_.calculateVelocity();
182// velocity_.addOutputVtkFields(this->resultWriter());
183 }
184
187 {
188 ScalarSolution *exactPressure = this->resultWriter().allocateManagedBuffer(this->gridView().size(0));
189
190 for(const auto& element : elements(this->gridView()))
191 {
192 (*exactPressure)[this->elementMapper().index(element)][0] = exact(element.geometry().center());
193 }
194
195 this->resultWriter().attachCellData(*exactPressure, "exact pressure");
196
197 this->spatialParams().addOutputVtkFields(this->resultWriter());
198
199 return;
200 }
201
207 Scalar temperatureAtPos(const GlobalPosition& globalPos) const
208 {
209 return 273.15 + 10; // -> 10°C
210 }
211
212 // \}
213
214
216 Scalar referencePressureAtPos(const GlobalPosition& globalPos) const
217 {
218 return 1e5; // -> 10°C
219 }
220
221 void source(PrimaryVariables &values,const Element& element) const
222 {
223 values = 0;
224
225 values[wPhaseIdx] = integratedSource_(element, 4);
226 }
227
234 void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition& globalPos) const
235 {
236
237 bcTypes.setAllDirichlet();
238 }
239
241 void dirichletAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
242 {
243 values[pwIdx] = exact(globalPos);
244 values[swIdx] = 1.0;
245 }
246
248 void neumannAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
249 {
250 values = 0;
251 }
252
253 Scalar exact (const GlobalPosition& globalPos) const
254 {
255 using std::sin;
256 using std::atan;
257
258 Scalar pi = 4.0*atan(1.0);
259
260 return (sin(pi*globalPos[0])*sin(pi*globalPos[1]));
261 }
262
263 Dune::FieldVector<Scalar,dim> exactGrad (const GlobalPosition& globalPos) const
264 {
265 Dune::FieldVector<Scalar,dim> grad(0);
266 using std::sin;
267 using std::cos;
268 using std::atan;
269
270 Scalar pi = 4.0*atan(1.0);
271 grad[0] = pi*cos(pi*globalPos[0])*sin(pi*globalPos[1]);
272 grad[1] = pi*cos(pi*globalPos[1])*sin(pi*globalPos[0]);
273
274 return grad;
275 }
276
277private:
278
279 Scalar integratedSource_(const Element& element, int integrationPoints) const
280 {
281 Scalar source = 0.;
282 LocalPosition localPos(0.0);
283 GlobalPosition globalPos(0.0);
284 Scalar halfInterval = 1.0/double(integrationPoints)/2.;
285 for (int i = 1; i <= integrationPoints; i++)
286 {
287 for (int j = 1; j <= integrationPoints; j++)
288 {
289 localPos[0] = double(i)/double(integrationPoints) - halfInterval;
290 localPos[1] = double(j)/double(integrationPoints) - halfInterval;
291 globalPos = element.geometry().global(localPos);
292 source += 1./(integrationPoints*integrationPoints) * evaluateSource_(globalPos);
293 }
294 }
295
296 return source;
297 }
298
299 Scalar evaluateSource_(const GlobalPosition& globalPos) const
300 {
301 Scalar temp = temperatureAtPos(globalPos);
302 Scalar referencePress = referencePressureAtPos(globalPos);
303
304 using std::sin;
305 using std::cos;
306 using std::atan;
307
308 Scalar pi = 4.0 * atan(1.0);
309 Scalar x = globalPos[0];
310 Scalar y = globalPos[1];
311
312 Scalar dpdx = pi * cos(pi * x) * sin(pi * y);
313 Scalar dpdy = pi * sin(pi * x) * cos(pi * y);
314 Scalar dppdxx = -pi * pi * sin(pi * x) * sin(pi * y);
315 Scalar dppdxy = pi * pi * cos(pi * x) * cos(pi * y);
316 Scalar dppdyx = dppdxy;
317 Scalar dppdyy = dppdxx;
318 Scalar kxx = (delta_* x*x + y*y)/(x*x + y*y);
319 Scalar kxy = -(1.0 - delta_) * x * y / (x*x + y*y);
320 Scalar kyy = (x*x + delta_*y*y)/(x*x + y*y);
321 Scalar dkxxdx = 2 * x * y*y * (delta_ - 1.0)/((x*x + y*y) * (x*x + y*y));
322 Scalar dkyydy = 2 * x*x * y * (delta_ - 1.0)/((x*x + y*y) * (x*x + y*y));
323 Scalar dkxydx = (1.0 - delta_) * y * (x*x - y*y) /((x*x + y*y) * (x*x + y*y));
324 Scalar dkxydy = (1.0 - delta_) * x * (y*y - x*x) /((x*x + y*y) * (x*x + y*y));
325
326 Scalar fx = dkxxdx * dpdx + kxx * dppdxx + dkxydx * dpdy + kxy * dppdyx;
327 Scalar fy = dkxydy * dpdx + kxy * dppdxy + dkyydy * dpdy + kyy * dppdyy;
328
329 return -(fx + fy) / WettingPhase::viscosity(temp, referencePress) * WettingPhase::density(temp, referencePress);
330 }
331
332 Scalar delta_;
333 FVVelocity<TypeTag, typename GET_PROP_TYPE(TypeTag, Velocity) > velocity_;
334};
335} //end namespace
336
337#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.
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
TODO: Remove this property as soon as the decoupled models are integrated.
Definition: common/properties.hh:95
The type of the fluid system to use.
Definition: common/properties.hh:223
Dune::BlockVector< Dune::FieldVector< Scalar, nComp > > * allocateManagedBuffer(int nEntities)
Allocate a managed buffer for a scalar field.
Definition: vtkmultiwriter.hh:140
void attachCellData(DataBuffer &buf, std::string name, int nComps=1)
Add a cell centered quantity to the output.
Definition: vtkmultiwriter.hh:205
Sequential SSOR-preconditioned BiCGSTAB solver.
Definition: seqsolverbackend.hh:246
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
Base class for stationary solution of a two-phase diffusion/pressure equation.
Definition: dumux/porousmediumflow/2p/sequential/diffusion/problem.hh:40
SpatialParams & spatialParams()
Returns the spatial parameters object.
Definition: dumux/porousmediumflow/2p/sequential/diffusion/problem.hh:213
void calculateVelocity()
Function which reconstructs a global velocity field.
Definition: sequential/cellcentered/velocity.hh:96
void initialize()
Initialize velocity implementation.
Definition: sequential/cellcentered/velocity.hh:56
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
VtkMultiWriter & resultWriter()
Definition: onemodelproblem.hh:639
Model & model()
Returns numerical model used for the problem.
Definition: onemodelproblem.hh:575
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 diffusion models from the FVCA5 benchmark.
Definition: test_diffusionproblem.hh:117
TestDiffusionProblem(Grid &grid)
Definition: test_diffusionproblem.hh:151
void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition &globalPos) const
Returns the type of boundary condition.
Definition: test_diffusionproblem.hh:234
void source(PrimaryVariables &values, const Element &element) const
Definition: test_diffusionproblem.hh:221
bool shouldWriteRestartFile() const
Definition: test_diffusionproblem.hh:175
Dune::FieldVector< Scalar, dim > exactGrad(const GlobalPosition &globalPos) const
Definition: test_diffusionproblem.hh:263
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature within the domain.
Definition: test_diffusionproblem.hh:207
void addOutputVtkFields()
Definition: test_diffusionproblem.hh:186
Scalar exact(const GlobalPosition &globalPos) const
Definition: test_diffusionproblem.hh:253
void init()
for this specific problem: initialize the saturation and afterwards the model
Definition: test_diffusionproblem.hh:158
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set dirichlet condition (saturation [-])
Definition: test_diffusionproblem.hh:241
Scalar referencePressureAtPos(const GlobalPosition &globalPos) const
Returns the reference pressure for evaluation of constitutive relations.
Definition: test_diffusionproblem.hh:216
typename SolutionTypes::ScalarSolution ScalarSolution
Definition: test_diffusionproblem.hh:149
void calculateFVVelocity()
Definition: test_diffusionproblem.hh:179
void neumannAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set neumann condition for phases (flux, [kg/(m^2 s)])
Definition: test_diffusionproblem.hh:248
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_diffusionproblem.hh:64
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_diffusionproblem.hh:82
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_diffusionproblem.hh:99
spatial parameters for the test problem for diffusion models.
Definition: test_diffusionspatialparams.hh:65
Base class for stationary solution of a two-phase diffusion/pressure equation.
Defines the properties required for finite volume pressure models in a two-phase sequential model.
Defines the properties required for (immiscible) twophase sequential models.
Finite volume velocity reconstruction.
Properties for the MPFA-O method.