3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/porousmediumflow/1p/implicit/convergence/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 *****************************************************************************/
24#ifndef DUMUX_INCOMPRESSIBLE_ONEP_CONVERGENCETEST_PROBLEM_HH
25#define DUMUX_INCOMPRESSIBLE_ONEP_CONVERGENCETEST_PROBLEM_HH
26
27#include <cmath>
28#include <dune/grid/yaspgrid.hh>
29#include <dune/geometry/quadraturerules.hh>
30
34
38
41
42#include "spatialparams.hh"
43
44namespace Dumux {
45// forward declarations
46template<class TypeTag> class OnePTestProblem;
47
48namespace Properties {
49
50// create the type tag nodes
51namespace TTag {
52struct OnePIncompressible { using InheritsFrom = std::tuple<OneP>; };
53struct OnePIncompressibleTpfa { using InheritsFrom = std::tuple<OnePIncompressible, CCTpfaModel>; };
54struct OnePIncompressibleMpfa { using InheritsFrom = std::tuple<OnePIncompressible, CCMpfaModel>; };
55struct OnePIncompressibleBox { using InheritsFrom = std::tuple<OnePIncompressible, BoxModel>; };
56} // end namespace TTag
57
58// Set the grid type
59template<class TypeTag>
60struct Grid<TypeTag, TTag::OnePIncompressible> { using type = Dune::YaspGrid<2>; };
61
62// Set the problem type
63template<class TypeTag>
64struct Problem<TypeTag, TTag::OnePIncompressible> { using type = OnePTestProblem<TypeTag>; };
65
66// set the spatial params
67template<class TypeTag>
68struct SpatialParams<TypeTag, TTag::OnePIncompressible>
69{
73};
74
75// use the incompressible local residual (provides analytic jacobian)
76template<class TypeTag>
77struct LocalResidual<TypeTag, TTag::OnePIncompressible> { using type = OnePIncompressibleLocalResidual<TypeTag>; };
78
79// the fluid system
80template<class TypeTag>
81struct FluidSystem<TypeTag, TTag::OnePIncompressible>
82{
85};
86
87} // end namespace Properties
88
93template<class TypeTag>
94class OnePTestProblem : public PorousMediumFlowProblem<TypeTag>
95{
96 using ParentType = PorousMediumFlowProblem<TypeTag>;
97
102
104 using FVElementGeometry = typename GridGeometry::LocalView;
105 using SubControlVolume = typename GridGeometry::SubControlVolume;
106 using GridView = typename GridGeometry::GridView;
107 using Element = typename GridView::template Codim<0>::Entity;
108 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
109
110public:
115 OnePTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
117 {}
118
124 BoundaryTypes boundaryTypesAtPos(const GlobalPosition& globalPos) const
125 {
126 BoundaryTypes values;
127 values.setAllDirichlet();
128 return values;
129 }
130
135 PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
136 { return exact(globalPos); }
137
145 template<class ElementVolumeVariables>
146 NumEqVector source(const Element& element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& elemVolVars,
149 const SubControlVolume& scv) const
150 {
151 static const auto order = getParam<Scalar>("Problem.SourceIntegrationOrder");
152 static const auto periodLength = getParam<Scalar>("Problem.ExactSolPeriodLength");
153 const auto& k = this->spatialParams().permeabilityAtPos(scv.center());
154
155 using std::sin;
156 using std::cos;
157
158 const auto eg = element.geometry();
159 const auto rule = Dune::QuadratureRules<Scalar, GridView::dimension>::rule(eg.type(), order);
160
161 Scalar source = 0.0;
162 for (auto qp : rule)
163 {
164 const auto p = eg.global(qp.position());
165 const auto x = p[0];
166 const auto y = p[1];
167
168 const auto preFactor = -1.0*periodLength*periodLength*M_PI*M_PI;
169 const auto preFactorArg = periodLength*M_PI;
170 const auto sineTerm = sin(preFactorArg*x);
171 const auto cosTerm = cos(preFactorArg*y);
172 const auto secondDeriv = preFactor*sineTerm*cosTerm;
173
174 // derivative in x and y are identical
175 source -= 2.0*k*secondDeriv*qp.weight()*eg.integrationElement(qp.position());
176 }
177
178 source /= eg.volume();
179 return NumEqVector(source);
180 }
181
185 Scalar temperature() const
186 { return 283.15; }
187
192 static PrimaryVariables exact(const GlobalPosition& globalPos)
193 {
194 const auto x = globalPos[0];
195 const auto y = globalPos[1];
196
197 using std::sin;
198 using std::cos;
199
200 static const auto periodLength = getParam<Scalar>("Problem.ExactSolPeriodLength");
201 const auto preFactorArg = periodLength*M_PI;
202 const auto u = sin(preFactorArg*x)*cos(preFactorArg*y);
203
204 return PrimaryVariables(u);
205 }
206};
207
208} // end namespace Dumux
209
210#endif
Defines a type tag and some properties for models using the box scheme.
Properties for all models using cell-centered finite volume scheme with mpfa.
Properties for all models using cell-centered finite volume scheme with TPFA.
Setting constant fluid properties via the input file.
A liquid phase consisting of a single component.
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
Class to specify the type of a boundary.
Definition: common/boundarytypes.hh:38
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
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
Definition: common/properties.hh:91
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
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...
Definition: 1p/incompressiblelocalresidual.hh:41
Base class for all fully implicit porous media problems.
Definition: dumux/porousmediumflow/problem.hh:39
SpatialParams & spatialParams()
Returns the spatial parameters object.
Definition: dumux/porousmediumflow/problem.hh:146
NumEqVector source(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolume &scv) const
Evaluates the source term within a sub-control volume.
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:146
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Evaluates the boundary conditions for a Dirichlet control volume.
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:135
OnePTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
The constructor.
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:115
static PrimaryVariables exact(const GlobalPosition &globalPos)
Returns the exact solution at a position.
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:192
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/porousmediumflow/1p/implicit/convergence/problem.hh:124
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:185
The spatial parameters class for the test problem using the incompressible 1p model.
Definition: multidomain/boundary/darcydarcy/1p_1p/spatialparams.hh:62
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:52
std::tuple< OneP > InheritsFrom
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:52
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:53
std::tuple< OnePIncompressible, CCTpfaModel > InheritsFrom
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:53
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:54
std::tuple< OnePIncompressible, CCMpfaModel > InheritsFrom
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:54
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:55
std::tuple< OnePIncompressible, BoxModel > InheritsFrom
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:55
Dune::YaspGrid< 2 > type
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:60
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:71
GetPropType< TypeTag, Properties::GridGeometry > GridGeometry
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:70
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/porousmediumflow/1p/implicit/convergence/problem.hh:83
A single-phase, isothermal flow model using the fully implicit scheme.
Base class for all porous media problems.
Element-wise calculation of the residual and its derivatives for a single-phase, incompressible,...
Definition of the spatial parameters for the MaxwellStefan problem.