24#ifndef DUMUX_DAM_BREAK_TEST_PROBLEM_HH
25#define DUMUX_DAM_BREAK_TEST_PROBLEM_HH
27#include <dune/grid/yaspgrid.hh>
43template <
class TypeTag>
55template<
class TypeTag>
56struct Grid<TypeTag, TTag::DamBreakWet>
57{
using type = Dune::YaspGrid<2, Dune::TensorProductCoordinates<GetPropType<TypeTag, Properties::Scalar>, 2> >; };
60template<
class TypeTag>
65template<
class TypeTag>
75template<
class TypeTag>
77{
static constexpr bool value =
true; };
79template<
class TypeTag>
81{
static constexpr bool value =
false; };
83template<
class TypeTag>
85{
static constexpr bool value =
false; };
108template <
class TypeTag>
120 using ElementVolumeVariables =
typename GridVariables::GridVolumeVariables::LocalView;
121 using ElementFluxVariablesCache =
typename GridVariables::GridFluxVariablesCache::LocalView;
123 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
125 using Element =
typename GridView::template Codim<0>::Entity;
126 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
132 name_ = getParam<std::string>(
"Problem.Name");
140 return exactWaterDepth_;
146 return exactVelocityX_;
150 template<
class SolutionVector,
class Gr
idVariables>
152 const GridVariables& gridVariables,
156 for (
const auto& element : elements(this->
gridGeometry().gridView()))
159 fvGeometry.bindElement(element);
161 auto elemVolVars =
localView(gridVariables.curGridVolVars());
162 elemVolVars.bindElement(element, fvGeometry, curSol);
164 const auto& globalPos = element.geometry().center();
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);
181 const auto eIdx = this->
gridGeometry().elementMapper().index(element);
182 exactWaterDepth_[eIdx] = riemannResult.waterDepth;
183 exactVelocityX_[eIdx] = riemannResult.velocityX;
197 const std::string&
name()
const
216 BoundaryTypes bcTypes;
217 bcTypes.setAllNeumann();
230 const FVElementGeometry& fvGeometry,
231 const ElementVolumeVariables& elemVolVars,
232 const ElementFluxVariablesCache& elemFluxVarsCache,
233 const SubControlVolumeFace& scvf)
const
235 NeumannFluxes values(0.0);
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());
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(),
256 values[Indices::massBalanceIdx] = riemannFlux[0];
257 values[Indices::velocityXIdx] = riemannFlux[1];
258 values[Indices::velocityYIdx] = riemannFlux[2];
281 PrimaryVariables values(0.0);
283 values[0] = initialWaterDepthRight_;
288 if (globalPos[0] < 10.0 + eps_)
290 values[0] = initialWaterDepthLeft_;
296 values[0] = initialWaterDepthRight_;
306 std::vector<Scalar> exactWaterDepth_;
307 std::vector<Scalar> exactVelocityX_;
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;
314 static constexpr Scalar eps_ = 1.0e-6;
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.