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