24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH
27#include <dune/common/fvector.hh>
41template<
class IvIndexSet,
class Scalar,
int dim,
int dimWorld>
50 using ctype =
typename GlobalCoordinate::value_type;
68 template<
class MpfaHelper,
class FVElementGeometry,
class SubControlVolume>
70 const FVElementGeometry& fvGeometry,
71 const SubControlVolume& scv,
73 const IvIndexSet& indexSet)
74 : indexSet_(&indexSet)
75 , globalScvIndex_(scv.dofIndex())
76 , localDofIndex_(localIndex)
79 const auto& center = scv.center();
83 for (
unsigned int coordIdx = 0; coordIdx <
myDimension; ++coordIdx)
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;
91 nus_ = helper.calculateInnerNormals(localBasis);
92 detX_ = helper.calculateDetX(localBasis);
101 {
return globalScvIndex_; }
105 {
return localDofIndex_; }
111 return indexSet_->localScvfIndex(localDofIndex_, coordDir);
118 return nus_[coordDir];
122 const IvIndexSet* indexSet_;
136template<
class IvIndexSet >
157 template<
class SubControlVolumeFace >
163 , scvfIdxGlobal_(scvf.index())
164 , localDofIndex_(localDofIdx)
165 , neighborScvIndicesLocal_(&localScvIndices)
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