26#ifndef DUMUX_ANGELI_TEST_PROBLEM_HH
27#define DUMUX_ANGELI_TEST_PROBLEM_HH
29#include <dune/grid/yaspgrid.hh>
37#include "../l2error.hh"
40template <
class TypeTag>
41class AngeliTestProblem;
50template<
class TypeTag>
60template<
class TypeTag>
61struct Grid<TypeTag, TTag::AngeliTest> {
using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
64template<
class TypeTag>
67template<
class TypeTag>
70template<
class TypeTag>
72template<
class TypeTag>
85template <
class TypeTag>
97 using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;
100 using Element =
typename GridGeometry::GridView::template Codim<0>::Entity;
101 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
102 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
108 : ParentType(gridGeometry)
110 kinematicViscosity_ = getParam<Scalar>(
"Component.LiquidKinematicViscosity", 1.0);
135 BoundaryTypes values;
138 values.setDirichlet(Indices::velocityXIdx);
139 values.setDirichlet(Indices::velocityYIdx);
153 const typename GridGeometry::LocalView& fvGeometry,
154 const typename GridGeometry::SubControlVolume& scv,
158 for (
const auto& scvf : scvfs(fvGeometry))
184 const Scalar x = globalPos[0];
185 const Scalar y = globalPos[1];
187 const Scalar t = time;
189 PrimaryVariables values;
191 values[Indices::pressureIdx] = - 0.25 * std::exp(-10.0 * kinematicViscosity_ * M_PI * M_PI * t) * M_PI * M_PI * (4.0 * std::cos(2.0 * M_PI * x) + std::cos(4.0 * M_PI * y));
192 values[Indices::velocityXIdx] = - 2.0 * M_PI * std::exp(- 5.0 * kinematicViscosity_ * M_PI * M_PI * t) * std::cos(M_PI * x) * std::sin(2.0 * M_PI * y);
193 values[Indices::velocityYIdx] = M_PI * std::exp(- 5.0 * kinematicViscosity_ * M_PI * M_PI * t) * std::sin(M_PI * x) * std::cos(2.0 * M_PI * y);
228 timeStepSize_ = timeStepSize;
232 static constexpr Scalar eps_ = 1e-6;
234 Scalar kinematicViscosity_;
236 Scalar timeStepSize_ = 0;
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
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:169
If disabled, the volume variables are not stored (reduces memory, but is slower)
Definition: common/properties.hh:178
specifies if data on flux vars should be saved (faster, but more memory consuming)
Definition: common/properties.hh:188
The type of the fluid system to use.
Definition: common/properties.hh:223
Navier-Stokes problem base class.
Definition: dumux/freeflow/navierstokes/problem.hh:63
A liquid phase consisting of a single component.
Definition: 1pliquid.hh:46
Test problem for the staggered grid (Angeli et al. 2017, ).
Definition: test/freeflow/navierstokes/angeli/problem.hh:87
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokes/angeli/problem.hh:210
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/freeflow/navierstokes/angeli/problem.hh:133
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Returns Dirichlet boundary values at a given position.
Definition: test/freeflow/navierstokes/angeli/problem.hh:170
bool isDirichletCell(const Element &element, const typename GridGeometry::LocalView &fvGeometry, const typename GridGeometry::SubControlVolume &scv, int pvIdx) const
Returns whether a fixed Dirichlet value shall be used at a given cell.
Definition: test/freeflow/navierstokes/angeli/problem.hh:152
void updateTime(const Scalar time)
Updates the time.
Definition: test/freeflow/navierstokes/angeli/problem.hh:218
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/navierstokes/angeli/problem.hh:118
void updateTimeStepSize(const Scalar timeStepSize)
Updates the time step size.
Definition: test/freeflow/navierstokes/angeli/problem.hh:226
AngeliTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokes/angeli/problem.hh:107
PrimaryVariables analyticalSolution(const GlobalPosition &globalPos, const Scalar time) const
Returns the analytical solution of the problem at a given position.
Definition: test/freeflow/navierstokes/angeli/problem.hh:182
Definition: test/freeflow/navierstokes/angeli/problem.hh:46
std::tuple< NavierStokes, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokes/angeli/problem.hh:46
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: test/freeflow/navierstokes/angeli/problem.hh:61
Defines a type tag and some properties for ree-flow models using the staggered scheme....
A single-phase, isothermal Navier-Stokes model.