24#ifndef DUMUX_TEST_2P_PROBLEM_HH
25#define DUMUX_TEST_2P_PROBLEM_HH
27#include <dune/grid/yaspgrid.hh>
28#include <dune/grid/utility/structuredgridfactory.hh>
45template<
class TypeTag>
46class TestDiffusionProblem;
61template<
class TypeTag>
79template<
class TypeTag>
80struct FluidSystem<TypeTag, TTag::FVMPFAOVelocity2PTestTypeTag>
96template<
class TypeTag>
97struct FluidSystem<TypeTag, TTag::MimeticPressure2PTestTypeTag>
115template<
class TypeTag>
120 using Grid =
typename GridView::Grid;
122 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
124 using WettingPhase =
typename GET_PROP(TypeTag, FluidSystem)::WettingPhase;
126 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
130 dim = GridView::dimension, dimWorld = GridView::dimensionworld
135 wPhaseIdx = Indices::wPhaseIdx,
136 pwIdx = Indices::pwIdx,
137 swIdx = Indices::swIdx
142 using Element =
typename GridView::Traits::template Codim<0>::Entity;
143 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
144 using LocalPosition = Dune::FieldVector<Scalar, dim>;
147 using SolutionTypes =
typename GET_PROP(TypeTag, SolutionTypes);
148 using PrimaryVariables =
typename SolutionTypes::PrimaryVariables;
154 delta_ = getParam<Scalar>(
"Problem.Delta", 1e-3);
162 for (
int i = 0; i < this->
gridView().size(0); i++)
164 this->
variables().cellData(i).setSaturation(wPhaseIdx, 1.0);
166 this->
model().initialize();
190 for(
const auto& element : elements(this->
gridView()))
192 (*exactPressure)[this->
elementMapper().index(element)][0] =
exact(element.geometry().center());
221 void source(PrimaryVariables &values,
const Element& element)
const
225 values[wPhaseIdx] = integratedSource_(element, 4);
237 bcTypes.setAllDirichlet();
241 void dirichletAtPos(PrimaryVariables &values,
const GlobalPosition& globalPos)
const
243 values[pwIdx] =
exact(globalPos);
248 void neumannAtPos(PrimaryVariables &values,
const GlobalPosition& globalPos)
const
253 Scalar
exact (
const GlobalPosition& globalPos)
const
258 Scalar pi = 4.0*atan(1.0);
260 return (sin(pi*globalPos[0])*sin(pi*globalPos[1]));
263 Dune::FieldVector<Scalar,dim>
exactGrad (
const GlobalPosition& globalPos)
const
265 Dune::FieldVector<Scalar,dim> grad(0);
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]);
279 Scalar integratedSource_(
const Element& element,
int integrationPoints)
const
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++)
287 for (
int j = 1; j <= integrationPoints; j++)
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);
299 Scalar evaluateSource_(
const GlobalPosition& globalPos)
const
308 Scalar pi = 4.0 * atan(1.0);
309 Scalar x = globalPos[0];
310 Scalar y = globalPos[1];
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));
326 Scalar fx = dkxxdx * dpdx + kxx * dppdxx + dkxydx * dpdy + kxy * dppdyx;
327 Scalar fy = dkxydy * dpdx + kxy * dppdxy + dkyydy * dpdy + kyy * dppdyy;
Setting constant fluid properties via the input file.
#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
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
Defines the properties required for finite volume pressure models in a two-phase sequential model.
Defines the properties required for (immiscible) twophase sequential models.
Properties for the MPFA-O method.
Base class for stationary solution of a two-phase diffusion/pressure equation.
Finite volume velocity reconstruction.