3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test_diffusionproblem3d.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2/*****************************************************************************
3 * See the file COPYING for full copying permissions. *
4 * *
5 * This program is free software: you can redistribute it and/or modify *
6 * it under the terms of the GNU General Public License as published by *
7 * the Free Software Foundation, either version 3 of the License, or *
8 * (at your option) any later version. *
9 * *
10 * This program is distributed in the hope that it will be useful, *
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
13 * GNU General Public License for more details. *
14 * *
15 * You should have received a copy of the GNU General Public License *
16 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
17*****************************************************************************/
22#ifndef DUMUX_TEST_DIFFUSION_3D_PROBLEM_HH
23#define DUMUX_TEST_DIFFUSION_3D_PROBLEM_HH
24
25#if HAVE_DUNE_ALUGRID
26#include <dune/alugrid/grid.hh>
27#endif
28#if HAVE_UG
29#include <dune/grid/uggrid.hh>
30#endif
31#include <dune/grid/yaspgrid.hh>
32
34
38
41
43
44namespace Dumux
45{
49template<class TypeTag>
50class TestDiffusion3DProblem;
51
53// Specify the properties
55namespace Properties
56{
58
59// Set the grid type
60#if HAVE_DUNE_ALUGRID
61SET_TYPE_PROP(DiffusionTest, Grid, Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming>);
62#elif HAVE_UG
63SET_TYPE_PROP(DiffusionTest, Grid, Dune::UGGrid<3>);
64#else
65SET_TYPE_PROP(DiffusionTest, Grid, Dune::YaspGrid<3>);
66#endif
67
69
70// Set the fluid system
71template<class TypeTag>
72struct FluidSystem<TypeTag, TTag::DiffusionTest>
73{
74 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
78};
79
80#if HAVE_SUPERLU
81SET_TYPE_PROP(DiffusionTest, LinearSolver, SuperLUBackend);
82#else
84#endif
85
86// set the types for the 2PFA FV method
87NEW_TYPE_TAG(FVTest, INHERITS_FROM(FVPressureTwoP, DiffusionTest));
88
89// set the types for the MPFA-L FV method
90NEW_TYPE_TAG(FVMPFAL3DTestTypeTag, INHERITS_FROM(FvMpfaL3dPressureTwoP, DiffusionTest));
91
92// set the types for the mimetic FD method
93NEW_TYPE_TAG(MimeticTest, INHERITS_FROM(MimeticPressureTwoP, DiffusionTest));
94}
95
101template<class TypeTag>
103{
105 using GridView = typename GET_PROP_TYPE(TypeTag, GridView);
106 using Grid = typename GridView::Grid;
107
108 using Indices = typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
109
110 using FluidSystem = typename GET_PROP_TYPE(TypeTag, FluidSystem);
111 using FluidState = typename GET_PROP_TYPE(TypeTag, FluidState);
112
113 using BoundaryTypes = typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
114
115 enum
116 {
117 dim = GridView::dimension, dimWorld = GridView::dimensionworld
118 };
119
120 enum
121 {
122 wPhaseIdx = Indices::wPhaseIdx,
123 nPhaseIdx = Indices::nPhaseIdx,
124 pWIdx = Indices::pressureIdx,
125 swIdx = Indices::swIdx,
126 pressureEqIdx = Indices::pressureEqIdx,
127 };
128
129 using Scalar = typename GET_PROP_TYPE(TypeTag, Scalar);
130
131 using Element = typename GridView::Traits::template Codim<0>::Entity;
132 using Intersection = typename GridView::Intersection;
133 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
134
135public:
136 using SolutionTypes = typename GET_PROP(TypeTag, SolutionTypes);
137 using PrimaryVariables = typename SolutionTypes::PrimaryVariables;
138 using ScalarSolution = typename SolutionTypes::ScalarSolution;
139
141 ParentType(grid), velocity_(*this)
142 { }
143
145 void init()
146 {
147 this->variables().initialize();
148 for (int i = 0; i < this->gridView().size(0); i++)
149 {
150 this->variables().cellData(i).setSaturation(wPhaseIdx, 1.0);
151 this->variables().cellData(i).setSaturation(nPhaseIdx, 0.0);
152 }
153 this->model().initialize();
154 }
155
157 {
158 velocity_.calculateVelocity();
159// velocity_.addOutputVtkFields(this->resultWriter());
160 }
161
164 {
165 ScalarSolution *exactPressure = this->resultWriter().allocateManagedBuffer(this->gridView().size(0));
166
167 for(const auto& element : elements(this->gridView()))
168 {
169 (*exactPressure)[this->elementMapper().index(element)][0] = exact(element.geometry().center());
170 }
171
172 this->resultWriter().attachCellData(*exactPressure, "exact pressure");
173
174 return;
175 }
176
180 // \{
181
183 { return false; }
184
190 Scalar temperatureAtPos(const GlobalPosition& globalPos) const
191 {
192 return 273.15 + 10; // -> 10°C
193 }
194
195 // \}
196
197
199 Scalar referencePressureAtPos(const GlobalPosition& globalPos) const
200 {
201 return 1e5; // -> 10°C
202 }
203
204 void sourceAtPos(PrimaryVariables &values,const GlobalPosition& globalPos) const
205 {
206 values = 0;
207 using std::sin;
208 using std::cos;
209 using std::atan;
210
211 double pi = 4.0*atan(1.0);
212
213 values[wPhaseIdx] = -pi*pi*cos(pi*globalPos[0])*cos(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0))+
214 3.0*pi*pi*sin(pi*globalPos[0])*sin(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0))-
215 pi*pi*sin(pi*globalPos[0])*cos(pi*(globalPos[1]+1.0/2.0))*cos(pi*(globalPos[2]+1.0/3.0));
216 }
217
224 void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition& globalPos) const
225 {
226 bcTypes.setAllDirichlet();
227 }
228
230 void dirichletAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
231 {
232 values[pWIdx] = exact(globalPos);
233 values[swIdx] = 1.0;
234 }
235
237 void neumannAtPos(PrimaryVariables &values, const GlobalPosition& globalPos) const
238 {
239 values = 0;
240 }
241
242 Scalar exact (const GlobalPosition& globalPos) const
243 {
244 using std::sin;
245 using std::atan;
246
247 double pi = 4.0*atan(1.0);
248
249 return (1.0+sin(pi*globalPos[0])*sin(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0)));
250 }
251
252 Dune::FieldVector<Scalar,dim> exactGrad (const GlobalPosition& globalPos) const
253 {
254 Dune::FieldVector<Scalar,dim> grad(0);
255 using std::sin;
256 using std::cos;
257 using std::atan;
258
259 double pi = 4.0*atan(1.0);
260 grad[0] = pi*cos(pi*globalPos[0])*sin(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0));
261 grad[1] = pi*sin(pi*globalPos[0])*cos(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0));
262 grad[2] = pi*sin(pi*globalPos[0])*sin(pi*(globalPos[1]+1.0/2.0))*cos(pi*(globalPos[2]+1.0/3.0));
263
264 return grad;
265 }
266
267private:
268 FVVelocity<TypeTag, typename GET_PROP_TYPE(TypeTag, Velocity) > velocity_;
269 static constexpr Scalar eps_ = 1e-4;
270
271};
272} //end namespace
273
274#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.
Properties for the MPFA-O method.
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.
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 ILU(n)-preconditioned GMRes solver.
Definition: seqsolverbackend.hh:691
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
Base class for finite volume velocity reconstruction.
Definition: sequential/cellcentered/velocity.hh:48
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
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 FVCA6 benchmark.
Definition: test_diffusionproblem3d.hh:103
void addOutputVtkFields()
Definition: test_diffusionproblem3d.hh:163
typename SolutionTypes::ScalarSolution ScalarSolution
Definition: test_diffusionproblem3d.hh:138
TestDiffusion3DProblem(Grid &grid)
Definition: test_diffusionproblem3d.hh:140
Dune::FieldVector< Scalar, dim > exactGrad(const GlobalPosition &globalPos) const
Definition: test_diffusionproblem3d.hh:252
void calculateFVVelocity()
Definition: test_diffusionproblem3d.hh:156
void dirichletAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set dirichlet condition (saturation [-])
Definition: test_diffusionproblem3d.hh:230
Scalar exact(const GlobalPosition &globalPos) const
Definition: test_diffusionproblem3d.hh:242
void sourceAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
Definition: test_diffusionproblem3d.hh:204
void init()
for this specific problem: initialize the saturation and afterwards the model
Definition: test_diffusionproblem3d.hh:145
bool shouldWriteRestartFile() const
Definition: test_diffusionproblem3d.hh:182
Scalar temperatureAtPos(const GlobalPosition &globalPos) const
Returns the temperature within the domain.
Definition: test_diffusionproblem3d.hh:190
Scalar referencePressureAtPos(const GlobalPosition &globalPos) const
Returns the reference pressure for evaluation of constitutive relations.
Definition: test_diffusionproblem3d.hh:199
void neumannAtPos(PrimaryVariables &values, const GlobalPosition &globalPos) const
set neumann condition for phases (flux, [kg/(m^2 s)])
Definition: test_diffusionproblem3d.hh:237
void boundaryTypesAtPos(BoundaryTypes &bcTypes, const GlobalPosition &globalPos) const
Returns the type of boundary condition.
Definition: test_diffusionproblem3d.hh:224
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition: test_diffusionproblem3d.hh:74
spatial parameters for the test problem for diffusion models.
Definition: test_diffusionspatialparams3d.hh:64
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.