version 3.8
multidomain/facet/cellcentered/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//
15#ifndef DUMUX_FACETCOUPLING_CC_LOCAL_RESIDUAL_HH
16#define DUMUX_FACETCOUPLING_CC_LOCAL_RESIDUAL_HH
17
21
22namespace Dumux {
23
31template<class TypeTag>
33{
35
37 using ElementVolumeVariables = typename GridVariables::GridVolumeVariables::LocalView;
38 using ElementFluxVariablesCache = typename GridVariables::GridFluxVariablesCache::LocalView;
39
40 using GridGeometry = typename GridVariables::GridGeometry;
41 using FVElementGeometry = typename GridGeometry::LocalView;
42 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
43 using Element = typename GridGeometry::GridView::template Codim<0>::Entity;
44
46
47public:
49 using ParentType::ParentType;
52
55 template< class Problem >
56 NumEqVector evalFlux(const Problem& problem,
57 const Element& element,
58 const FVElementGeometry& fvGeometry,
59 const ElementVolumeVariables& elemVolVars,
60 const ElementFluxVariablesCache& elemFluxVarsCache,
61 const SubControlVolumeFace& scvf) const
62 {
63 // Even if scvf.boundary=true, compute flux on interior boundaries
64 if (problem.couplingManager().isOnInteriorBoundary(element, scvf))
65 return this->asImp().computeFlux(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache);
66 else
67 return ParentType::evalFlux(problem, element, fvGeometry, elemVolVars, elemFluxVarsCache, scvf);
68 }
69};
70
71} // end namespace Dumux
72
73#endif
Calculates the element-wise residual for cell-centered discretization schemes.
Calculates the element-wise residual for cell-centered discretization schemes in models where couplin...
Definition: multidomain/facet/cellcentered/localresidual.hh:33
NumEqVector evalFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Definition: multidomain/facet/cellcentered/localresidual.hh:56
Calculates the element-wise residual for the cell-centered discretization schemes.
Definition: cclocalresidual.hh:31
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:49
typename ParentType::ElementResidualVector ElementResidualVector
Definition: cclocalresidual.hh:45
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
Defines all properties used in Dumux.
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
A helper to deduce a vector with the same size as numbers of equations.