22#ifndef DUMUX_TEST_DIFFUSION_3D_PROBLEM_HH
23#define DUMUX_TEST_DIFFUSION_3D_PROBLEM_HH
26#include <dune/alugrid/grid.hh>
29#include <dune/grid/uggrid.hh>
31#include <dune/grid/yaspgrid.hh>
49template<
class TypeTag>
61SET_TYPE_PROP(DiffusionTest,
Grid, Dune::ALUGrid<3, 3, Dune::cube, Dune::nonconforming>);
71template<
class TypeTag>
101template<
class TypeTag>
104 using ParentType = DiffusionProblem2P<TypeTag>;
106 using Grid =
typename GridView::Grid;
108 using Indices =
typename GET_PROP_TYPE(TypeTag, ModelTraits)::Indices;
110 using FluidSystem =
typename GET_PROP_TYPE(TypeTag, FluidSystem);
111 using FluidState =
typename GET_PROP_TYPE(TypeTag, FluidState);
113 using BoundaryTypes =
typename GET_PROP_TYPE(TypeTag, BoundaryTypes);
117 dim = GridView::dimension, dimWorld = GridView::dimensionworld
122 wPhaseIdx = Indices::wPhaseIdx,
123 nPhaseIdx = Indices::nPhaseIdx,
124 pWIdx = Indices::pressureIdx,
125 swIdx = Indices::swIdx,
126 pressureEqIdx = Indices::pressureEqIdx,
131 using Element =
typename GridView::Traits::template Codim<0>::Entity;
132 using Intersection =
typename GridView::Intersection;
133 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
141 ParentType(grid), velocity_(*this)
148 for (
int i = 0; i < this->
gridView().size(0); i++)
150 this->
variables().cellData(i).setSaturation(wPhaseIdx, 1.0);
151 this->
variables().cellData(i).setSaturation(nPhaseIdx, 0.0);
153 this->
model().initialize();
158 velocity_.calculateVelocity();
167 for(
const auto& element : elements(this->
gridView()))
169 (*exactPressure)[this->
elementMapper().index(element)][0] =
exact(element.geometry().center());
211 double pi = 4.0*atan(1.0);
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));
226 bcTypes.setAllDirichlet();
232 values[pWIdx] =
exact(globalPos);
242 Scalar
exact (
const GlobalPosition& globalPos)
const
247 double pi = 4.0*atan(1.0);
249 return (1.0+sin(pi*globalPos[0])*sin(pi*(globalPos[1]+1.0/2.0))*sin(pi*(globalPos[2]+1.0/3.0)));
252 Dune::FieldVector<Scalar,dim>
exactGrad (
const GlobalPosition& globalPos)
const
254 Dune::FieldVector<Scalar,dim> grad(0);
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));
269 static constexpr Scalar eps_ = 1e-4;
#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.
#define INHERITS_FROM(...)
Syntactic sugar for NEW_TYPE_TAG.
Definition propertysystemmacros.hh:142
#define SET_TYPE_PROP(EffTypeTagName, PropTagName,...)
Set a property which defines a type.
Definition propertysystemmacros.hh:232
make the local view function available whenever we use the grid geometry
Definition adapt.hh:29
Definition common/properties.hh:47
Type tag TestDiffusionSpatialParams3d
Definition test_diffusionspatialparams3d.hh:41
Type tag for numeric models.
Definition grid.hh:35
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 finite volume velocity reconstruction.
Definition sequential/cellcentered/velocity.hh:48
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
typename SolutionTypes::PrimaryVariables PrimaryVariables
Definition test_diffusionproblem3d.hh:137
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
typename GET_PROP(TypeTag, SolutionTypes) SolutionTypes
Definition test_diffusionproblem3d.hh:136
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
FluidSystems::OnePLiquid< Scalar, Components::Constant< 1, Scalar > > NonwettingPhase
Definition test_diffusionproblem3d.hh:76
FluidSystems::OnePLiquid< Scalar, Components::Constant< 1, Scalar > > WettingPhase
Definition test_diffusionproblem3d.hh:75
typename GET_PROP_TYPE(TypeTag, Scalar) Scalar
Definition test_diffusionproblem3d.hh:74
FluidSystems::TwoPImmiscible< Scalar, WettingPhase, NonwettingPhase > type
Definition test_diffusionproblem3d.hh:77
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.