12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH 
   13#define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH 
   15#include <dune/common/fvector.hh> 
   29template< 
class IvIndexSet, 
class Scalar, 
int dim, 
int dimWorld>
 
   38    using ctype = 
typename GlobalCoordinate::value_type;
 
   56    template<
class MpfaHelper, 
class FVElementGeometry, 
class SubControlVolume>
 
   58                                     const FVElementGeometry& fvGeometry,
 
   59                                     const SubControlVolume& scv,
 
   61                                     const IvIndexSet& indexSet)
 
   62    : indexSet_(&indexSet)
 
   63    , globalScvIndex_(scv.dofIndex())
 
   64    , localDofIndex_(localIndex)
 
   67        const auto& 
center = scv.center();
 
   71        for (
unsigned int coordIdx = 0; coordIdx < 
myDimension; ++coordIdx)
 
   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;
 
   79        nus_ = helper.calculateInnerNormals(localBasis);
 
   80        detX_ = helper.calculateDetX(localBasis);
 
 
   89    { 
return globalScvIndex_; }
 
 
   93    { 
return localDofIndex_; }
 
 
   99        return indexSet_->localScvfIndex(localDofIndex_, coordDir);
 
 
  106        return nus_[coordDir];
 
 
  110    const IvIndexSet* indexSet_;
 
 
  124template< 
class IvIndexSet >
 
  145    template< 
class SubControlVolumeFace >
 
  151    , scvfIdxGlobal_(scvf.index())
 
  152    , localDofIndex_(localDofIdx)
 
  153    , neighborScvIndicesLocal_(&localScvIndices)
 
 
 
typename GlobalCoordinate::value_type ctype
Definition discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:38
typename IndexSet::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 IndexSet::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
const ScvfNeighborLocalIndexSet & neighboringLocalScvIndices() const
Returns the local indices of the scvs neighboring this scvf.
Definition discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:164
typename IndexSet::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
Definition discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:167
typename IndexSet::ScvfNeighborLocalIndexSet ScvfNeighborLocalIndexSet
Definition discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:127
CCMpfaOInteractionVolumeLocalScvf()=default
The default constructor.
typename IndexSet::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