3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 3 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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
26
27#include <dune/common/fvector.hh>
28
29namespace Dumux
30{
31
41template< class IvIndexSet, class Scalar, int dim, int dimWorld>
43{
44
45public:
46 // export some types
47 using GridIndexType = typename IvIndexSet::GridIndexType;
48 using LocalIndexType = typename IvIndexSet::LocalIndexType;
49 using GlobalCoordinate = Dune::FieldVector<Scalar, dimWorld>;
50 using ctype = typename GlobalCoordinate::value_type;
51 using LocalBasis = std::array< GlobalCoordinate, dim >;
52
53 static constexpr int myDimension = dim;
54 static constexpr int worldDimension = dimWorld;
55
58
68 template<class MpfaHelper, class FVElementGeometry, class SubControlVolume>
69 CCMpfaOInteractionVolumeLocalScv(const MpfaHelper& helper,
70 const FVElementGeometry& fvGeometry,
71 const SubControlVolume& scv,
72 const LocalIndexType localIndex,
73 const IvIndexSet& indexSet)
74 : indexSet_(&indexSet)
75 , globalScvIndex_(scv.dofIndex())
76 , localDofIndex_(localIndex)
77 {
78 // center of the global scv
79 const auto& center = scv.center();
80
81 // set up local basis
82 LocalBasis localBasis;
83 for (unsigned int coordIdx = 0; coordIdx < myDimension; ++coordIdx)
84 {
85 const auto scvfIdx = indexSet.nodalIndexSet().gridScvfIndex(localDofIndex_, coordIdx);
86 const auto& scvf = fvGeometry.scvf(scvfIdx);
87 localBasis[coordIdx] = scvf.ipGlobal();
88 localBasis[coordIdx] -= center;
89 }
90
91 nus_ = helper.calculateInnerNormals(localBasis);
92 detX_ = helper.calculateDetX(localBasis);
93 }
94
96 ctype detX() const
97 { return detX_; }
98
101 { return globalScvIndex_; }
102
105 { return localDofIndex_; }
106
108 LocalIndexType localScvfIndex(unsigned int coordDir) const
109 {
110 assert(coordDir < myDimension);
111 return indexSet_->localScvfIndex(localDofIndex_, coordDir);
112 }
113
115 const GlobalCoordinate& nu(unsigned int coordDir) const
116 {
117 assert(coordDir < myDimension);
118 return nus_[coordDir];
119 }
120
121private:
122 const IvIndexSet* indexSet_;
123 GridIndexType globalScvIndex_;
124 LocalIndexType localDofIndex_;
125 LocalBasis nus_;
126 ctype detX_;
127};
128
136template< class IvIndexSet >
138{
139 using ScvfNeighborLocalIndexSet = typename IvIndexSet::ScvfNeighborLocalIndexSet;
140
141public:
142 // export index types
143 using GridIndexType = typename IvIndexSet::GridIndexType;
144 using LocalIndexType = typename IvIndexSet::LocalIndexType;
145
148
157 template< class SubControlVolumeFace >
158 CCMpfaOInteractionVolumeLocalScvf(const SubControlVolumeFace& scvf,
159 const ScvfNeighborLocalIndexSet& localScvIndices,
160 const LocalIndexType localDofIdx,
161 const bool isDirichlet)
162 : isDirichlet_(isDirichlet)
163 , scvfIdxGlobal_(scvf.index())
164 , localDofIndex_(localDofIdx)
165 , neighborScvIndicesLocal_(&localScvIndices)
166 {}
167
170 LocalIndexType localDofIndex() const { return localDofIndex_; }
171
173 GridIndexType gridScvfIndex() const { return scvfIdxGlobal_; }
174
176 const ScvfNeighborLocalIndexSet& neighboringLocalScvIndices() const { return *neighborScvIndicesLocal_; }
177
179 bool isDirichlet() const { return isDirichlet_; }
180
181private:
182 bool isDirichlet_;
183 GridIndexType scvfIdxGlobal_;
184 LocalIndexType localDofIndex_;
185 const ScvfNeighborLocalIndexSet* neighborScvIndicesLocal_;
186};
187
188} // end namespace Dumux
189
190#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 GlobalCoordinate::value_type ctype
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:50
typename IvIndexSet::GridIndexType GridIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:47
static constexpr int worldDimension
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:54
static constexpr int myDimension
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:53
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:96
std::array< GlobalCoordinate, dim > LocalBasis
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:51
GridIndexType gridScvIndex() const
grid index related to this scv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:100
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:115
LocalIndexType localDofIndex() const
returns the index in the set of cell unknowns of the iv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:104
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:69
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:48
Dune::FieldVector< Scalar, dimWorld > GlobalCoordinate
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:49
LocalIndexType localScvfIndex(unsigned int coordDir) const
iv-local index of the coordir's scvf in this scv
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:108
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
GridIndexType gridScvfIndex() const
returns the grid view-global index of this scvf
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:173
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
CCMpfaOInteractionVolumeLocalScvf()=default
The default constructor.
typename IvIndexSet::LocalIndexType LocalIndexType
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:144
LocalIndexType localDofIndex() const
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:170
CCMpfaOInteractionVolumeLocalScvf(const SubControlVolumeFace &scvf, const ScvfNeighborLocalIndexSet &localScvIndices, const LocalIndexType localDofIdx, const bool isDirichlet)
The constructor.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:158