3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
freeflow/shallowwater/localresidual.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_FREEFLOW_SHALLOW_WATER_LOCAL_RESIDUAL_HH
25#define DUMUX_FREEFLOW_SHALLOW_WATER_LOCAL_RESIDUAL_HH
26
30
31namespace Dumux{
32
37template<class TypeTag>
39: public GetPropType<TypeTag, Properties::BaseLocalResidual>
40{
46 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
47 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
48 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
49 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
50 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
51 using Element = typename GridView::template Codim<0>::Entity;
54
55public:
56
57 using ParentType::ParentType;
58
70 NumEqVector computeStorage(const Problem& problem,
71 const SubControlVolume& scv,
72 const VolumeVariables& volVars) const
73 {
74 // partial time derivative of the phase mass
75 NumEqVector storage(0.0);
76 storage[Indices::massBalanceIdx] = volVars.waterDepth();
77 storage[Indices::momentumXBalanceIdx] = volVars.waterDepth() * volVars.velocity(0);
78 storage[Indices::momentumYBalanceIdx] = volVars.waterDepth() * volVars.velocity(1);
79 return storage;
80 }
81
92 NumEqVector computeFlux(const Problem& problem,
93 const Element& element,
94 const FVElementGeometry& fvGeometry,
95 const ElementVolumeVariables& elemVolVars,
96 const SubControlVolumeFace& scvf,
97 const ElementFluxVariablesCache& elemFluxVarsCache) const
98 {
99 NumEqVector flux(0.0);
100 FluxVariables fluxVars;
101 flux += fluxVars.advectiveFlux(problem, element, fvGeometry, elemVolVars, scvf);
102
103 // Compute viscous momentum flux contribution if required
104 static const bool enableViscousFlux = getParamFromGroup<bool>(problem.paramGroup(), "ShallowWater.EnableViscousFlux", false);
105 if (enableViscousFlux)
106 flux += fluxVars.viscousFlux(problem, element, fvGeometry, elemVolVars, scvf);
107 return flux;
108 }
109};
110} // end namespace Dumux
111
112#endif
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
typename NumEqVectorTraits< PrimaryVariables >::type NumEqVector
A vector with the same size as numbers of equations This is the default implementation and has to be ...
Definition: numeqvector.hh:46
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Element-wise calculation of the residual for the shallow water equations.
Definition: freeflow/shallowwater/localresidual.hh:40
NumEqVector computeFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const ElementFluxVariablesCache &elemFluxVarsCache) const
Evaluate the mass/momentum flux over a face of a sub control volume.
Definition: freeflow/shallowwater/localresidual.hh:92
NumEqVector computeStorage(const Problem &problem, const SubControlVolume &scv, const VolumeVariables &volVars) const
Evaluate the rate of change of all conservation quantites (e.g. mass, momentum) within a sub-control ...
Definition: freeflow/shallowwater/localresidual.hh:70
Declares all properties used in Dumux.