3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
test/freeflow/shallowwater/dambreak/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 2 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_DAM_BREAK_TEST_PROBLEM_HH
25#define DUMUX_DAM_BREAK_TEST_PROBLEM_HH
26
27#include <dune/grid/yaspgrid.hh>
29#include "spatialparams.hh"
30
35
36
37namespace Dumux {
38
43template <class TypeTag>
44class DamBreakProblem;
45
46
47// Specify the properties for the problem
48namespace Properties {
49
50// Create new type tags
51namespace TTag {
52struct DamBreakWet { using InheritsFrom = std::tuple<ShallowWater, CCTpfaModel>; };
53} // end namespace TTag
54
55template<class TypeTag>
56struct Grid<TypeTag, TTag::DamBreakWet>
57{ using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
58
59// Set the problem property
60template<class TypeTag>
61struct Problem<TypeTag, TTag::DamBreakWet>
63
64// Set the spatial parameters
65template<class TypeTag>
66struct SpatialParams<TypeTag, TTag::DamBreakWet>
67{
68private:
71public:
73};
74
75template<class TypeTag>
76struct EnableGridGeometryCache<TypeTag, TTag::DamBreakWet>
77{ static constexpr bool value = true; };
78
79template<class TypeTag>
80struct EnableGridVolumeVariablesCache<TypeTag, TTag::DamBreakWet>
81{ static constexpr bool value = false; };
82
83template<class TypeTag>
84struct EnableGridFluxVariablesCache<TypeTag, TTag::DamBreakWet>
85{ static constexpr bool value = false; };
86
87} // end namespace Properties
88
89
108template <class TypeTag>
110{
118
120 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
121 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
122 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
123 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
125 using Element = typename GridView::template Codim<0>::Entity;
126 using GlobalPosition = typename Element::Geometry::GlobalCoordinate;
127
128public:
129 DamBreakProblem(std::shared_ptr<const GridGeometry> gridGeometry)
131 {
132 name_ = getParam<std::string>("Problem.Name");
133 exactWaterDepth_.resize(gridGeometry->numDofs(), 0.0);
134 exactVelocityX_.resize(gridGeometry->numDofs(), 0.0);
135 }
136
138 const std::vector<Scalar>& getExactWaterDepth()
139 {
140 return exactWaterDepth_;
141 }
142
144 const std::vector<Scalar>& getExactVelocityX()
145 {
146 return exactVelocityX_;
147 }
148
150 template<class SolutionVector, class GridVariables>
151 void updateAnalyticalSolution(const SolutionVector& curSol,
152 const GridVariables& gridVariables,
153 const Scalar time)
154 {
155 //compute solution for all elements
156 for (const auto& element : elements(this->gridGeometry().gridView()))
157 {
158 auto fvGeometry = localView(this->gridGeometry());
159 fvGeometry.bindElement(element);
160
161 auto elemVolVars = localView(gridVariables.curGridVolVars());
162 elemVolVars.bindElement(element, fvGeometry, curSol);
163
164 const auto& globalPos = element.geometry().center();
165
166 //compute the position s and the initial water depth at the gate, velocities are zero
167 const Scalar s = (globalPos[0] - gatePosition_)/time;
168 const Scalar waterDepthLeft = initialWaterDepthLeft_;
169 const Scalar waterDepthRight = initialWaterDepthRight_;
170 const auto gravity = this->spatialParams().gravity(globalPos);
171
172 auto riemannResult = ShallowWater::exactRiemann(waterDepthLeft,
173 waterDepthRight,
174 0.0,
175 0.0,
176 0.0,
177 0.0,
178 gravity,
179 s);
180
181 const auto eIdx = this->gridGeometry().elementMapper().index(element);
182 exactWaterDepth_[eIdx] = riemannResult.waterDepth;
183 exactVelocityX_[eIdx] = riemannResult.velocityX;
184 }
185 }
186
190 // \{
191
197 const std::string& name() const
198 { return name_; }
199
200
201 // \}
202
206 // \{
207
214 BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
215 {
216 BoundaryTypes bcTypes;
217 bcTypes.setAllNeumann();
218 return bcTypes;
219 }
220
229 NeumannFluxes neumann(const Element& element,
230 const FVElementGeometry& fvGeometry,
231 const ElementVolumeVariables& elemVolVars,
232 const ElementFluxVariablesCache& elemFluxVarsCache,
233 const SubControlVolumeFace& scvf) const
234 {
235 NeumannFluxes values(0.0);
236
237 // we need the Riemann invariants to compute the values depending of the boundary type
238 // since we use a weak imposition we do not have a dirichlet value. We impose fluxes
239 // based on q,h, etc. computed with the Riemann invariants
240 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
241 const auto& insideVolVars = elemVolVars[insideScv];
242 const auto& nxy = scvf.unitOuterNormal();
243 const auto gravity = this->spatialParams().gravity(scvf.center());
244
245 auto riemannFlux = ShallowWater::riemannProblem(insideVolVars.waterDepth(),
246 insideVolVars.waterDepth(),
247 insideVolVars.velocity(0),
248 -insideVolVars.velocity(0),
249 insideVolVars.velocity(1),
250 -insideVolVars.velocity(1),
251 insideVolVars.bedSurface(),
252 insideVolVars.bedSurface(),
253 gravity,
254 nxy);
255
256 values[Indices::massBalanceIdx] = riemannFlux[0];
257 values[Indices::velocityXIdx] = riemannFlux[1];
258 values[Indices::velocityYIdx] = riemannFlux[2];
259
260 return values;
261 }
262
263 // \}
264
268 // \{
269
278 PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
279 {
280
281 PrimaryVariables values(0.0);
282
283 values[0] = initialWaterDepthRight_;
284 values[1] = 0.0;
285 values[2] = 0.0;
286
287 // water level on the left side of the gate
288 if (globalPos[0] < 10.0 + eps_)
289 {
290 values[0] = initialWaterDepthLeft_;
291 }
292
293 //water level on the right side of the gate
294 else
295 {
296 values[0] = initialWaterDepthRight_;
297 }
298
299 return values;
300 };
301
302 // \}
303
304private:
305
306 std::vector<Scalar> exactWaterDepth_;
307 std::vector<Scalar> exactVelocityX_;
308
309 static constexpr Scalar initialWaterDepthLeft_ = 4.0;
310 static constexpr Scalar initialWaterDepthRight_ = 1.0;
311 static constexpr Scalar channelLenght_ = 20.0;
312 static constexpr Scalar gatePosition_ = 10.0;
313
314 static constexpr Scalar eps_ = 1.0e-6;
315 std::string name_;
316};
317
318} //end namespace Dumux
319
320#endif
Properties for all models using cell-centered finite volume scheme with TPFA.
Function to compute the Riemann flux at the interface.
This file includes a function to construct the Riemann problem.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
RiemannSolution< Scalar > exactRiemann(const Scalar dl, const Scalar dr, const Scalar ul, const Scalar ur, const Scalar vl, const Scalar vr, const Scalar grav, const Scalar s=0.0)
Exact Riemann solver for Shallow water equations.
Definition: exactriemann.hh:61
std::array< Scalar, 3 > riemannProblem(const Scalar waterDepthLeft, const Scalar waterDepthRight, Scalar velocityXLeft, Scalar velocityXRight, Scalar velocityYLeft, Scalar velocityYRight, const Scalar bedSurfaceLeft, const Scalar bedSurfaceRight, const Scalar gravity, const GlobalPosition &nxy)
Construct Riemann Problem and solve it.
Definition: riemannproblem.hh:66
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
Base class for all finite-volume problems.
Definition: common/fvproblem.hh:50
const GridGeometry & gridGeometry() const
The finite volume grid geometry.
Definition: common/fvproblem.hh:588
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 spatial parameters object.
Definition: common/properties.hh:221
Shallow water problem base class.
Definition: dumux/freeflow/shallowwater/problem.hh:39
const SpatialParams & spatialParams() const
Returns the spatial parameters object.
Definition: dumux/freeflow/shallowwater/problem.hh:78
A simple dam break test for the shallow water equations.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:110
DamBreakProblem(std::shared_ptr< const GridGeometry > gridGeometry)
Definition: test/freeflow/shallowwater/dambreak/problem.hh:129
NeumannFluxes neumann(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Specifies the neumann bounday.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:229
const std::vector< Scalar > & getExactVelocityX()
Get the analytical velocity.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:144
const std::vector< Scalar > & getExactWaterDepth()
Get the analytical water depth.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:138
BoundaryTypes boundaryTypesAtPos(const GlobalPosition &globalPos) const
Specifies which kind of boundary condition should be used for which equation on a given boundary segm...
Definition: test/freeflow/shallowwater/dambreak/problem.hh:214
PrimaryVariables initialAtPos(const GlobalPosition &globalPos) const
Evaluate the initial values for a control volume.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:278
const std::string & name() const
The problem name.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:197
void updateAnalyticalSolution(const SolutionVector &curSol, const GridVariables &gridVariables, const Scalar time)
Udpate the analytical solution.
Definition: test/freeflow/shallowwater/dambreak/problem.hh:151
Definition: test/freeflow/shallowwater/dambreak/problem.hh:52
std::tuple< ShallowWater, CCTpfaModel > InheritsFrom
Definition: test/freeflow/shallowwater/dambreak/problem.hh:52
Dune::YaspGrid< 2, Dune::TensorProductCoordinates< GetPropType< TypeTag, Properties::Scalar >, 2 > > type
Definition: test/freeflow/shallowwater/dambreak/problem.hh:57
The spatial parameters class for the dam break test.
Definition: freeflow/shallowwater/dambreak/spatialparams.hh:40
A two-dimensional shallow water equations model The two-dimensonal shallow water equations (SWEs) can...
Definition of the spatial parameters for the MaxwellStefan problem.