3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
cclocalresidual.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_CC_LOCAL_RESIDUAL_HH
26#define DUMUX_CC_LOCAL_RESIDUAL_HH
27
32
33namespace Dumux {
34
40template<class TypeTag>
41class CCLocalResidual : public FVLocalResidual<TypeTag>
42{
45 using Element = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView::template Codim<0>::Entity;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
50 using FVElementGeometry = typename GetPropType<TypeTag, Properties::GridGeometry>::LocalView;
51 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
52
53public:
55 using ParentType::ParentType;
56
59 const Problem& problem,
60 const Element& element,
61 const FVElementGeometry& fvGeometry,
62 const ElementVolumeVariables& elemVolVars,
63 const ElementBoundaryTypes& elemBcTypes,
64 const ElementFluxVariablesCache& elemFluxVarsCache,
65 const SubControlVolumeFace& scvf) const
66 {
67 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
68 const auto localScvIdx = scv.localDofIndex();
69 residual[localScvIdx] += this->asImp().evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
70 }
71
73 NumEqVector evalFlux(const Problem& problem,
74 const Element& element,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementFluxVariablesCache& elemFluxVarsCache,
78 const SubControlVolumeFace& scvf) const
79 {
80 NumEqVector flux(0.0);
81
82 // inner faces
83 if (!scvf.boundary())
84 {
85 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
86 }
87
88 // boundary faces
89 else
90 {
91 const auto& bcTypes = problem.boundaryTypes(element, scvf);
92
93 // Dirichlet boundaries
94 if (bcTypes.hasDirichlet() && !bcTypes.hasNeumann())
95 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
96
97 // Neumann and Robin ("solution dependent Neumann") boundary conditions
98 else if (bcTypes.hasNeumann() && !bcTypes.hasDirichlet())
99 {
100 auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
101
102 // multiply neumann fluxes with the area and the extrusion factor
103 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
104 neumannFluxes *= scvf.area()*elemVolVars[scv].extrusionFactor();
105
106 flux += neumannFluxes;
107 }
108
109 else
110 DUNE_THROW(Dune::NotImplemented, "Mixed boundary conditions. Use pure boundary conditions by converting Dirichlet BCs to Robin BCs");
111 }
112
113 return flux;
114 }
115};
116
117} // end namespace Dumux
118
119#endif
The element-wise residual for finite volume schemes.
A arithmetic block vector type based on DUNE's reserved vector.
Helpers for deprecation.
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
Calculates the element-wise residual for the cell-centered discretization schemes.
Definition: cclocalresidual.hh:42
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
evaluate the flux residual for a sub control volume face
Definition: cclocalresidual.hh:73
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 the flux residual for a sub control volume face and add to residual
Definition: cclocalresidual.hh:58
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cclocalresidual.hh:54
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:45
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:66
Declares all properties used in Dumux.