3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokes/kovasznay/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 *****************************************************************************/
25#ifndef DUMUX_KOVASZNAY_TEST_PROBLEM_HH
26#define DUMUX_KOVASZNAY_TEST_PROBLEM_HH
27
28#ifndef UPWINDSCHEMEORDER
29#define UPWINDSCHEMEORDER 0
30#endif
31
32#include <dune/grid/yaspgrid.hh>
33
36
40#include "../l2error.hh"
41
42namespace Dumux {
43template <class TypeTag>
44class KovasznayTestProblem;
45
46namespace Properties {
47// Create new type tags
48namespace TTag {
49struct KovasznayTest { using InheritsFrom = std::tuple<NavierStokes, StaggeredFreeFlowModel>; };
50} // end namespace TTag
51
52// the fluid system
53template<class TypeTag>
54struct FluidSystem<TypeTag, TTag::KovasznayTest>
55{
58};
59
60// Set the grid type
61template<class TypeTag>
62struct Grid<TypeTag, TTag::KovasznayTest> { using type = Dune::YaspGrid<2, Dune::EquidistantOffsetCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
63
64// Set the problem property
65template<class TypeTag>
66struct Problem<TypeTag, TTag::KovasznayTest> { using type = Dumux::KovasznayTestProblem<TypeTag> ; };
67
68template<class TypeTag>
69struct EnableGridGeometryCache<TypeTag, TTag::KovasznayTest> { static constexpr bool value = true; };
70
71template<class TypeTag>
72struct EnableGridFluxVariablesCache<TypeTag, TTag::KovasznayTest> { static constexpr bool value = true; };
73template<class TypeTag>
74struct EnableGridVolumeVariablesCache<TypeTag, TTag::KovasznayTest> { static constexpr bool value = true; };
75
76template<class TypeTag>
77struct UpwindSchemeOrder<TypeTag, TTag::KovasznayTest> { static constexpr int value = UPWINDSCHEMEORDER; };
78} // end namespace Properties
79
88template <class TypeTag>
90{
91 using ParentType = NavierStokesProblem<TypeTag>;
92
95 using FVElementGeometry = typename GridGeometry::LocalView;
96 using SubControlVolume = typename GridGeometry::SubControlVolume;
103
104 static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
105 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
106 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
107 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
108
109 static constexpr auto upwindSchemeOrder = getPropValue<TypeTag, Properties::UpwindSchemeOrder>();
110
111public:
112 KovasznayTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
113 : ParentType(gridGeometry), eps_(1e-6)
114 {
115 printL2Error_ = getParam<bool>("Problem.PrintL2Error");
116 std::cout<< "upwindSchemeOrder is: " << GridGeometry::upwindStencilOrder() << "\n";
117 kinematicViscosity_ = getParam<Scalar>("Component.LiquidKinematicViscosity", 1.0);
118 Scalar reynoldsNumber = 1.0 / kinematicViscosity_;
119 lambda_ = 0.5 * reynoldsNumber
120 - std::sqrt(reynoldsNumber * reynoldsNumber * 0.25 + 4.0 * M_PI * M_PI);
121
122 createAnalyticalSolution_();
123 }
124
128 // \{
129
131 {
132 return false;
133 }
134
135 void printL2Error(const SolutionVector& curSol) const
136 {
137 if(printL2Error_)
138 {
140 const auto l2error = L2Error::calculateL2Error(*this, curSol);
141 const int numCellCenterDofs = this->gridGeometry().numCellCenterDofs();
142 const int numFaceDofs = this->gridGeometry().numFaceDofs();
143 std::cout << std::setprecision(8) << "** L2 error (abs/rel) for "
144 << std::setw(6) << numCellCenterDofs << " cc dofs and " << numFaceDofs << " face dofs (total: " << numCellCenterDofs + numFaceDofs << "): "
145 << std::scientific
146 << "L2(p) = " << l2error.first[Indices::pressureIdx] << " / " << l2error.second[Indices::pressureIdx]
147 << " , L2(vx) = " << l2error.first[Indices::velocityXIdx] << " / " << l2error.second[Indices::velocityXIdx]
148 << " , L2(vy) = " << l2error.first[Indices::velocityYIdx] << " / " << l2error.second[Indices::velocityYIdx]
149 << std::endl;
150 }
151 }
152
158 Scalar temperature() const
159 { return 298.0; }
160
161
167 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
168 {
169 return NumEqVector(0.0);
170 }
171
172 // \}
176 // \{
177
184 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
185 {
186 BoundaryTypes values;
187
188 // set Dirichlet values for the velocity everywhere
189 values.setDirichlet(Indices::velocityXIdx);
190 values.setDirichlet(Indices::velocityYIdx);
191
192 return values;
193 }
194
203 bool isDirichletCell(const Element& element,
204 const FVElementGeometry& fvGeometry,
205 const SubControlVolume& scv,
206 int pvIdx) const
207 {
208 // set fixed pressure in all cells at the left boundary
209 auto isAtLeftBoundary = [&](const FVElementGeometry& fvGeometry)
210 {
211 if (fvGeometry.hasBoundaryScvf())
212 {
213 for (const auto& scvf : scvfs(fvGeometry))
214 if (scvf.boundary() && scvf.center()[0] < this->gridGeometry().bBoxMin()[0] + eps_)
215 return true;
216 }
217 return false;
218 };
219 return (isAtLeftBoundary(fvGeometry) && pvIdx == Indices::pressureIdx);
220 }
221
227 PrimaryVariables dirichletAtPos(const GlobalPosition & globalPos) const
228 {
229 // use the values of the analytical solution
230 return analyticalSolution(globalPos);
231 }
232
238 PrimaryVariables analyticalSolution(const GlobalPosition& globalPos) const
239 {
240 Scalar x = globalPos[0];
241 Scalar y = globalPos[1];
242
243 PrimaryVariables values;
244 values[Indices::pressureIdx] = 0.5 * (1.0 - std::exp(2.0 * lambda_ * x));
245 values[Indices::velocityXIdx] = 1.0 - std::exp(lambda_ * x) * std::cos(2.0 * M_PI * y);
246 values[Indices::velocityYIdx] = 0.5 * lambda_ / M_PI * std::exp(lambda_ * x) * std::sin(2.0 * M_PI * y);
247
248 return values;
249 }
250
251 // \}
252
256 // \{
257
263 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
264 {
265 PrimaryVariables values;
266 values[Indices::pressureIdx] = 0.0;
267 values[Indices::velocityXIdx] = 0.0;
268 values[Indices::velocityYIdx] = 0.0;
269
270 return values;
271 }
272
277 {
278 return analyticalPressure_;
279 }
280
285 {
286 return analyticalVelocity_;
287 }
288
293 {
294 return analyticalVelocityOnFace_;
295 }
296
297private:
298
304 void createAnalyticalSolution_()
305 {
306 analyticalPressure_.resize(this->gridGeometry().numCellCenterDofs());
307 analyticalVelocity_.resize(this->gridGeometry().numCellCenterDofs());
308 analyticalVelocityOnFace_.resize(this->gridGeometry().numFaceDofs());
309
310 for (const auto& element : elements(this->gridGeometry().gridView()))
311 {
312 auto fvGeometry = localView(this->gridGeometry());
313 fvGeometry.bindElement(element);
314 for (auto&& scv : scvs(fvGeometry))
315 {
316 auto ccDofIdx = scv.dofIndex();
317 auto ccDofPosition = scv.dofPosition();
318 auto analyticalSolutionAtCc = analyticalSolution(ccDofPosition);
319
320 // velocities on faces
321 for (auto&& scvf : scvfs(fvGeometry))
322 {
323 const auto faceDofIdx = scvf.dofIndex();
324 const auto faceDofPosition = scvf.center();
325 const auto dirIdx = scvf.directionIndex();
326 const auto analyticalSolutionAtFace = analyticalSolution(faceDofPosition);
327 analyticalVelocityOnFace_[faceDofIdx][dirIdx] = analyticalSolutionAtFace[Indices::velocity(dirIdx)];
328 }
329
330 analyticalPressure_[ccDofIdx] = analyticalSolutionAtCc[Indices::pressureIdx];
331
332 for(int dirIdx = 0; dirIdx < ModelTraits::dim(); ++dirIdx)
333 analyticalVelocity_[ccDofIdx][dirIdx] = analyticalSolutionAtCc[Indices::velocity(dirIdx)];
334 }
335 }
336 }
337
338 Scalar eps_;
339
340 Scalar kinematicViscosity_;
341 Scalar lambda_;
342 bool printL2Error_;
343 std::vector<Scalar> analyticalPressure_;
344 std::vector<VelocityVector> analyticalVelocity_;
345 std::vector<VelocityVector> analyticalVelocityOnFace_;
346};
347} // end namespace Dumux
348
349#endif
Setting constant fluid properties via the input file.
A liquid phase consisting of a single component.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
Specifies the order of the upwinding scheme (1 == first order, 2 == second order(tvd methods))
Definition: common/properties.hh:305
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 (Kovasznay 1948, )
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:90
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/kovasznay/problem.hh:184
Scalar temperature() const
Returns the temperature within the domain in [K].
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:158
bool shouldWriteRestartFile() const
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:130
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Returns Dirichlet boundary values at a given position.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:227
void printL2Error(const SolutionVector &curSol) const
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:135
auto & getAnalyticalVelocitySolution() const
Returns the analytical solution for the velocity.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:284
bool isDirichletCell(const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, int pvIdx) const
Returns whether a fixed Dirichlet value shall be used at a given cell.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:203
auto & getAnalyticalPressureSolution() const
Returns the analytical solution for the pressure.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:276
KovasznayTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:112
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Returns the sources within the domain.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:167
auto & getAnalyticalVelocitySolutionOnFace() const
Returns the analytical solution for the velocity at the faces.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:292
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:263
PrimaryVariables analyticalSolution(const GlobalPosition &globalPos) const
Returns the analytical solution of the problem at a given position.
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:238
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:49
std::tuple< NavierStokes, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:49
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:56
Dune::YaspGrid< 2, Dune::EquidistantOffsetCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:62
Routines to calculate the discrete L2 error.
Definition: l2error.hh:38
Defines a type tag and some properties for ree-flow models using the staggered scheme....
A single-phase, isothermal Navier-Stokes model.
#define UPWINDSCHEMEORDER
Definition: test/freeflow/navierstokes/kovasznay/problem.hh:29