12#ifndef DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
13#define DUMUX_MULTIDOMAIN_FACET_CC_MPFA_O_INTERACTIONVOLUME_HH
27template<
class Traits >
class CCMpfaOFacetCouplingInteractionVolume;
38template<
class NodalIndexSet,
class Scalar >
46 static constexpr int dim = NodalIndexSet::Traits::GridView::dimension;
47 static constexpr int dimWorld = NodalIndexSet::Traits::GridView::dimensionworld;
56 template<
class Problem,
class FVElementGeometry,
class ElemVolVars>
66template<
class Traits >
70 using GridView =
typename Traits::GridView;
71 using Element =
typename GridView::template Codim<0>::Entity;
73 using IndexSet =
typename Traits::IndexSet;
74 using GridIndexType =
typename IndexSet::GridIndexType;
75 using LocalIndexType =
typename IndexSet::LocalIndexType;
76 using Stencil =
typename IndexSet::NodalGridStencilType;
78 using LocalScvType =
typename Traits::LocalScvType;
79 using LocalScvfType =
typename Traits::LocalScvfType;
80 using LocalFaceData =
typename Traits::LocalFaceData;
89 GridIndexType scvfIdx_;
96 GridIndexType
scvfIndex()
const {
return scvfIdx_; }
103 template<
class Problem,
class FVElementGeometry >
104 void bind(
const IndexSet& indexSet,
105 const Problem& problem,
106 const FVElementGeometry& fvGeometry)
110 stencil_ = &indexSet.nodalIndexSet().gridScvIndices();
113 std::size_t numFacetElems = 0;
114 std::size_t numOutsideFaces = 0;
115 std::vector<bool> isOnInteriorBoundary(indexSet.numFaces(),
false);
116 for (LocalIndexType fIdx = 0; fIdx < indexSet.numFaces(); ++fIdx)
118 const auto& scvf = fvGeometry.scvf(indexSet.gridScvfIndex(fIdx));
119 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
120 if (problem.couplingManager().isOnInteriorBoundary(
element, scvf))
123 if (!scvf.boundary())
124 numOutsideFaces += scvf.numOutsideScvs();
125 isOnInteriorBoundary[fIdx] =
true;
126 interiorBoundaryData_.emplace_back( scvf.index() );
131 numFaces_ = indexSet.numFaces() + numOutsideFaces;
132 const auto numLocalScvs = indexSet.numScvs();
133 const auto numGlobalScvfs = indexSet.nodalIndexSet().numScvfs();
136 elements_.clear(); elements_.reserve(numLocalScvs);
137 scvs_.clear(); scvs_.reserve(numLocalScvs);
138 scvfs_.clear(); scvfs_.reserve(numFaces_);
139 localFaceData_.clear(); localFaceData_.reserve(numGlobalScvfs);
140 dirichletData_.clear(); dirichletData_.reserve(numFaces_);
144 numKnowns_ = numLocalScvs + numFacetElems;
147 std::unordered_map<GridIndexType, LocalIndexType> scvfIndexMap;
150 LocalIndexType facetElementCounter = 0;
151 for (LocalIndexType faceIdxLocal = 0; faceIdxLocal < indexSet.numFaces(); ++faceIdxLocal)
153 const auto gridScvfIdx = indexSet.gridScvfIndex(faceIdxLocal);
154 const auto& flipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[gridScvfIdx];
155 const auto& scvf = fvGeometry.scvf(gridScvfIdx);
156 const auto element = fvGeometry.gridGeometry().element(scvf.insideScvIdx());
159 const auto& neighborScvIndicesLocal = indexSet.neighboringLocalScvIndices(faceIdxLocal);
160 const auto numNeighborScvs = neighborScvIndicesLocal.size();
163 const auto curLocalScvfIdx = scvfs_.size();
164 scvfIndexMap[gridScvfIdx] = curLocalScvfIdx;
165 localFaceData_.emplace_back(curLocalScvfIdx, neighborScvIndicesLocal[0], scvf.index());
168 if (isOnInteriorBoundary[faceIdxLocal])
170 const LocalIndexType facetLocalDofIdx = numLocalScvs + facetElementCounter++;
171 const bool isDirichlet = problem.interiorBoundaryTypes(
element, scvf).hasOnlyDirichlet();
175 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, facetLocalDofIdx,
true, facetLocalDofIdx);
177 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false, facetLocalDofIdx);
180 if (!scvf.boundary())
182 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
184 const auto outsideGridScvfIdx = flipScvfIdxSet[i-1];
185 const auto& flipScvf = fvGeometry.scvf(outsideGridScvfIdx);
186 const auto& outsideFlipScvfIdxSet = fvGeometry.gridGeometry().flipScvfIndexSet()[outsideGridScvfIdx];
190 auto outsideNeighborScvIdxSet = neighborScvIndicesLocal;
191 outsideNeighborScvIdxSet[0] = outsideNeighborScvIdxSet[i];
192 for (LocalIndexType j = 0; j < outsideFlipScvfIdxSet.size(); ++j)
194 const auto flipScvfIdx = outsideFlipScvfIdxSet[j];
195 auto it = std::find(flipScvfIdxSet.begin(), flipScvfIdxSet.end(), flipScvfIdx);
198 if (it != flipScvfIdxSet.end())
199 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[
std::distance(flipScvfIdxSet.begin(), it)+1];
204 assert(flipScvfIdx == gridScvfIdx);
205 outsideNeighborScvIdxSet[j+1] = neighborScvIndicesLocal[0];
209 scvfIndexMap[outsideGridScvfIdx] = curLocalScvfIdx+i;
210 localFaceData_.emplace_back(curLocalScvfIdx+i, outsideNeighborScvIdxSet[0], flipScvf.index());
212 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, facetLocalDofIdx,
true, facetLocalDofIdx);
214 scvfs_.emplace_back(flipScvf, outsideNeighborScvIdxSet, numUnknowns_++,
false, facetLocalDofIdx);
220 else if (scvf.boundary())
222 if (problem.boundaryTypes(
element, scvf).hasOnlyDirichlet())
224 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numKnowns_++,
true);
225 dirichletData_.emplace_back(scvf.outsideScvIdx());
228 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false);
234 scvfs_.emplace_back(scvf, neighborScvIndicesLocal, numUnknowns_++,
false);
237 for (LocalIndexType i = 1; i < numNeighborScvs; ++i)
240 const auto outsideLocalScvIdx = neighborScvIndicesLocal[i];
241 const auto& flipScvfIndex = fvGeometry.gridGeometry().flipScvfIndexSet()[scvf.index()][i-1];
242 const auto& flipScvf = fvGeometry.scvf(flipScvfIndex);
243 scvfIndexMap[flipScvfIndex] = curLocalScvfIdx;
244 localFaceData_.emplace_back(curLocalScvfIdx,
253 for (LocalIndexType scvIdxLocal = 0; scvIdxLocal < numLocalScvs; scvIdxLocal++)
255 elements_.emplace_back(fvGeometry.gridGeometry().element(
stencil()[scvIdxLocal] ));
256 scvs_.emplace_back(fvGeometry.gridGeometry().mpfaHelper(),
258 fvGeometry.scv(
stencil()[scvIdxLocal] ),
267 {
return numFaces_; }
271 {
return numUnknowns_; }
275 {
return numKnowns_; }
279 {
return scvs_.size(); }
283 {
return *stencil_; }
286 const Element&
element(LocalIndexType ivLocalScvIdx)
const
287 {
return elements_[ivLocalScvIdx]; }
290 const LocalScvfType&
localScvf(LocalIndexType ivLocalScvfIdx)
const
291 {
return scvfs_[ivLocalScvfIdx]; }
294 const LocalScvType&
localScv(LocalIndexType ivLocalScvIdx)
const
295 {
return scvs_[ivLocalScvIdx]; }
299 {
return localFaceData_; }
303 {
return dirichletData_; }
307 {
return interiorBoundaryData_; }
310 template<
class FVElementGeometry >
311 auto getScvGeometry(LocalIndexType ivLocalScvIdx,
const FVElementGeometry& fvGeometry)
const
321 template<
class IvIndexSetContainer,
324 class FlipScvfIndexSet >
326 ScvfIndexMap& scvfIndexMap,
327 const NodalIndexSet& nodalIndexSet,
328 const FlipScvfIndexSet& flipScvfIndexSet)
339 const Stencil* stencil_;
342 std::vector<Element> elements_;
343 std::vector<LocalScvType> scvs_;
344 std::vector<LocalScvfType> scvfs_;
345 std::vector<LocalFaceData> localFaceData_;
346 std::vector<DirichletData> dirichletData_;
347 std::vector<InteriorBoundaryData> interiorBoundaryData_;
350 std::size_t numFaces_;
351 std::size_t numUnknowns_;
352 std::size_t numKnowns_;
Base class for the interaction volumes of mpfa methods. It defines the interface and actual implement...
Definition: interactionvolumebase.hh:56
Define data structure to store which scvfs lie on interior boundaries.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:88
GridIndexType scvfIndex() const
Return corresponding scvf index.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:96
InteriorBoundaryData(GridIndexType scvfIdx)
Constructor.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:93
Forward declaration of the facet coupling MPFA-O interaction volume.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:69
std::size_t numScvs() const
returns the number of scvs embedded in this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:278
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:290
const std::vector< LocalFaceData > & localFaceData() const
returns a reference to the container with the local face data
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:298
auto getScvGeometry(LocalIndexType ivLocalScvIdx, const FVElementGeometry &fvGeometry) const
returns the geometry of the i-th local scv
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:311
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:302
typename CCMpfaOInteractionVolume< Traits >::DirichletData DirichletData
Reuse standard o-scheme's Dirichlet Data class.
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:84
const std::vector< InteriorBoundaryData > & interiorBoundaryData() const
returns a reference to the data container on interior boundaries
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:306
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:274
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:104
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:286
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:316
std::size_t numUnknowns() const
returns the number of intermediate unknowns within this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:270
const Stencil & stencil() const
returns the cell-stencil of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:282
std::size_t numFaces() const
returns the number of primary scvfs of this interaction volume
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:266
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:294
static void addIVIndexSets(IvIndexSetContainer &ivIndexSetContainer, ScvfIndexMap &scvfIndexMap, const NodalIndexSet &nodalIndexSet, const FlipScvfIndexSet &flipScvfIndexSet)
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:325
static constexpr MpfaMethods MpfaMethod
publicly state the mpfa-scheme this interaction volume is associated with
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:100
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:35
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:120
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
Specialization of the interaction volume-local assembler class for the schemes using an mpfa-o type a...
Definition: multidomain/facet/cellcentered/mpfa/localassembler.hh:43
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 interaction volume of the mpfa-o scheme.
MpfaMethods
The available mpfa schemes in Dumux.
Definition: methods.hh:22
static ctype distance(const Dune::FieldVector< ctype, dimWorld > &a, const Dune::FieldVector< ctype, dimWorld > &b)
Compute the shortest distance between two points.
Definition: distance.hh:282
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.
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. This uses dynamic types types for ...
Definition: discretization/cellcentered/mpfa/omethod/interactionvolume.hh:48
The default interaction volume traits class for the mpfa-o method in the context of facet coupling....
Definition: multidomain/facet/cellcentered/mpfa/interactionvolume.hh:41
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: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