12#ifndef DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH 
   13#define DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH 
   17#include <dune/common/exceptions.hh> 
   18#include <dune/geometry/multilineargeometry.hh> 
   57    using GridView = 
typename T::GridView;
 
   58    using Element = 
typename GridView::template Codim<0>::Entity;
 
   60    using NodalStencilType = 
typename T::IndexSet::NodalGridStencilType;
 
   61    using LocalIndexType = 
typename T::IndexSet::LocalIndexType;
 
   62    using LocalScvType = 
typename T::LocalScvType;
 
   63    using LocalScvfType = 
typename T::LocalScvfType;
 
   65    using ScvGeometry = Dune::MultiLinearGeometry<
typename LocalScvType::ctype,
 
   66                                                  LocalScvType::myDimension,
 
   67                                                  LocalScvType::worldDimension>;
 
   74    template< 
class Problem, 
class FVElementGeometry >
 
   75    void bind(
const typename Traits::IndexSet& indexSet,
 
   76              const Problem& problem,
 
   77              const FVElementGeometry& fvGeometry)
 
   78    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a bind() function"); }
 
 
   82    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a numFaces() function"); }
 
 
   86    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a numUnknowns() function"); }
 
 
   90    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a numKnowns() function"); }
 
 
   94    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a numScvs() function"); }
 
 
   97    template< 
class FVElementGeometry >
 
   98    ScvGeometry 
computeScvGeometry(LocalIndexType ivLocalScvIdx, 
const FVElementGeometry& fvGeometry)
 
   99    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a computeScvGeometry() function"); }
 
 
  105    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a localFaceData() function"); }
 
 
  109    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a stencil() function"); }
 
 
  112    const LocalScvfType& 
localScvf(LocalIndexType ivLocalScvfIdx)
 const 
  113    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a localScvf() function"); }
 
 
  116    const LocalScvType& 
localScv(LocalIndexType ivLocalScvIdx)
 const 
  117    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a localScv() function"); }
 
 
  120    const Element& 
element(LocalIndexType ivLocalScvIdx)
 const 
  121    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide an element() function"); }
 
 
  124    template< 
class NodalIndexSet >
 
  126    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide a numIVAtVertex() function"); }
 
 
  130    template< 
class IvIndexSetContainer,
 
  133              class FlipScvfIndexSet >
 
  135                               ScvfIndexMap& scvfIndexMap,
 
  136                               const NodalIndexSet& nodalIndexSet,
 
  137                               const FlipScvfIndexSet& flipScvfIndexSet)
 
  138    { DUNE_THROW(Dune::NotImplemented, 
"Interaction volume does not provide an addIVIndexSets() function"); }
 
 
 
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition interactionvolumebase.hh:56
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition interactionvolumebase.hh:116
std::size_t numFaces() const
returns the number of "primary" scvfs of this interaction volume
Definition interactionvolumebase.hh:81
void bind(const typename Traits::IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Prepares everything for the assembly.
Definition interactionvolumebase.hh:75
T Traits
state the traits type publicly
Definition interactionvolumebase.hh:71
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the element in which the scv with the given local idx is embedded in
Definition interactionvolumebase.hh:120
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition interactionvolumebase.hh:93
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition interactionvolumebase.hh:89
static std::size_t numIVAtVertex(const NodalIndexSet &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition interactionvolumebase.hh:125
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition interactionvolumebase.hh:85
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition interactionvolumebase.hh:134
const std::vector< typename Traits::LocalFaceData > & localFaceData() const
Definition interactionvolumebase.hh:104
const NodalStencilType & stencil() const
returns the cell-stencil of this interaction volume
Definition interactionvolumebase.hh:108
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition interactionvolumebase.hh:112
ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition interactionvolumebase.hh:98