version 3.8
multidomain/facet/box/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//
13#ifndef DUMUX_FACETCOUPLING_BOX_LOCAL_RESIDUAL_HH
14#define DUMUX_FACETCOUPLING_BOX_LOCAL_RESIDUAL_HH
15
16#include <dune/geometry/type.hh>
17
22
23namespace Dumux {
24
30template<class TypeTag>
32{
37 using FVElementGeometry = typename GridGeometry::LocalView;
38 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
39 using Extrusion = Extrusion_t<GridGeometry>;
40 using GridView = typename GridGeometry::GridView;
41 using Element = typename GridView::template Codim<0>::Entity;
43 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
44 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
46
47public:
49 using ParentType::ParentType;
50
53 const Problem& problem,
54 const Element& element,
55 const FVElementGeometry& fvGeometry,
56 const ElementVolumeVariables& elemVolVars,
57 const ElementBoundaryTypes& elemBcTypes,
58 const ElementFluxVariablesCache& elemFluxVarsCache,
59 const SubControlVolumeFace& scvf) const
60 {
61 const auto flux = evalFlux(problem, element, fvGeometry, elemVolVars, elemBcTypes, elemFluxVarsCache, scvf);
62
63 // inner faces
64 if (!scvf.boundary() && !scvf.interiorBoundary())
65 {
66 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
67 const auto& outsideScv = fvGeometry.scv(scvf.outsideScvIdx());
68 residual[insideScv.localDofIndex()] += flux;
69 residual[outsideScv.localDofIndex()] -= flux;
70 }
71
72 // interior/domain boundary faces
73 else
74 {
75 const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx());
76 residual[insideScv.localDofIndex()] += flux;
77 }
78 }
79
81 NumEqVector evalFlux(const Problem& problem,
82 const Element& element,
83 const FVElementGeometry& fvGeometry,
84 const ElementVolumeVariables& elemVolVars,
85 const ElementBoundaryTypes& elemBcTypes,
86 const ElementFluxVariablesCache& elemFluxVarsCache,
87 const SubControlVolumeFace& scvf) const
88 {
89 NumEqVector flux(0.0);
90
91 // inner or interior boundary faces
92 if (!scvf.boundary() || scvf.interiorBoundary())
93 flux += this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
94
95 // boundary faces
96 else
97 {
98 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
99 const auto& bcTypes = elemBcTypes.get(fvGeometry, scv);
100
101 // Neumann and Robin ("solution dependent Neumann") boundary conditions
102 if (bcTypes.hasNeumann())
103 {
104 auto neumannFluxes = problem.neumann(element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
105
106 // multiply neumann fluxes with the area and the extrusion factor
107 neumannFluxes *= Extrusion::area(fvGeometry, scvf)*elemVolVars[scv].extrusionFactor();
108
109 // only add fluxes to equations for which Neumann is set
110 for (int eqIdx = 0; eqIdx < NumEqVector::dimension; ++eqIdx)
111 if (bcTypes.isNeumann(eqIdx))
112 flux[eqIdx] += neumannFluxes[eqIdx];
113 }
114 }
115
116 return flux;
117 }
118};
119
120} // end namespace Dumux
121
122#endif
The element-wise residual for the box scheme.
Definition: multidomain/facet/box/localresidual.hh:32
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: multidomain/facet/box/localresidual.hh:81
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: multidomain/facet/box/localresidual.hh:52
The element-wise residual for finite volume schemes.
Definition: fvlocalresidual.hh:35
Implementation & asImp()
Definition: fvlocalresidual.hh:487
const Problem & problem() const
the problem
Definition: fvlocalresidual.hh:472
ReservedBlockVector< NumEqVector, FVElementGeometry::maxNumElementScvs > ElementResidualVector
the container storing all element residuals
Definition: fvlocalresidual.hh:57
A arithmetic block vector type based on DUNE's reserved vector.
Definition: reservedblockvector.hh:26
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
The element-wise residual for finite volume schemes.
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
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
A helper to deduce a vector with the same size as numbers of equations.