14#ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
15#define DUMUX_DISCRETIZATION_CC_MPFA_O_STATIC_INTERACTIONVOLUME_HH
17#include <dune/common/fmatrix.hh>
18#include <dune/common/fvector.hh>
19#include <dune/common/exceptions.hh>
37template<
class Traits >
class CCMpfaOStaticInteractionVolume;
51template<
class NI,
class S,
int C,
int F >
55 using GridIndexType =
typename NI::GridIndexType;
56 using LocalIndexType =
typename NI::LocalIndexType;
58 static constexpr int dim = NI::Traits::GridView::dimension;
59 static constexpr int dimWorld = NI::Traits::GridView::dimensionworld;
61 using DimVector = Dune::FieldVector<S, dim>;
62 using FaceOmegas = Dune::ReservedVector<DimVector, 2>;
67 using OmegaStorage = std::array< FaceOmegas, F >;
69 using AMatrix = Dune::FieldMatrix< S, F, F >;
70 using BMatrix = Dune::FieldMatrix< S, F, C >;
71 using CMatrix = Dune::FieldMatrix< S, F, F >;
72 using DMatrix = Dune::FieldMatrix< S, F, C >;
73 using TMatrix = Dune::FieldMatrix< S, F, C >;
74 using CellVector = Dune::FieldVector< S, C >;
75 using FaceVector = Dune::FieldVector< S, F >;
80 using GridView =
typename NI::Traits::GridView;
97 template<
class Problem,
class FVElementGeometry,
class ElemVolVars>
111template<
class Traits >
115 using GridView =
typename Traits::GridView;
116 using Element =
typename GridView::template Codim<0>::Entity;
118 using IndexSet =
typename Traits::IndexSet;
119 using GridIndexType =
typename IndexSet::GridIndexType;
120 using LocalIndexType =
typename IndexSet::LocalIndexType;
121 using Stencil =
typename IndexSet::NodalGridStencilType;
123 using LocalScvType =
typename Traits::LocalScvType;
124 using LocalScvfType =
typename Traits::LocalScvfType;
125 using LocalFaceData =
typename Traits::LocalFaceData;
127 static constexpr int numScvf = Traits::numScvfs;
128 static constexpr int numScv = Traits::numScvs;
132 static_assert(int(GridView::dimension)==int(GridView::dimensionworld),
"static iv does not work on surface grids");
139 { DUNE_THROW(Dune::InvalidStateException,
"Static interaction volume cannot be used on boundaries!"); }
146 template<
class Problem,
class FVElementGeometry >
147 void bind(
const IndexSet& indexSet,
148 const Problem& problem,
149 const FVElementGeometry& fvGeometry)
153 assert(indexSet.numScvs() == numScv);
154 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
157 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numScv; scvIdxLocal++)
159 elements_[scvIdxLocal] = fvGeometry.gridGeometry().element(
stencil()[scvIdxLocal] );
160 scvs_[scvIdxLocal] = LocalScvType(fvGeometry.gridGeometry().mpfaHelper(),
162 fvGeometry.scv(
stencil()[scvIdxLocal] ),
168 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < numScvf; ++faceIdxLocal)
170 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(faceIdxLocal));
171 assert(!scvf.boundary());
174 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
177 scvfs_[faceIdxLocal] = LocalScvfType(scvf, neighborScvIndicesLocal, faceIdxLocal,
false);
178 localFaceData_[faceIdxLocal*2] = LocalFaceData(faceIdxLocal, neighborScvIndicesLocal[0], scvf.index());
181 const auto outsideLocalScvIdx = neighborScvIndicesLocal[1];
182 for (
int coord = 0; coord < GridView::dimension; ++coord)
184 if (indexSet.localScvfIndex(outsideLocalScvIdx, coord) == faceIdxLocal)
186 const auto globalScvfIdx = indexSet.nodalIndexSet().gridScvfIndex(outsideLocalScvIdx, coord);
187 const auto& flipScvf = fvGeometry.scvf(globalScvfIdx);
188 localFaceData_[faceIdxLocal*2+1] = LocalFaceData(faceIdxLocal,
197 assert(localFaceData_[faceIdxLocal*2+1].ivLocalInsideScvIndex() == outsideLocalScvIdx);
219 {
return *stencil_; }
222 const Element&
element(LocalIndexType ivLocalScvIdx)
const
223 {
return elements_[ivLocalScvIdx]; }
226 const LocalScvfType&
localScvf(LocalIndexType ivLocalScvfIdx)
const
227 {
return scvfs_[ivLocalScvfIdx]; }
230 const LocalScvType&
localScv(LocalIndexType ivLocalScvIdx)
const
231 {
return scvs_[ivLocalScvIdx]; }
235 {
return localFaceData_; }
240 {
return dirichletData_; }
243 template<
class FVElementGeometry >
244 auto getScvGeometry(LocalIndexType ivLocalScvIdx,
const FVElementGeometry& fvGeometry)
const
254 template<
class IvIndexSetContainer,
257 class FlipScvfIndexSet >
259 ScvfIndexMap& scvfIndexMap,
260 const NodalIndexSet& nodalIndexSet,
261 const FlipScvfIndexSet& flipScvfIndexSet)
272 const Stencil * stencil_ =
nullptr;
275 std::array<Element, numScv> elements_;
276 std::array<LocalScvType, numScv> scvs_;
277 std::array<LocalScvfType, numScvf> scvfs_;
278 std::array<LocalFaceData, numScvf*2> localFaceData_;
281 std::array<DirichletData, 0> dirichletData_;
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:56
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:267
The interaction volume index set class for the mpfa-o scheme.
Definition: interactionvolumeindexset.hh:29
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:31
static ScvGeometry computeScvGeometry(LocalIndexType ivLocalScvIdx, const InteractionVolume &iv, const FVElementGeometry &fvGeometry)
returns the geometry of the i-th local scv
Definition: scvgeometryhelper.hh:61
Forward declaration of the o-method's static interaction volume.
Definition: staticinteractionvolume.hh:114
static constexpr std::size_t numFaces()
returns the number of primary scvfs of this interaction volume
Definition: staticinteractionvolume.hh:202
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:147
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: staticinteractionvolume.hh:143
static constexpr std::size_t numKnowns()
returns the number of (in this context) known solution values within this interaction volume
Definition: staticinteractionvolume.hh:210
static constexpr std::size_t numUnknowns()
returns the number of intermediate unknowns within this interaction volume
Definition: staticinteractionvolume.hh:206
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: staticinteractionvolume.hh:218
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:222
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: staticinteractionvolume.hh:249
const std::array< LocalFaceData, numScvf *2 > & localFaceData() const
returns a reference to the container with the local face data
Definition: staticinteractionvolume.hh:234
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: staticinteractionvolume.hh:258
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: staticinteractionvolume.hh:244
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: staticinteractionvolume.hh:230
const std::array< DirichletData, 0 > & dirichletData() const
Definition: staticinteractionvolume.hh:239
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: staticinteractionvolume.hh:226
static constexpr std::size_t numScvs()
returns the number of scvs embedded in this interaction volume
Definition: staticinteractionvolume.hh:214
General implementation of a data structure holding interaction volume-local information for a grid su...
Definition: localfacedata.hh:30
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: discretization/cellcentered/mpfa/omethod/localassembler.hh:36
Class for the index sets of the dual grid in mpfa schemes.
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:22
Base class for interaction volumes of mpfa methods.
Class for the index set within an interaction volume of the mpfa-o scheme.
Data structure holding interaction volume-local information for a grid subb-control volume face embed...
Define some often used mathematical functions.
The available mpfa schemes in Dumux.
Classes for sub control entities of the mpfa-o method in the context of facet coupling.
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
The default interaction volume traits class for the mpfa-o method with known size of the interaction ...
Definition: staticinteractionvolume.hh:53
static constexpr int numScvs
export the number of scvs in the interaction volumes
Definition: staticinteractionvolume.hh:92
typename NI::Traits::GridView GridView
export the type of grid view
Definition: staticinteractionvolume.hh:80
MVTraits MatVecTraits
export the matrix/vector traits to be used by the iv
Definition: staticinteractionvolume.hh:90
static constexpr int numScvfs
export the number of scvfs in the interaction volumes
Definition: staticinteractionvolume.hh:94
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme.
Definition: discretization/cellcentered/mpfa/omethod/localsubcontrolentities.hh:126
This does not work on surface grids.
Definition: staticinteractionvolume.hh:137
GridIndexType volVarIndex() const
Definition: staticinteractionvolume.hh:138