3.5-git
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
35
36namespace Dumux {
37
44template<class TypeTag>
45class BoxLocalResidual : public FVLocalResidual<TypeTag>
46{
51 using Element = typename GridView::template Codim<0>::Entity;
54 using FVElementGeometry = typename GridGeometry::LocalView;
55 using Extrusion = Extrusion_t<GridGeometry>;
56 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
57 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
58 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
60
61public:
63 using ParentType::ParentType;
64
67 const Problem& problem,
68 const Element& element,
69 const FVElementGeometry& fvGeometry,
70 const ElementVolumeVariables& elemVolVars,
71 const ElementBoundaryTypes& elemBcTypes,
72 const ElementFluxVariablesCache& elemFluxVarsCache,
73 const SubControlVolumeFace& scvf) const
74 {
75 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
76 if (!scvf.boundary())
77 {
78 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
79 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
80 residual[insideScv.localDofIndex()] += flux;
81 residual[outsideScv.localDofIndex()] -= flux;
82 }
83 else
84 {
85 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
86 residual[insideScv.localDofIndex()] += flux;
87 }
88 }
89
91 NumEqVector evalFlux(const Problem& problem,
92 const Element& element,
93 const FVElementGeometry& fvGeometry,
94 const ElementVolumeVariables& elemVolVars,
95 const ElementBoundaryTypes& elemBcTypes,
96 const ElementFluxVariablesCache& elemFluxVarsCache,
97 const SubControlVolumeFace& scvf) const
98 {
99 NumEqVector flux(0.0);
100
101 // inner faces
102 if (!scvf.boundary())
103 {
104 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
105 }
106
107 // boundary faces
108 else
109 {
110 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
111 const auto& bcTypes = elemBcTypes[scv.localDofIndex()];
112
113 // Treat Neumann and Robin ("solution dependent Neumann") boundary conditions.
114 // For Dirichlet there is no addition to the residual here but they
115 // are enforced strongly by replacing the residual entry afterwards.
116 if (bcTypes.hasNeumann())
117 {
118 auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
119
120 // multiply neumann fluxes with the area and the extrusion factor
121 neumannFluxes *= Extrusion::area(scvf)*elemVolVars[scv].extrusionFactor();
122
123 // only add fluxes to equations for which Neumann is set
124 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
125 if (bcTypes.isNeumann(eqIdx))
126 flux[eqIdx] += neumannFluxes[eqIdx];
127 }
128 }
129
130 return flux;
131 }
132};
133
134} // end namespace Dumux
135
136#endif
The element-wise residual for finite volume schemes.
A helper to deduce a vector with the same size as numbers of equations.
Helper classes to compute the integration elements.
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 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
Definition: propertysystem.hh:150
The element-wise residual for the box scheme.
Definition: boxlocalresidual.hh:46
typename ParentType::ElementResidualVector ElementResidualVector
Definition: boxlocalresidual.hh:62
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:91
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:66
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:47
Implementation & asImp()
Definition: fvlocalresidual.hh:501
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:486
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:69
Declares all properties used in Dumux.