version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_FREEFLOW_SHALLOW_WATER_LOCAL_RESIDUAL_HH
13#define DUMUX_FREEFLOW_SHALLOW_WATER_LOCAL_RESIDUAL_HH
14
18
19namespace Dumux{
20
25template<class TypeTag>
27: public GetPropType<TypeTag, Properties::BaseLocalResidual>
28{
34 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
35 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
36 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
37 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
38 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
39 using Element = typename GridView::template Codim<0>::Entity;
42
43public:
44
45 using ParentType::ParentType;
46
58 NumEqVector computeStorage(const Problem& problem,
59 const SubControlVolume& scv,
60 const VolumeVariables& volVars) const
61 {
62 // partial time derivative of the phase mass
63 NumEqVector storage(0.0);
64 storage[Indices::massBalanceIdx] = volVars.waterDepth();
65 storage[Indices::momentumXBalanceIdx] = volVars.waterDepth() * volVars.velocity(0);
66 storage[Indices::momentumYBalanceIdx] = volVars.waterDepth() * volVars.velocity(1);
67 return storage;
68 }
69
80 NumEqVector computeFlux(const Problem& problem,
81 const Element& element,
82 const FVElementGeometry& fvGeometry,
83 const ElementVolumeVariables& elemVolVars,
84 const SubControlVolumeFace& scvf,
85 const ElementFluxVariablesCache& elemFluxVarsCache) const
86 {
87 NumEqVector flux(0.0);
88 FluxVariables fluxVars;
89 flux += fluxVars.advectiveFlux(problem, element, fvGeometry, elemVolVars, scvf);
90
91 // Compute viscous momentum flux contribution if required
92 static const bool enableViscousFlux = getParamFromGroup<bool>(problem.paramGroup(), "ShallowWater.EnableViscousFlux", false);
93 if (enableViscousFlux)
94 flux += fluxVars.viscousFlux(problem, element, fvGeometry, elemVolVars, scvf);
95 return flux;
96 }
97};
98} // end namespace Dumux
99
100#endif
Element-wise calculation of the residual for the shallow water equations.
Definition: freeflow/shallowwater/localresidual.hh:28
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:80
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:58
Defines all properties used in Dumux.
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:34
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Definition: adapt.hh:17
A helper to deduce a vector with the same size as numbers of equations.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.