version 3.8
multidomain/facet/cellcentered/mpfa/localsubcontrolentities.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_MULTDOMAIN_FACET_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
14#define DUMUX_MULTDOMAIN_FACET_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
15
16#include <array>
17
18#include <dune/common/fvector.hh>
20
21namespace Dumux {
22
32template< class IvIndexSet, class Scalar, int dim, int dimWorld>
34: public CCMpfaOInteractionVolumeLocalScv<IvIndexSet, Scalar, dim, dimWorld>
35{
37
38public:
39 // export some types
40 using typename ParentType::LocalIndexType;
41
42 // export dimension
43 static constexpr int myDimension = dim;
44
47
58 template<class MpfaHelper, class FVElementGeometry, class SubControlVolume, class IndexMap>
60 const FVElementGeometry& fvGeometry,
61 const SubControlVolume& scv,
62 const LocalIndexType localIndex,
63 const IvIndexSet& indexSet,
64 const IndexMap& scvfGridToLocalIndexMap)
65 : ParentType(helper, fvGeometry, scv, localIndex, indexSet)
66 {
67 // set up local scvf indices
68 const auto& nis = indexSet.nodalIndexSet();
69 for (unsigned int dir = 0; dir < myDimension; ++dir)
70 localScvfIndices_[dir] = scvfGridToLocalIndexMap.at(nis.gridScvfIndex(this->localDofIndex(), dir));
71 }
72
74 LocalIndexType localScvfIndex(unsigned int coordDir) const
75 {
76 assert(coordDir < myDimension);
77 return localScvfIndices_[coordDir];
78 }
79
80private:
81 std::array<LocalIndexType, dim> localScvfIndices_;
82};
83
91template< class IvIndexSet >
93: public CCMpfaOInteractionVolumeLocalScvf< IvIndexSet >
94{
96 using ScvfNeighborLocalIndexSet = typename IvIndexSet::ScvfNeighborLocalIndexSet;
97
98public:
99 // pull up index types
100 using typename ParentType::LocalIndexType;
101 using typename ParentType::GridIndexType;
102
104 using ParentType::ParentType;
105
116 template< class SubControlVolumeFace >
117 CCMpfaOFacetCouplingInteractionVolumeLocalScvf(const SubControlVolumeFace& scvf,
118 const ScvfNeighborLocalIndexSet& localScvIndices,
119 const LocalIndexType localDofIdx,
120 const bool isDirichlet,
121 const LocalIndexType coupledFacetLocalDof)
122 : ParentType(scvf, neighborLocalScvIndices_, localDofIdx, isDirichlet)
123 , isInteriorBoundary_(true)
124 , coupledFacetLocalDofIndex_(coupledFacetLocalDof)
125 , neighborLocalScvIndices_(localScvIndices)
126 {}
127
130 { assert(isInteriorBoundary_); return coupledFacetLocalDofIndex_; }
131
134 { return isOnInteriorBoundary() ? neighborLocalScvIndices_ : ParentType::neighboringLocalScvIndices(); }
135
137 bool isOnInteriorBoundary() const { return isInteriorBoundary_; }
138
139private:
140 bool isInteriorBoundary_{false};
141 LocalIndexType coupledFacetLocalDofIndex_{0};
142 ScvfNeighborLocalIndexSet neighborLocalScvIndices_;
143};
144
145} // end namespace Dumux
146
147#endif
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme in the context of...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:35
static constexpr int myDimension
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:43
CCMpfaOFacetCouplingInteractionVolumeLocalScv(const MpfaHelper &helper, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const LocalIndexType localIndex, const IvIndexSet &indexSet, const IndexMap &scvfGridToLocalIndexMap)
The constructor.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:59
LocalIndexType localScvfIndex(unsigned int coordDir) const
iv-local index of the coordir's scvf in this scv
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:74
CCMpfaOFacetCouplingInteractionVolumeLocalScv()=default
The default constructor.
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:31
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:36
Classes for sub control entities of the mpfa-o method.
Definition: adapt.hh:17
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme in the conte...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:94
bool isOnInteriorBoundary() const
Returns true if this face is on an interior boundary.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:137
CCMpfaOFacetCouplingInteractionVolumeLocalScvf(const SubControlVolumeFace &scvf, const ScvfNeighborLocalIndexSet &localScvIndices, const LocalIndexType localDofIdx, const bool isDirichlet, const LocalIndexType coupledFacetLocalDof)
The constructor for interior boundary faces.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:117
LocalIndexType coupledFacetLocalDofIndex() const
Returns the iv-local dof index of the coupled facet element.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:129
const ScvfNeighborLocalIndexSet & neighboringLocalScvIndices() const
Returns the local indices of the scvs neighboring this scvf.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:133
typename IvIndexSet::ScvfNeighborLocalIndexSet ScvfNeighborLocalIndexSet
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:96
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:132
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:126
const ScvfNeighborLocalIndexSet & neighboringLocalScvIndices() const
Returns the local indices of the scvs neighboring this scvf.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:164
typename IvIndexSet::GridIndexType GridIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:131
bool isDirichlet() const
states if this is scvf is on a Dirichlet boundary
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:167
typename IvIndexSet::ScvfNeighborLocalIndexSet ScvfNeighborLocalIndexSet
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:127
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:132