3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boxlocalresidual.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 *****************************************************************************/
25#ifndef DUMUX_BOX_LOCAL_RESIDUAL_HH
26#define DUMUX_BOX_LOCAL_RESIDUAL_HH
27
28#include <dune/geometry/type.hh>
29#include <dune/istl/matrix.hh>
30
34
35namespace Dumux {
36
43template<class TypeTag>
44class BoxLocalResidual : public FVLocalResidual<TypeTag>
45{
50 using Element = typename GridView::template Codim<0>::Entity;
53 using FVElementGeometry = typename GridGeometry::LocalView;
54 using Extrusion = Extrusion_t<GridGeometry>;
55 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
56 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
57 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
59
60public:
62 using ParentType::ParentType;
63
66 const Problem& problem,
67 const Element& element,
68 const FVElementGeometry& fvGeometry,
69 const ElementVolumeVariables& elemVolVars,
70 const ElementBoundaryTypes& elemBcTypes,
71 const ElementFluxVariablesCache& elemFluxVarsCache,
72 const SubControlVolumeFace& scvf) const
73 {
74 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
75 if (!scvf.boundary())
76 {
77 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
78 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
79 residual[insideScv.localDofIndex()] += flux;
80 residual[outsideScv.localDofIndex()] -= flux;
81 }
82 else
83 {
84 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
85 residual[insideScv.localDofIndex()] += flux;
86 }
87 }
88
90 NumEqVector evalFlux(const Problem& problem,
91 const Element& element,
92 const FVElementGeometry& fvGeometry,
93 const ElementVolumeVariables& elemVolVars,
94 const ElementBoundaryTypes& elemBcTypes,
95 const ElementFluxVariablesCache& elemFluxVarsCache,
96 const SubControlVolumeFace& scvf) const
97 {
98 NumEqVector flux(0.0);
99
100 // inner faces
101 if (!scvf.boundary())
102 {
103 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
104 }
105
106 // boundary faces
107 else
108 {
109 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
110 const auto& bcTypes = elemBcTypes[scv.localDofIndex()];
111
112 // Treat Neumann and Robin ("solution dependent Neumann") boundary conditions.
113 // For Dirichlet there is no addition to the residual here but they
114 // are enforced strongly by replacing the residual entry afterwards.
115 if (bcTypes.hasNeumann())
116 {
117 auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
118
119 // multiply neumann fluxes with the area and the extrusion factor
120 neumannFluxes *= Extrusion::area(scvf)*elemVolVars[scv].extrusionFactor();
121
122 // only add fluxes to equations for which Neumann is set
123 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
124 if (bcTypes.isNeumann(eqIdx))
125 flux[eqIdx] += neumannFluxes[eqIdx];
126 }
127 }
128
129 return flux;
130 }
131};
132
133} // end namespace Dumux
134
135#endif
The element-wise residual for finite volume schemes.
Helper classes to compute the integration elements.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
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 element-wise residual for the box scheme.
Definition: boxlocalresidual.hh:45
typename ParentType::ElementResidualVector ElementResidualVector
Definition: boxlocalresidual.hh:61
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate flux residuals for one sub control volume face
Definition: boxlocalresidual.hh:90
void evalFlux(ElementResidualVector &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate flux residuals for one sub control volume face and add to residual
Definition: boxlocalresidual.hh:65
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:46
Implementation & asImp()
Definition: fvlocalresidual.hh:503
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:488
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:68
Declares all properties used in Dumux.