24#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
25#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
39template<
class Traits >
class CCMpfaOFacetCouplingInteractionVolume;
50template<
class NodalIndexSet,
class Scalar >
58 static constexpr int dim = NodalIndexSet::Traits::GridView::dimension;
59 static constexpr int dimWorld = NodalIndexSet::Traits::GridView::dimensionworld;
68 template<
class Problem,
class FVElementGeometry,
class ElemVolVars>
78template<
class Traits >
82 using GridView =
typename Traits::GridView;
83 using Element =
typename GridView::template Codim<0>::Entity;
85 using IndexSet =
typename Traits::IndexSet;
86 using GridIndexType =
typename IndexSet::GridIndexType;
87 using LocalIndexType =
typename IndexSet::LocalIndexType;
88 using Stencil =
typename IndexSet::NodalGridStencilType;
90 using LocalScvType =
typename Traits::LocalScvType;
91 using LocalScvfType =
typename Traits::LocalScvfType;
92 using LocalFaceData =
typename Traits::LocalFaceData;
101 GridIndexType scvfIdx_;
115 template<
class Problem,
class FVElementGeometry >
116 void bind(
const IndexSet& indexSet,
117 const Problem& problem,
118 const FVElementGeometry& fvGeometry)
122 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
125 std::size_t numFacetElems = 0;
126 std::size_t numOutsideFaces = 0;
127 std::vector<bool> isOnInteriorBoundary(indexSet.numFaces(),
false);
128 for (LocalIndexType fIdx = 0; fIdx < indexSet.numFaces(); ++fIdx)
130 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(fIdx));
131 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
132 if (problem.couplingManager().isOnInteriorBoundary(
element, scvf))
135 numOutsideFaces += scvf.numOutsideScvs();
136 isOnInteriorBoundary[fIdx] =
true;
137 interiorBoundaryData_.emplace_back( scvf.index() );
142 numFaces_ = indexSet.numFaces() + numOutsideFaces;
143 const auto numLocalScvs = indexSet.numScvs();
144 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
147 elements_.clear(); elements_.reserve(numLocalScvs);
148 scvs_.clear(); scvs_.reserve(numLocalScvs);
149 scvfs_.clear(); scvfs_.reserve(numFaces_);
150 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
151 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
155 numKnowns_ = numLocalScvs + numFacetElems;
158 std::unordered_map<GridIndexType, LocalIndexType> scvfIndexMap;
161 LocalIndexType facetElementCounter = 0;
162 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < indexSet.numFaces(); ++faceIdxLocal)
164 const auto gridScvfIdx = indexSet.gridScvfIndex(faceIdxLocal);
165 const auto& flipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[gridScvfIdx];
166 const auto& scvf = fvGeometry.scvf(gridScvfIdx);
167 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
170 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
171 const auto numNeighborScvs = neighborScvIndicesLocal.size();
174 const auto curLocalScvfIdx = scvfs_.size();
175 scvfIndexMap[gridScvfIdx] = curLocalScvfIdx;
176 localFaceData_.emplace_back(curLocalScvfIdx, neighborScvIndicesLocal[0], scvf.index());
179 if (isOnInteriorBoundary[faceIdxLocal])
181 const LocalIndexType facetLocalDofIdx = numLocalScvs + facetElementCounter++;
182 const bool isDirichlet = problem.interiorBoundaryTypes(
element, scvf).hasOnlyDirichlet();
186 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, facetLocalDofIdx,
true, facetLocalDofIdx);
188 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false, facetLocalDofIdx);
191 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
193 const auto outsideGridScvfIdx = flipScvfIdxSet[i-1];
194 const auto& flipScvf = fvGeometry.scvf(outsideGridScvfIdx);
195 const auto& outsideFlipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[outsideGridScvfIdx];
199 auto outsideNeighborScvIdxSet = neighborScvIndicesLocal;
200 outsideNeighborScvIdxSet[0] = outsideNeighborScvIdxSet[i];
201 for (LocalIndexType j = 0; j < outsideFlipScvfIdxSet.size(); ++j)
203 const auto flipScvfIdx = outsideFlipScvfIdxSet[j];
204 auto it = std::find(flipScvfIdxSet.begin(), flipScvfIdxSet.end(), flipScvfIdx);
207 if (it != flipScvfIdxSet.end())
208 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[
std::distance(flipScvfIdxSet.begin(), it)+1];
213 assert(flipScvfIdx == gridScvfIdx);
214 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[0];
218 scvfIndexMap[outsideGridScvfIdx] = curLocalScvfIdx+i;
219 localFaceData_.emplace_back(curLocalScvfIdx+i, outsideNeighborScvIdxSet[0], flipScvf.index());
221 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, facetLocalDofIdx,
true, facetLocalDofIdx);
223 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, numUnknowns_++,
false, facetLocalDofIdx);
228 else if (scvf.boundary())
230 if (problem.boundaryTypes(
element, scvf).hasOnlyDirichlet())
232 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++,
true);
233 dirichletData_.emplace_back(scvf.outsideScvIdx());
236 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false);
242 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false);
245 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
248 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
249 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
250 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
251 scvfIndexMap[flipScvfIndex] = curLocalScvfIdx;
252 localFaceData_.emplace_back(curLocalScvfIdx,
261 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
263 elements_.emplace_back(fvGeometry.gridGeometry().element(
stencil()[scvIdxLocal] ));
264 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
266 fvGeometry.scv(
stencil()[scvIdxLocal] ),
275 {
return numFaces_; }
279 {
return numUnknowns_; }
283 {
return numKnowns_; }
287 {
return scvs_.size(); }
291 {
return *stencil_; }
294 const Element&
element(LocalIndexType ivLocalScvIdx)
const
295 {
return elements_[ivLocalScvIdx]; }
298 const LocalScvfType&
localScvf(LocalIndexType ivLocalScvfIdx)
const
299 {
return scvfs_[ivLocalScvfIdx]; }
302 const LocalScvType&
localScv(LocalIndexType ivLocalScvIdx)
const
303 {
return scvs_[ivLocalScvIdx]; }
307 {
return localFaceData_; }
311 {
return dirichletData_; }
315 {
return interiorBoundaryData_; }
318 template<
class FVElementGeometry >
319 auto getScvGeometry(LocalIndexType ivLocalScvIdx,
const FVElementGeometry& fvGeometry)
const
329 template<
class IvIndexSetContainer,
332 class FlipScvfIndexSet >
334 ScvfIndexMap& scvfIndexMap,
335 const NodalIndexSet& nodalIndexSet,
336 const FlipScvfIndexSet& flipScvfIndexSet)
347 const Stencil* stencil_;
350 std::vector<Element> elements_;
351 std::vector<LocalScvType> scvs_;
352 std::vector<LocalScvfType> scvfs_;
353 std::vector<LocalFaceData> localFaceData_;
354 std::vector<DirichletData> dirichletData_;
355 std::vector<InteriorBoundaryData> interiorBoundaryData_;
358 std::size_t numFaces_;
359 std::size_t numUnknowns_;
360 std::size_t numKnowns_;
Helper class providing functionality to compute the geometry of the interaction-volume local sub-cont...
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.
ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: geometry/distance.hh:138
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
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:279
The default interaction volume traits class for the mpfa-o method. This uses dynamic types types for ...
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:60
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:132
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 facet coupling MPFA-O interaction volume.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:81
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:286
const LocalScvfType & localScvf(LocalIndexType ivLocalScvfIdx) const
returns the local scvf entity corresponding to a given iv-local scvf idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:298
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:306
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:319
const std::vector< DirichletData > & dirichletData() const
returns a reference to the information container on Dirichlet BCs within this iv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:310
typename CCMpfaOInteractionVolume< Traits >::DirichletData DirichletData
Reuse standard o-scheme's Dirichlet Data class.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:96
const std::vector< InteriorBoundaryData > & interiorBoundaryData() const
returns a reference to the data container on interior boundaries
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:314
std::size_t numKnowns() const
returns the number of (in this context) known solution values within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:282
void bind(const IndexSet &indexSet, const Problem &problem, const FVElementGeometry &fvGeometry)
Sets up the local scope for a given iv index set.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:116
const Element & element(LocalIndexType ivLocalScvIdx) const
returns the grid element corresponding to a given iv-local scv idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:294
static constexpr std::size_t numIVAtVertex(const NI &nodalIndexSet)
returns the number of interaction volumes living around a vertex
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:324
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:278
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:290
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:274
const LocalScvType & localScv(LocalIndexType ivLocalScvIdx) const
returns the local scv entity corresponding to a given iv-local scv idx
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:302
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:333
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:112
The default interaction volume traits class for the mpfa-o method in the context of facet coupling....
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:53
Define data structure to store which scvfs lie on interior boundaries.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:100
GridIndexType scvfIndex() const
Return corresponding scvf index.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:108
InteriorBoundaryData(GridIndexType scvfIdx)
Constructor.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:105
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:55
Class for the interaction volume-local sub-control volume used in the mpfa-o scheme in the context of...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:47
Class for the interaction volume-local sub-control volume face used in the mpfa-o scheme in the conte...
Definition: multidomain/facet/cellcentered/mpfa/localsubcontrolentities.hh:106
Class for the interaction volume of the mpfa-o scheme.
Classes for sub control entities of the mpfa-o method in the context of facet coupling.