24#ifndef DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
25#define DUMUX_DISCRETIZATION_CC_MPFA_INTERACTIONVOLUMEBASE_HH
29#include <dune/common/exceptions.hh>
30#include <dune/geometry/multilineargeometry.hh>
69 using GridView =
typename T::GridView;
70 using Element =
typename GridView::template Codim<0>::Entity;
72 using NodalStencilType =
typename T::IndexSet::NodalGridStencilType;
73 using LocalIndexType =
typename T::IndexSet::LocalIndexType;
74 using LocalScvType =
typename T::LocalScvType;
75 using LocalScvfType =
typename T::LocalScvfType;
77 using ScvGeometry = Dune::MultiLinearGeometry<
typename LocalScvType::ctype,
78 LocalScvType::myDimension,
79 LocalScvType::worldDimension>;
86 template<
class Problem,
class FVElementGeometry >
87 void bind(
const typename Traits::IndexSet& indexSet,
88 const Problem& problem,
89 const FVElementGeometry& fvGeometry)
90 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a bind() function"); }
94 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a numFaces() function"); }
98 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a numUnknowns() function"); }
102 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a numKnowns() function"); }
106 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a numScvs() function"); }
109 template<
class FVElementGeometry >
111 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a computeScvGeometry() function"); }
117 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a localFaceData() function"); }
121 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a stencil() function"); }
124 const LocalScvfType&
localScvf(LocalIndexType ivLocalScvfIdx)
const
125 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a localScvf() function"); }
128 const LocalScvType&
localScv(LocalIndexType ivLocalScvIdx)
const
129 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a localScv() function"); }
132 const Element&
element(LocalIndexType ivLocalScvIdx)
const
133 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide an element() function"); }
136 template<
class NodalIndexSet >
138 { DUNE_THROW(Dune::NotImplemented,
"Interaction volume does not provide a numIVAtVertex() function"); }
142 template<
class IvIndexSetContainer,
145 class FlipScvfIndexSet >
147 ScvfIndexMap& scvfIndexMap,
148 const NodalIndexSet& nodalIndexSet,
149 const FlipScvfIndexSet& flipScvfIndexSet)
150 { 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:68
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: interactionvolumebase.hh:128
std::size_t numFaces() const
returns the number of "primary" scvfs of this interaction volume
Definition: interactionvolumebase.hh:93
void bind(const typename Traits::IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Prepares everything for the assembly.
Definition: interactionvolumebase.hh:87
T Traits
state the traits type publicly
Definition: interactionvolumebase.hh:83
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the element in which the scv with the given local idx is embedded in
Definition: interactionvolumebase.hh:132
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: interactionvolumebase.hh:105
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: interactionvolumebase.hh:101
static std::size_t numIVAtVertex(const NodalIndexSet &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: interactionvolumebase.hh:137
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: interactionvolumebase.hh:97
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: interactionvolumebase.hh:146
const std::vector< typename Traits::LocalFaceData > & localFaceData() const
Definition: interactionvolumebase.hh:116
const NodalStencilType & stencil() const
returns the cell-stencil of this interaction volume
Definition: interactionvolumebase.hh:120
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: interactionvolumebase.hh:124
ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: interactionvolumebase.hh:110