3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/navierstokes/donea/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 *****************************************************************************/
24#ifndef DUMUX_DONEA_TEST_PROBLEM_HH
25#define DUMUX_DONEA_TEST_PROBLEM_HH
26
27#ifndef ENABLECACHING
28#define ENABLECACHING 0
29#endif
30
31#include <dune/grid/yaspgrid.hh>
32
35
39#include "../l2error.hh"
40
41
42namespace Dumux
43{
44template <class TypeTag>
45class DoneaTestProblem;
46
47namespace Properties
48{
49// Create new type tags
50namespace TTag {
51struct DoneaTest { using InheritsFrom = std::tuple<NavierStokes, StaggeredFreeFlowModel>; };
52} // end namespace TTag
53
54// the fluid system
55template<class TypeTag>
56struct FluidSystem<TypeTag, TTag::DoneaTest>
57{
60};
61
62// Set the grid type
63template<class TypeTag>
64struct Grid<TypeTag, TTag::DoneaTest> { using type = Dune::YaspGrid<2>; };
65
66// Set the problem property
67template<class TypeTag>
68struct Problem<TypeTag, TTag::DoneaTest> { using type = Dumux::DoneaTestProblem<TypeTag> ; };
69
70template<class TypeTag>
71struct EnableGridGeometryCache<TypeTag, TTag::DoneaTest> { static constexpr bool value = ENABLECACHING; };
72template<class TypeTag>
73struct EnableGridFluxVariablesCache<TypeTag, TTag::DoneaTest> { static constexpr bool value = ENABLECACHING; };
74template<class TypeTag>
75struct EnableGridVolumeVariablesCache<TypeTag, TTag::DoneaTest> { static constexpr bool value = ENABLECACHING; };
76template<class TypeTag>
77struct EnableGridFaceVariablesCache<TypeTag, TTag::DoneaTest> { static constexpr bool value = ENABLECACHING; };
78}
79
88template <class TypeTag>
90{
91 using ParentType = NavierStokesProblem<TypeTag>;
92
101
102 static constexpr auto dimWorld = GetPropType<TypeTag, Properties::GridView>::dimensionworld;
103 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
104 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
105 using VelocityVector = Dune::FieldVector<Scalar, dimWorld>;
106
107public:
108 DoneaTestProblem(std::shared_ptr<const GridGeometry> gridGeometry)
109 : ParentType(gridGeometry), eps_(1e-6)
110 {
111 printL2Error_ = getParam<bool>("Problem.PrintL2Error");
112 createAnalyticalSolution_();
113 }
114
118 // \{
119
121 {
122 return false;
123 }
124
125 void printL2Error(const SolutionVector& curSol) const
126 {
127 if(printL2Error_)
128 {
130 const auto l2error = L2Error::calculateL2Error(*this, curSol);
131 const int numCellCenterDofs = this->gridGeometry().numCellCenterDofs();
132 const int numFaceDofs = this->gridGeometry().numFaceDofs();
133 std::cout << std::setprecision(8) << "** L2 error (abs/rel) for "
134 << std::setw(6) << numCellCenterDofs << " cc dofs and " << numFaceDofs << " face dofs (total: " << numCellCenterDofs + numFaceDofs << "): "
135 << std::scientific
136 << "L2(p) = " << l2error.first[Indices::pressureIdx] << " / " << l2error.second[Indices::pressureIdx]
137 << " , L2(vx) = " << l2error.first[Indices::velocityXIdx] << " / " << l2error.second[Indices::velocityXIdx]
138 << " , L2(vy) = " << l2error.first[Indices::velocityYIdx] << " / " << l2error.second[Indices::velocityYIdx]
139 << std::endl;
140 }
141 }
142
148 Scalar temperature() const
149 { return 298.0; }
150
156 NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
157 {
158 NumEqVector source(0.0);
159 Scalar x = globalPos[0];
160 Scalar y = globalPos[1];
161
162 source[Indices::momentumXBalanceIdx] = (12.0-24.0*y) * x*x*x*x + (-24.0 + 48.0*y)* x*x*x
163 + (-48.0*y + 72.0*y*y - 48.0*y*y*y + 12.0)* x*x
164 + (-2.0 + 24.0*y - 72.0*y*y + 48.0*y*y*y)*x
165 + 1.0 - 4.0*y + 12.0*y*y - 8.0*y*y*y;
166 source[Indices::momentumYBalanceIdx] = (8.0 - 48.0*y + 48.0*y*y)*x*x*x + (-12.0 + 72.0*y - 72.0*y*y)*x*x
167 + (4.0 - 24.0*y + 48.0*y*y - 48.0*y*y*y + 24.0*y*y*y*y)*x - 12.0*y*y
168 + 24.0*y*y*y - 12.0*y*y*y*y;
169 return source;
170 }
171 // \}
175 // \{
176
183 BoundaryTypes boundaryTypesAtPos(const GlobalPosition& globalPos) const
184 {
185 BoundaryTypes values;
186
187 // set Dirichlet values for the velocity and pressure everywhere
188 values.setDirichlet(Indices::velocityXIdx);
189 values.setDirichlet(Indices::velocityYIdx);
190
191 return values;
192 }
193
202 bool isDirichletCell(const Element& element,
203 const typename GridGeometry::LocalView& fvGeometry,
204 const typename GridGeometry::SubControlVolume& scv,
205 int pvIdx) const
206 {
207 bool onBoundary = false;
208 for (const auto& scvf : scvfs(fvGeometry))
209 onBoundary = std::max(onBoundary, scvf.boundary());
210 return onBoundary;
211 }
212
218 PrimaryVariables dirichletAtPos(const GlobalPosition& globalPos) const
219 {
220 // use the values of the analytical solution
221 return analyticalSolution(globalPos);
222 }
223
229 PrimaryVariables analyticalSolution(const GlobalPosition& globalPos) const
230 {
231 Scalar x = globalPos[0];
232 Scalar y = globalPos[1];
233
234 PrimaryVariables values;
235 values[Indices::pressureIdx] = x * (1.0-x); // p(x,y) = x(1-x) [Donea2003]
236 values[Indices::velocityXIdx] = x*x * (1.0 - x)*(1.0 - x) * (2.0*y - 6.0*y*y + 4.0*y*y*y);
237 values[Indices::velocityYIdx] = -1.0*y*y * (1.0 - y)*(1.0 - y) * (2.0*x - 6.0*x*x + 4.0*x*x*x);
238
239 return values;
240 }
241
242 // \}
243
247 // \{
248
254 PrimaryVariables initialAtPos(const GlobalPosition& globalPos) const
255 {
256 PrimaryVariables values;
257 values[Indices::pressureIdx] = 0.0;
258 values[Indices::velocityXIdx] = 0.0;
259 values[Indices::velocityYIdx] = 0.0;
260
261 return values;
262 }
263
268 {
269 return analyticalPressure_;
270 }
271
276 {
277 return analyticalVelocity_;
278 }
279
284 {
285 return analyticalVelocityOnFace_;
286 }
287
288private:
289
293 void createAnalyticalSolution_()
294 {
295 analyticalPressure_.resize(this->gridGeometry().numCellCenterDofs());
296 analyticalVelocity_.resize(this->gridGeometry().numCellCenterDofs());
297 analyticalVelocityOnFace_.resize(this->gridGeometry().numFaceDofs());
298
299 for (const auto& element : elements(this->gridGeometry().gridView()))
300 {
301 auto fvGeometry = localView(this->gridGeometry());
302 fvGeometry.bindElement(element);
303 for (auto&& scv : scvs(fvGeometry))
304 {
305 auto ccDofIdx = scv.dofIndex();
306 auto ccDofPosition = scv.dofPosition();
307 auto analyticalSolutionAtCc = analyticalSolution(ccDofPosition);
308
309 // velocities on faces
310 for (auto&& scvf : scvfs(fvGeometry))
311 {
312 const auto faceDofIdx = scvf.dofIndex();
313 const auto faceDofPosition = scvf.center();
314 const auto dirIdx = scvf.directionIndex();
315 const auto analyticalSolutionAtFace = analyticalSolution(faceDofPosition);
316 analyticalVelocityOnFace_[faceDofIdx][dirIdx] = analyticalSolutionAtFace[Indices::velocity(dirIdx)];
317 }
318
319 analyticalPressure_[ccDofIdx] = analyticalSolutionAtCc[Indices::pressureIdx];
320
321 for(int dirIdx = 0; dirIdx < ModelTraits::dim(); ++dirIdx)
322 analyticalVelocity_[ccDofIdx][dirIdx] = analyticalSolutionAtCc[Indices::velocity(dirIdx)];
323 }
324 }
325 }
326
327 Scalar eps_;
328 bool printL2Error_;
329 std::vector<Scalar> analyticalPressure_;
330 std::vector<VelocityVector> analyticalVelocity_;
331 std::vector<VelocityVector> analyticalVelocityOnFace_;
332};
333} // end namespace Dumux
334
335#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
Switch on/off caching of face variables.
Definition: common/properties.hh:303
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 (Donea 2003, ).
Definition: test/freeflow/navierstokes/donea/problem.hh:90
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluates the initial value for a control volume.
Definition: test/freeflow/navierstokes/donea/problem.hh:254
PrimaryVariables analyticalSolution(const GlobalPosition &globalPos) const
Return the analytical solution of the problem at a given position.
Definition: test/freeflow/navierstokes/donea/problem.hh:229
Scalar temperature() const
Return the temperature within the domain in [K].
Definition: test/freeflow/navierstokes/donea/problem.hh:148
auto & getAnalyticalVelocitySolutionOnFace() const
Returns the analytical solution for the velocity at the faces.
Definition: test/freeflow/navierstokes/donea/problem.hh:283
void printL2Error(const SolutionVector &curSol) const
Definition: test/freeflow/navierstokes/donea/problem.hh:125
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/donea/problem.hh:202
auto & getAnalyticalVelocitySolution() const
Returns the analytical solution for the velocity.
Definition: test/freeflow/navierstokes/donea/problem.hh:275
bool shouldWriteRestartFile() const
Definition: test/freeflow/navierstokes/donea/problem.hh:120
auto & getAnalyticalPressureSolution() const
Returns the analytical solution for the pressure.
Definition: test/freeflow/navierstokes/donea/problem.hh:267
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/donea/problem.hh:183
PrimaryVariables dirichletAtPos(const GlobalPosition &globalPos) const
Return dirichlet boundary values at a given position.
Definition: test/freeflow/navierstokes/donea/problem.hh:218
DoneaTestProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/navierstokes/donea/problem.hh:108
NumEqVector sourceAtPos(const GlobalPosition &globalPos) const
Return the sources within the domain.
Definition: test/freeflow/navierstokes/donea/problem.hh:156
Definition: test/freeflow/navierstokes/donea/problem.hh:51
std::tuple< NavierStokes, StaggeredFreeFlowModel > InheritsFrom
Definition: test/freeflow/navierstokes/donea/problem.hh:51
GetPropType< TypeTag, Properties::Scalar > Scalar
Definition: test/freeflow/navierstokes/donea/problem.hh:58
Dune::YaspGrid< 2 > type
Definition: test/freeflow/navierstokes/donea/problem.hh:64
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 ENABLECACHING
Definition: test/freeflow/navierstokes/donea/problem.hh:28