3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 2 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_MULTDOMAIN_FACET_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
26#define DUMUX_MULTDOMAIN_FACET_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
27
28#include <array>
29
30#include <dune/common/fvector.hh>
32
33namespace Dumux {
34
44template< class IvIndexSet, class Scalar, int dim, int dimWorld>
46: public CCMpfaOInteractionVolumeLocalScv<IvIndexSet, Scalar, dim, dimWorld>
47{
49
50public:
51 // export some types
52 using typename ParentType::LocalIndexType;
53
54 // export dimension
55 static constexpr int myDimension = dim;
56
59
70 template<class MpfaHelper, class FVElementGeometry, class SubControlVolume, class IndexMap>
72 const FVElementGeometry& fvGeometry,
73 const SubControlVolume& scv,
74 const LocalIndexType localIndex,
75 const IvIndexSet& indexSet,
76 const IndexMap& scvfGridToLocalIndexMap)
77 : ParentType(helper, fvGeometry, scv, localIndex, indexSet)
78 {
79 // set up local scvf indices
80 const auto& nis = indexSet.nodalIndexSet();
81 for (unsigned int dir = 0; dir < myDimension; ++dir)
82 localScvfIndices_[dir] = scvfGridToLocalIndexMap.at(nis.gridScvfIndex(this->localDofIndex(), dir));
83 }
84
86 LocalIndexType localScvfIndex(unsigned int coordDir) const
87 {
88 assert(coordDir < myDimension);
89 return localScvfIndices_[coordDir];
90 }
91
92private:
93 std::array<LocalIndexType, dim> localScvfIndices_;
94};
95
103template< class IvIndexSet >
105: public CCMpfaOInteractionVolumeLocalScvf< IvIndexSet >
106{
108 using ScvfNeighborLocalIndexSet = typename IvIndexSet::ScvfNeighborLocalIndexSet;
109
110public:
111 // pull up index types
112 using typename ParentType::LocalIndexType;
113 using typename ParentType::GridIndexType;
114
116 using ParentType::ParentType;
117
128 template< class SubControlVolumeFace >
129 CCMpfaOFacetCouplingInteractionVolumeLocalScvf(const SubControlVolumeFace& scvf,
130 const ScvfNeighborLocalIndexSet& localScvIndices,
131 const LocalIndexType localDofIdx,
132 const bool isDirichlet,
133 const LocalIndexType coupledFacetLocalDof)
134 : ParentType(scvf, neighborLocalScvIndices_, localDofIdx, isDirichlet)
135 , isInteriorBoundary_(true)
136 , coupledFacetLocalDofIndex_(coupledFacetLocalDof)
137 , neighborLocalScvIndices_(localScvIndices)
138 {}
139
142 { assert(isInteriorBoundary_); return coupledFacetLocalDofIndex_; }
143
146 { return isOnInteriorBoundary() ? neighborLocalScvIndices_ : ParentType::neighboringLocalScvIndices(); }
147
149 bool isOnInteriorBoundary() const { return isInteriorBoundary_; }
150
151private:
152 bool isInteriorBoundary_{false};
153 LocalIndexType coupledFacetLocalDofIndex_{0};
154 ScvfNeighborLocalIndexSet neighborLocalScvIndices_;
155};
156
157} // end namespace Dumux
158
159#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:43
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:48
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:138
const ScvfNeighborLocalIndexSet & neighboringLocalScvIndices() const
Returns the local indices of the scvs neighboring this scvf.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:176
typename IvIndexSet::GridIndexType GridIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:143
bool isDirichlet() const
states if this is scvf is on a Dirichlet boundary
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:179
typename IvIndexSet::ScvfNeighborLocalIndexSet ScvfNeighborLocalIndexSet
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:139
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:144
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:47
static constexpr int myDimension
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:55
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:71
LocalIndexType localScvfIndex(unsigned int coordDir) const
iv-local index of the coordir's scvf in this scv
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:86
CCMpfaOFacetCouplingInteractionVolumeLocalScv()=default
The default constructor.
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:106
bool isOnInteriorBoundary() const
Returns true if this face is on an interior boundary.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:149
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:129
LocalIndexType coupledFacetLocalDofIndex() const
Returns the iv-local dof index of the coupled facet element.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:141
const ScvfNeighborLocalIndexSet & neighboringLocalScvIndices() const
Returns the local indices of the scvs neighboring this scvf.
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:145
typename IvIndexSet::ScvfNeighborLocalIndexSet ScvfNeighborLocalIndexSet
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:108
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:144
Classes for sub control entities of the mpfa-o method.