24#ifndef DUMUX_BOX_FACETCOUPLING_MAPPER_HH
25#define DUMUX_BOX_FACETCOUPLING_MAPPER_HH
27#include <dune/common/indices.hh>
48template<
class BulkFVG,
class LowDimFVG, std::
size_t bulkId, std::
size_t lowDimId>
53 using LowDimElement =
typename LowDimFVG::GridView::template Codim<0>::Entity;
56 using BulkGridView =
typename BulkFVG::GridView;
57 using LowDimGridView =
typename LowDimFVG::GridView;
58 static constexpr int bulkDim = BulkGridView::dimension;
59 static constexpr int lowDimDim = LowDimGridView::dimension;
63 using ParentType::bulkGridId;
64 using ParentType::facetGridId;
74 template<
class Embeddings >
75 void update(
const BulkFVG& bulkFvGridGeometry,
76 const LowDimFVG& lowDimFvGridGeometry,
77 std::shared_ptr<const Embeddings> embeddings)
81 update(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, GridAdapter(embeddings));
92 template<
class Embeddings,
class CodimOneGr
idAdapter >
93 void update(
const BulkFVG& bulkFvGridGeometry,
94 const LowDimFVG& lowDimFvGridGeometry,
95 std::shared_ptr<const Embeddings> embeddings,
99 auto addCouplingEntryPolicy = [&] (
auto&& adjoinedEntityIndices,
100 const LowDimElement& lowDimElement,
101 const LowDimFVG& lowDimFvGridGeometry,
102 const BulkFVG& bulkFvGridGeometry)
107 const auto lowDimElemIdx = lowDimFvGridGeometry.elementMapper().index(lowDimElement);
108 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
111 const auto& eg = lowDimElement.geometry();
112 const auto numElementCorners = eg.corners();
113 std::vector<BulkIndexType> elementCorners(numElementCorners);
114 for (
int i = 0; i < numElementCorners; ++i)
115 elementCorners[i] = codimOneGridAdapter.
bulkGridVertexIndex(lowDimElement.template subEntity<lowDimDim>(i));
118 const auto unsortedElemCorners = elementCorners;
119 std::sort(elementCorners.begin(), elementCorners.end());
121 auto fvGeometry =
localView(bulkFvGridGeometry);
122 for (
auto bulkElemIdx : adjoinedEntityIndices)
124 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
125 const auto bulkRefElem = referenceElement(bulkElement);
128 const auto coupledFacetIndex = [&]
131 unsigned int coupledFacetIndex = 0;
132 std::vector<unsigned int> handledFacets;
133 for (
const auto& is : intersections(bulkFvGridGeometry.gridView(), bulkElement))
136 if (std::count(handledFacets.begin(), handledFacets.end(), is.indexInInside()))
139 handledFacets.push_back(is.indexInInside());
142 const auto numCorners = is.geometry().corners();
143 std::vector<BulkIndexType> facetIndices(
numCorners);
146 const auto vIdxLocal = bulkRefElem.subEntity(is.indexInInside(), 1, i, bulkDim);
147 facetIndices[i] = bulkFvGridGeometry.vertexMapper().vertexIndex(bulkElement.template subEntity<bulkDim>(vIdxLocal));
150 std::sort(facetIndices.begin(), facetIndices.end());
151 if ( std::equal(facetIndices.begin(), facetIndices.end(), elementCorners.begin(), elementCorners.end()) )
153 coupledFacetIndex = is.indexInInside();
160 DUNE_THROW(Dune::InvalidStateException,
"Could not find the bulk element coupling facet!");
162 return coupledFacetIndex;
166 fvGeometry.bindElement(bulkElement);
168 unsigned int foundCounter = 0;
169 std::vector<BulkIndexType> embeddedScvfIndices(numElementCorners);
170 for (
const auto& scvf : scvfs(fvGeometry))
172 if (scvf.interiorBoundary() && scvf.facetIndexInElement() == coupledFacetIndex)
177 const auto vIdxLocal = bulkRefElem.subEntity(coupledFacetIndex, 1, scvf.indexInElementFacet(), bulkDim);
178 const auto vIdxGlobal = bulkFvGridGeometry.vertexMapper().vertexIndex(bulkElement, vIdxLocal, bulkDim);
179 const auto it = std::find(unsortedElemCorners.begin(), unsortedElemCorners.end(), vIdxGlobal);
180 assert(it != unsortedElemCorners.end());
181 const auto lowDimElemLocalCornerIdx =
std::distance(unsortedElemCorners.begin(), it);
182 embeddedScvfIndices[lowDimElemLocalCornerIdx] = scvf.index();
188 if (foundCounter != numElementCorners)
189 DUNE_THROW(Dune::InvalidStateException,
"Found " << foundCounter <<
" instead of " << numElementCorners <<
" coupling scvfs in the bulk element");
192 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
194 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
195 : std::vector<LowDimIndexType>({lowDimElemIdx});
197 for (
auto lowDimDofIdx : lowDimElementDofs)
199 bulkData.couplingStencil.push_back(lowDimDofIdx);
200 auto& couplingScvfs = bulkData.dofToCouplingScvfMap[lowDimDofIdx];
201 couplingScvfs.insert(couplingScvfs.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
205 bulkData.couplingElementStencil.push_back(lowDimElemIdx);
206 auto& elemToScvfMap = bulkData.elementToScvfMap[lowDimElemIdx];
207 elemToScvfMap.insert(elemToScvfMap.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
210 lowDimData.embedments.emplace_back(bulkElemIdx, std::move(embeddedScvfIndices));
212 const auto bulkElementDofs = this->extractNodalDofs_(bulkElement, bulkFvGridGeometry);
213 for (
auto bulkDofIdx : bulkElementDofs)
214 lowDimData.couplingStencil.push_back(bulkDofIdx);
219 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
222 auto makeStencilUnique = [] (
auto& data)
224 auto& cs = data.second.couplingStencil;
225 std::sort(cs.begin(), cs.end());
226 cs.erase( std::unique(cs.begin(), cs.end()), cs.end() );
229 auto& lowDimCouplingData = this->couplingMap_(facetGridId, bulkGridId);
230 std::for_each(lowDimCouplingData.begin(), lowDimCouplingData.end(), makeStencilUnique);
235 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
236 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeStencilUnique);
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
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:294
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr Box box
Definition: method.hh:136
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:83
std::size_t numCorners(Shape shape)
Returns the number of corners of a given geometry.
Definition: throatproperties.hh:232
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:39
void update(const BulkFVG &bulkFvGridGeometry, const LowDimFVG &lowDimFvGridGeometry, std::shared_ptr< const Embeddings > embeddings, const CodimOneGridAdapter &codimOneGridAdapter)
Update coupling maps with a given grid adapter.
Definition: facet/box/couplingmapper.hh:93
void update(const BulkFVG &bulkFvGridGeometry, const LowDimFVG &lowDimFvGridGeometry, std::shared_ptr< const Embeddings > embeddings)
Update coupling maps. This is the standard interface required by any mapper implementation.
Definition: facet/box/couplingmapper.hh:75
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:52
BulkIndexType bulkGridVertexIndex(const FacetGridVertex &v) const
Returns the index within the d-dimensional grid of a vertex of the (d-1)-dimensional grid.
Definition: codimonegridadapter.hh:145
Implementation for the coupling mapper that sets up and stores the coupling maps between two domains ...
Definition: facet/couplingmapper.hh:54
Base class for the coupling mapper that sets up and stores the coupling maps between two domains of d...
Definition: couplingmapperbase.hh:53