version 3.10-dev
discretization/cellcentered/mpfa/omethod/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//
12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
13#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
14
15#include <dune/common/fvector.hh>
16
17namespace Dumux
18{
19
29template< class IvIndexSet, class Scalar, int dim, int dimWorld>
31{
32
33public:
34 // export some types
35 using GridIndexType = typename IvIndexSet::GridIndexType;
36 using LocalIndexType = typename IvIndexSet::LocalIndexType;
37 using GlobalCoordinate = Dune::FieldVector<Scalar, dimWorld>;
38 using ctype = typename GlobalCoordinate::value_type;
39 using LocalBasis = std::array< GlobalCoordinate, dim >;
40
41 static constexpr int myDimension = dim;
42 static constexpr int worldDimension = dimWorld;
43
46
56 template<class MpfaHelper, class FVElementGeometry, class SubControlVolume>
57 CCMpfaOInteractionVolumeLocalScv(const MpfaHelper& helper,
58 const FVElementGeometry& fvGeometry,
59 const SubControlVolume& scv,
60 const LocalIndexType localIndex,
61 const IvIndexSet& indexSet)
62 : indexSet_(&indexSet)
63 , globalScvIndex_(scv.dofIndex())
64 , localDofIndex_(localIndex)
65 {
66 // center of the global scv
67 const auto& center = scv.center();
68
69 // set up local basis
70 LocalBasis localBasis;
71 for (unsigned int coordIdx = 0; coordIdx < myDimension; ++coordIdx)
72 {
73 const auto scvfIdx = indexSet.nodalIndexSet().gridScvfIndex(localDofIndex_, coordIdx);
74 const auto& scvf = fvGeometry.scvf(scvfIdx);
75 localBasis[coordIdx] = scvf.ipGlobal();
76 localBasis[coordIdx] -= center;
77 }
78
79 nus_ = helper.calculateInnerNormals(localBasis);
80 detX_ = helper.calculateDetX(localBasis);
81 }
82
84 ctype detX() const
85 { return detX_; }
86
89 { return globalScvIndex_; }
90
93 { return localDofIndex_; }
94
96 LocalIndexType localScvfIndex(unsigned int coordDir) const
97 {
98 assert(coordDir < myDimension);
99 return indexSet_->localScvfIndex(localDofIndex_, coordDir);
100 }
101
103 const GlobalCoordinate& nu(unsigned int coordDir) const
104 {
105 assert(coordDir < myDimension);
106 return nus_[coordDir];
107 }
108
109private:
110 const IvIndexSet* indexSet_;
111 GridIndexType globalScvIndex_;
112 LocalIndexType localDofIndex_;
113 LocalBasis nus_;
114 ctype detX_;
115};
116
124template< class IvIndexSet >
126{
127 using ScvfNeighborLocalIndexSet = typename IvIndexSet::ScvfNeighborLocalIndexSet;
128
129public:
130 // export index types
131 using GridIndexType = typename IvIndexSet::GridIndexType;
132 using LocalIndexType = typename IvIndexSet::LocalIndexType;
133
136
145 template< class SubControlVolumeFace >
146 CCMpfaOInteractionVolumeLocalScvf(const SubControlVolumeFace& scvf,
147 const ScvfNeighborLocalIndexSet& localScvIndices,
148 const LocalIndexType localDofIdx,
149 const bool isDirichlet)
150 : isDirichlet_(isDirichlet)
151 , scvfIdxGlobal_(scvf.index())
152 , localDofIndex_(localDofIdx)
153 , neighborScvIndicesLocal_(&localScvIndices)
154 {}
155
158 LocalIndexType localDofIndex() const { return localDofIndex_; }
159
161 GridIndexType gridScvfIndex() const { return scvfIdxGlobal_; }
162
164 const ScvfNeighborLocalIndexSet& neighboringLocalScvIndices() const { return *neighborScvIndicesLocal_; }
165
167 bool isDirichlet() const { return isDirichlet_; }
168
169private:
170 bool isDirichlet_;
171 GridIndexType scvfIdxGlobal_;
172 LocalIndexType localDofIndex_;
173 const ScvfNeighborLocalIndexSet* neighborScvIndicesLocal_;
174};
175
176} // end namespace Dumux
177
178#endif
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:31
typename GlobalCoordinate::value_type ctype
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:38
typename IvIndexSet::GridIndexType GridIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:35
static constexpr int worldDimension
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:42
static constexpr int myDimension
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:41
CCMpfaOInteractionVolumeLocalScv()=default
The default constructor.
ctype detX() const
detX is needed for setting up the omegas in the interaction volumes
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:84
std::array< GlobalCoordinate, dim > LocalBasis
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:39
GridIndexType gridScvIndex() const
grid index related to this scv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:88
const GlobalCoordinate & nu(unsigned int coordDir) const
the nu vectors are needed for setting up the omegas of the iv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:103
LocalIndexType localDofIndex() const
returns the index in the set of cell unknowns of the iv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:92
CCMpfaOInteractionVolumeLocalScv(const MpfaHelper &helper, const FVElementGeometry &fvGeometry, const SubControlVolume &scv, const LocalIndexType localIndex, const IvIndexSet &indexSet)
The constructor.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:57
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:36
Dune::FieldVector< Scalar, dimWorld > GlobalCoordinate
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:37
LocalIndexType localScvfIndex(unsigned int coordDir) const
iv-local index of the coordir's scvf in this scv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:96
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Definition: adapt.hh:17
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
GridIndexType gridScvfIndex() const
returns the grid view-global index of this scvf
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:161
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
CCMpfaOInteractionVolumeLocalScvf()=default
The default constructor.
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:132
LocalIndexType localDofIndex() const
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:158
CCMpfaOInteractionVolumeLocalScvf(const SubControlVolumeFace &scvf, const ScvfNeighborLocalIndexSet &localScvIndices, const LocalIndexType localDofIdx, const bool isDirichlet)
The constructor.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:146