3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokes/angeli/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 *****************************************************************************/
26#ifndef DUMUX_ANGELI_TEST_PROBLEM_HH
27#define DUMUX_ANGELI_TEST_PROBLEM_HH
28
29#include <dune/grid/yaspgrid.hh>
30
33
37#include "../l2error.hh"
38
39namespace Dumux {
40template <class TypeTag>
41class AngeliTestProblem;
42
43namespace Properties {
44// Create new type tags
45namespace TTag {
46struct AngeliTest { using InheritsFrom = std::tuple<NavierStokes, StaggeredFreeFlowModel>; };
47} // end namespace TTag
48
49// the fluid system
50template<class TypeTag>
51struct FluidSystem<TypeTag, TTag::AngeliTest>
52{
53private:
55public:
57};
58
59// Set the grid type
60template<class TypeTag>
61struct Grid<TypeTag, TTag::AngeliTest> { using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
62
63// Set the problem property
64template<class TypeTag>
65struct Problem<TypeTag, TTag::AngeliTest> { using type = Dumux::AngeliTestProblem<TypeTag> ; };
66
67template<class TypeTag>
68struct EnableGridGeometryCache<TypeTag, TTag::AngeliTest> { static constexpr bool value = true; };
69
70template<class TypeTag>
71struct EnableGridFluxVariablesCache<TypeTag, TTag::AngeliTest> { static constexpr bool value = true; };
72template<class TypeTag>
73struct EnableGridVolumeVariablesCache<TypeTag, TTag::AngeliTest> { static constexpr bool value = true; };
74} // end namespace Properties
75
85template <class TypeTag>
87{
88 using ParentType = NavierStokesProblem<TypeTag>;
89
97 using TimeLoopPtr = std::shared_ptr<TimeLoop<Scalar>>;
98
100 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
101 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
102 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
103
104public:
106
107 AngeliTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
108 : ParentType(gridGeometry)
109 {
110 kinematicViscosity_ = getParam<Scalar>("Component.LiquidKinematicViscosity", 1.0);
111 }
112
118 Scalar temperature() const
119 { return 298.0; }
120
121 // \}
125 // \{
126
133 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
134 {
135 BoundaryTypes values;
136
137 // set Dirichlet values for the velocity everywhere
138 values.setDirichlet(Indices::velocityXIdx);
139 values.setDirichlet(Indices::velocityYIdx);
140
141 return values;
142 }
143
152 bool isDirichletCell(const Element& element,
153 const typename GridGeometry::LocalView& fvGeometry,
154 const typename GridGeometry::SubControlVolume& scv,
155 int pvIdx) const
156 {
157 // set a fixed pressure in all cells at the boundary
158 for (const auto& scvf : scvfs(fvGeometry))
159 if (scvf.boundary())
160 return true;
161
162 return false;
163 }
164
170 PrimaryVariables dirichletAtPos(const GlobalPosition & globalPos) const
171 {
172 // use the values of the analytical solution
173 return analyticalSolution(globalPos, time_+timeStepSize_);
174 }
175
182 PrimaryVariables analyticalSolution(const GlobalPosition& globalPos, const Scalar time) const
183 {
184 const Scalar x = globalPos[0];
185 const Scalar y = globalPos[1];
186
187 const Scalar t = time;
188
189 PrimaryVariables values;
190
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);
194
195 return values;
196 }
197
198 // \}
199
203 // \{
204
210 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
211 {
212 return analyticalSolution(globalPos, 0);
213 }
214
218 void updateTime(const Scalar time)
219 {
220 time_ = time;
221 }
222
226 void updateTimeStepSize(const Scalar timeStepSize)
227 {
228 timeStepSize_ = timeStepSize;
229 }
230
231private:
232 static constexpr Scalar eps_ = 1e-6;
233
234 Scalar kinematicViscosity_;
235 Scalar time_ = 0;
236 Scalar timeStepSize_ = 0;
237};
238} // end namespace Dumux
239
240#endif
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.