26#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
27#define DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
29#include <dune/common/fmatrix.hh>
30#include <dune/common/fvector.hh>
31#include <dune/common/exceptions.hh>
49template<
class Traits >
class CCMpfaOStaticInteractionVolume;
63template<
class NI,
class S,
int C,
int F >
67 using GridIndexType =
typename NI::GridIndexType;
68 using LocalIndexType =
typename NI::LocalIndexType;
70 static constexpr int dim = NI::Traits::GridView::dimension;
71 static constexpr int dimWorld = NI::Traits::GridView::dimensionworld;
73 using DimVector = Dune::FieldVector<S, dim>;
74 using FaceOmegas = Dune::ReservedVector<DimVector, 2>;
79 using OmegaStorage = std::array< FaceOmegas, F >;
81 using AMatrix = Dune::FieldMatrix< S, F, F >;
82 using BMatrix = Dune::FieldMatrix< S, F, C >;
83 using CMatrix = Dune::FieldMatrix< S, F, F >;
84 using DMatrix = Dune::FieldMatrix< S, F, C >;
85 using TMatrix = Dune::FieldMatrix< S, F, C >;
86 using CellVector = Dune::FieldVector< S, C >;
87 using FaceVector = Dune::FieldVector< S, F >;
92 using GridView =
typename NI::Traits::GridView;
109 template<
class Problem,
class FVElementGeometry,
class ElemVolVars>
123template<
class Traits >
127 using GridView =
typename Traits::GridView;
128 using Element =
typename GridView::template Codim<0>::Entity;
130 using IndexSet =
typename Traits::IndexSet;
131 using GridIndexType =
typename IndexSet::GridIndexType;
132 using LocalIndexType =
typename IndexSet::LocalIndexType;
133 using Stencil =
typename IndexSet::NodalGridStencilType;
135 using LocalScvType =
typename Traits::LocalScvType;
136 using LocalScvfType =
typename Traits::LocalScvfType;
137 using LocalFaceData =
typename Traits::LocalFaceData;
139 static constexpr int numScvf = Traits::numScvfs;
140 static constexpr int numScv = Traits::numScvs;
144 static_assert(int(GridView::dimension)==int(GridView::dimensionworld),
"static iv does not work on surface grids");
151 { DUNE_THROW(Dune::InvalidStateException,
"Static interaction volume cannot be used on bounaries!"); }
158 template<
class Problem,
class FVElementGeometry >
159 void bind(
const IndexSet& indexSet,
160 const Problem& problem,
161 const FVElementGeometry& fvGeometry)
165 assert(indexSet.numScvs() == numScv);
166 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
169 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numScv; scvIdxLocal++)
171 elements_[scvIdxLocal] = fvGeometry.gridGeometry().element(
stencil()[scvIdxLocal] );
172 scvs_[scvIdxLocal] = LocalScvType(fvGeometry.gridGeometry().mpfaHelper(),
174 fvGeometry.scv(
stencil()[scvIdxLocal] ),
180 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < numScvf; ++faceIdxLocal)
182 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(faceIdxLocal));
183 assert(!scvf.boundary());
186 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
189 scvfs_[faceIdxLocal] = LocalScvfType(scvf, neighborScvIndicesLocal, faceIdxLocal,
false);
190 localFaceData_[faceIdxLocal*2] = LocalFaceData(faceIdxLocal, neighborScvIndicesLocal[0], scvf.index());
193 const auto outsideLocalScvIdx = neighborScvIndicesLocal[1];
194 for (
int coord = 0; coord < GridView::dimension; ++coord)
196 if (indexSet.localScvfIndex(outsideLocalScvIdx, coord) == faceIdxLocal)
198 const auto globalScvfIdx = indexSet.nodalIndexSet().gridScvfIndex(outsideLocalScvIdx, coord);
199 const auto& flipScvf = fvGeometry.scvf(globalScvfIdx);
200 localFaceData_[faceIdxLocal*2+1] = LocalFaceData(faceIdxLocal,
209 assert(localFaceData_[faceIdxLocal*2+1].ivLocalInsideScvIndex() == outsideLocalScvIdx);
231 {
return *stencil_; }
234 const Element&
element(LocalIndexType ivLocalScvIdx)
const
235 {
return elements_[ivLocalScvIdx]; }
238 const LocalScvfType&
localScvf(LocalIndexType ivLocalScvfIdx)
const
239 {
return scvfs_[ivLocalScvfIdx]; }
242 const LocalScvType&
localScv(LocalIndexType ivLocalScvIdx)
const
243 {
return scvs_[ivLocalScvIdx]; }
247 {
return localFaceData_; }
252 {
return dirichletData_; }
255 template<
class FVElementGeometry >
256 auto getScvGeometry(LocalIndexType ivLocalScvIdx,
const FVElementGeometry& fvGeometry)
const
266 template<
class IvIndexSetContainer,
269 class FlipScvfIndexSet >
271 ScvfIndexMap& scvfIndexMap,
272 const NodalIndexSet& nodalIndexSet,
273 const FlipScvfIndexSet& flipScvfIndexSet)
284 const Stencil * stencil_ =
nullptr;
287 std::array<Element, numScv> elements_;
288 std::array<LocalScvType, numScv> scvs_;
289 std::array<LocalScvfType, numScvf> scvfs_;
290 std::array<LocalFaceData, numScvf*2> localFaceData_;
293 std::array<DirichletData, 0> dirichletData_;
Define some often used mathematical functions.
Class for the index set within an interaction volume of the mpfa-o scheme.
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
Class for the index sets of the dual grid in mpfa schemes.
Base class for interaction volumes of mpfa methods.
Data structure holding interaction volume-local information for a grid subb-control volume face embed...
The available mpfa schemes in Dumux.
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:34
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:68
General implementation of a data structure holding interaction volume-local information for a grid su...
Definition: localfacedata.hh:42
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:279
The interaction volume index set class for the mpfa-o scheme.
Definition: interactionvolumeindexset.hh:41
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: discretization/cellcentered/mpfa/omethod/localassembler.hh:48
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:43
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:138
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:73
Forward declaration of the o-method's static interaction volume.
Definition: staticinteractionvolume.hh:126
static constexpr std::size_t numFaces()
returns the number of primary scvfs of this interaction volume
Definition: staticinteractionvolume.hh:214
void bind(const IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Sets up the local scope for a given iv index set.
Definition: staticinteractionvolume.hh:159
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: staticinteractionvolume.hh:155
static constexpr std::size_t numKnowns()
returns the number of (in this context) known solution values within this interaction volume
Definition: staticinteractionvolume.hh:222
static constexpr std::size_t numUnknowns()
returns the number of intermediate unknowns within this interaction volume
Definition: staticinteractionvolume.hh:218
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: staticinteractionvolume.hh:230
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:234
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: staticinteractionvolume.hh:261
const std::array< LocalFaceData, numScvf *2 > & localFaceData() const
returns a reference to the container with the local face data
Definition: staticinteractionvolume.hh:246
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: staticinteractionvolume.hh:270
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: staticinteractionvolume.hh:256
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:242
const std::array< DirichletData, 0 > & dirichletData() const
Definition: staticinteractionvolume.hh:251
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: staticinteractionvolume.hh:238
static constexpr std::size_t numScvs()
returns the number of scvs embedded in this interaction volume
Definition: staticinteractionvolume.hh:226
The default interaction volume traits class for the mpfa-o method with known size of the interaction ...
Definition: staticinteractionvolume.hh:65
static constexpr int numScvs
export the number of scvs in the interaction volumes
Definition: staticinteractionvolume.hh:104
typename NI::Traits::GridView GridView
export the type of grid view
Definition: staticinteractionvolume.hh:92
MVTraits MatVecTraits
export the matrix/vector traits to be used by the iv
Definition: staticinteractionvolume.hh:102
static constexpr int numScvfs
export the number of scvfs in the interaction volumes
Definition: staticinteractionvolume.hh:106
This does not work on surface grids.
Definition: staticinteractionvolume.hh:149
GridIndexType volVarIndex() const
Definition: staticinteractionvolume.hh:150
Classes for sub control entities of the mpfa-o method in the context of facet coupling.