3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
Loading...
Searching...
No Matches
facet/box/couplingmapper.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
22#ifndef DUMUX_BOX_FACETCOUPLING_MAPPER_HH
23#define DUMUX_BOX_FACETCOUPLING_MAPPER_HH
24
25#include <dune/common/indices.hh>
26#include <dune/geometry/referenceelements.hh>
27
33
34namespace Dumux {
35
47template<class BulkFVG, class LowDimFVG, std::size_t bulkId, std::size_t lowDimId>
48class FacetCouplingMapper<BulkFVG, LowDimFVG, bulkId, lowDimId, DiscretizationMethod::box>
49: public virtual FacetCouplingMapperBase<BulkFVG, LowDimFVG, bulkId, lowDimId>
50{
52 using LowDimElement = typename LowDimFVG::GridView::template Codim<0>::Entity;
53
54 // dimensions of the two grids
55 using BulkGridView = typename BulkFVG::GridView;
56 using LowDimGridView = typename LowDimFVG::GridView;
57 static constexpr int bulkDim = BulkGridView::dimension;
58 static constexpr int lowDimDim = LowDimGridView::dimension;
59
60public:
64
73 template< class Embeddings >
74 void update(const BulkFVG& bulkFvGridGeometry,
75 const LowDimFVG& lowDimFvGridGeometry,
76 std::shared_ptr<const Embeddings> embeddings)
77 {
78 // forward to update function with instantiated vertex adapter
80 update(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, GridAdapter(embeddings));
81 }
82
91 template< class Embeddings, class CodimOneGridAdapter >
92 void update(const BulkFVG& bulkFvGridGeometry,
93 const LowDimFVG& lowDimFvGridGeometry,
94 std::shared_ptr<const Embeddings> embeddings,
95 const CodimOneGridAdapter& codimOneGridAdapter)
96 {
97 // define the policy how to add map entries for given lowdim element and adjoined entity indices
98 auto addCouplingEntryPolicy = [&] (auto&& adjoinedEntityIndices,
99 const LowDimElement& lowDimElement,
100 const LowDimFVG& lowDimFvGridGeometry,
101 const BulkFVG& bulkFvGridGeometry)
102 {
103 using LowDimIndexType = typename IndexTraits<LowDimGridView>::GridIndex;
104 using BulkIndexType = typename IndexTraits<BulkGridView>::GridIndex;
105 using BulkReferenceElements = Dune::ReferenceElements<typename BulkGridView::ctype, bulkDim>;
106
107 const auto lowDimElemIdx = lowDimFvGridGeometry.elementMapper().index(lowDimElement);
108 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
109
110 // determine corner indices (in bulk grid indices)
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));
116
117 // save unsorted set of corner indices and search scvfs in adjoined entities
118 const auto unsortedElemCorners = elementCorners;
119 std::sort(elementCorners.begin(), elementCorners.end());
120 for (auto bulkElemIdx : adjoinedEntityIndices)
121 {
122 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
123 const auto bulkRefElem = BulkReferenceElements::general(bulkElement.type());
124
125 // find the bulk element facet that lies on this low dim element (assumes conformity!)
126 bool found = false;
127 unsigned int coupledFacetIndex;
128 std::vector<unsigned int> handledFacets;
129 for (const auto& is : intersections(bulkFvGridGeometry.gridView(), bulkElement))
130 {
131 // skip already handled facets (necessary for e.g. Dune::FoamGrid)
132 if (std::count(handledFacets.begin(), handledFacets.end(), is.indexInInside()))
133 continue;
134
135 handledFacets.push_back(is.indexInInside());
136
137 // determine if it lies on low dim element by comparing corner indices
138 const auto numCorners = is.geometry().corners();
139 std::vector<BulkIndexType> facetIndices(numCorners);
140 for (int i = 0; i < numCorners; ++i)
141 {
142 const auto vIdxLocal = bulkRefElem.subEntity(is.indexInInside(), 1, i, bulkDim);
143 facetIndices[i] = bulkFvGridGeometry.vertexMapper().vertexIndex(bulkElement.template subEntity<bulkDim>(vIdxLocal));
144 }
145
146 std::sort(facetIndices.begin(), facetIndices.end());
147 if ( std::equal(facetIndices.begin(), facetIndices.end(), elementCorners.begin(), elementCorners.end()) )
148 {
149 coupledFacetIndex = is.indexInInside();
150 found = true; break;
151 }
152 }
153
154 // ensure that we found the facet!
155 if (!found)
156 DUNE_THROW(Dune::InvalidStateException, "Could not find the bulk element coupling facet!");
157
158 // we should always find numElementCorners coupling scvfs
159 auto fvGeometry = localView(bulkFvGridGeometry);
160 fvGeometry.bindElement(bulkElement);
161
162 unsigned int foundCounter = 0;
163 std::vector<BulkIndexType> embeddedScvfIndices(numElementCorners);
164 for (const auto& scvf : scvfs(fvGeometry))
165 {
166 if (scvf.interiorBoundary() && scvf.facetIndexInElement() == coupledFacetIndex)
167 {
168 // we want to order the set of scvfs lying on the lower-dimensional element such that the i-th scvf
169 // coincides with the i-th low dim element corner. This will help later to identify which scvf fluxes
170 // enter which scv of the low dim element if the lower-dimensional domain uses the box scheme
171 const auto vIdxLocal = bulkRefElem.subEntity(coupledFacetIndex, 1, scvf.indexInElementFacet(), bulkDim);
172 const auto vIdxGlobal = bulkFvGridGeometry.vertexMapper().vertexIndex(bulkElement, vIdxLocal, bulkDim);
173 const auto it = std::find(unsortedElemCorners.begin(), unsortedElemCorners.end(), vIdxGlobal);
174 assert(it != unsortedElemCorners.end());
175 const auto lowDimElemLocalCornerIdx = std::distance(unsortedElemCorners.begin(), it);
176 embeddedScvfIndices[lowDimElemLocalCornerIdx] = scvf.index();
177 foundCounter++;
178 }
179 }
180
181 // ensure we found all scvfs
182 if (foundCounter != numElementCorners)
183 DUNE_THROW(Dune::InvalidStateException, "Found " << foundCounter << " instead of " << numElementCorners << " coupling scvfs in the bulk element");
184
185 // add each dof in the low dim element to coupling stencil of the bulk element
186 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
187 const auto lowDimElementDofs = LowDimFVG::discMethod == DiscretizationMethod::box
188 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
189 : std::vector<LowDimIndexType>({lowDimElemIdx});
190
191 for (auto lowDimDofIdx : lowDimElementDofs)
192 {
193 bulkData.couplingStencil.push_back(lowDimDofIdx);
194 auto& couplingScvfs = bulkData.dofToCouplingScvfMap[lowDimDofIdx];
195 couplingScvfs.insert(couplingScvfs.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
196 }
197
198 // add info on which scvfs coincide with which low dim element
199 bulkData.couplingElementStencil.push_back(lowDimElemIdx);
200 auto& elemToScvfMap = bulkData.elementToScvfMap[lowDimElemIdx];
201 elemToScvfMap.insert(elemToScvfMap.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
202
203 // add embedment and the bulk cell dofs to coupling stencil of low dim element
204 lowDimData.embedments.emplace_back(bulkElemIdx, std::move(embeddedScvfIndices));
205
206 const auto bulkElementDofs = this->extractNodalDofs_(bulkElement, bulkFvGridGeometry);
207 for (auto bulkDofIdx : bulkElementDofs)
208 lowDimData.couplingStencil.push_back(bulkDofIdx);
209 }
210 };
211
212 // let the parent do the update subject to the execution policy defined above
213 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
214
215 // coupling stencils might not be unique with the policy above
216 auto makeStencilUnique = [] (auto& data)
217 {
218 auto& cs = data.second.couplingStencil;
219 std::sort(cs.begin(), cs.end());
220 cs.erase( std::unique(cs.begin(), cs.end()), cs.end() );
221 };
222
223 auto& lowDimCouplingData = this->couplingMap_(facetGridId, bulkGridId);
224 std::for_each(lowDimCouplingData.begin(), lowDimCouplingData.end(), makeStencilUnique);
225
226 // bulk coupling stencil is only non-unique if box is used
227 if (LowDimFVG::discMethod == DiscretizationMethod::box)
228 {
229 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
230 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeStencilUnique);
231 }
232 }
233};
234
235} // end namespace Dumux
236
237#endif // DUMUX_FACETCOUPLING_CCTPFA_COUPLING_MAPPER_HH
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
Base class for the coupling mapper that sets up and stores the coupling maps between two domains of d...
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:38
DiscretizationMethod
The available discretization methods in Dumux.
Definition method.hh:37
@ box
Definition method.hh:38
Definition adapt.hh:29
typename GridView::IndexSet::IndexType GridIndex
Definition indextraits.hh:39
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:74
static constexpr auto bulkGridId
export domain ids
Definition couplingmapperbase.hh:129
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:92
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition codimonegridadapter.hh:53
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:147
Implementation for the coupling mapper that sets up and stores the coupling maps between two domains ...
Definition facet/couplingmapper.hh:52
Base class for the coupling mapper that sets up and stores the coupling maps between two domains of d...
Definition couplingmapperbase.hh:51
static constexpr auto facetGridId
Definition couplingmapperbase.hh:130
std::vector< typename IndexTraits< typename GridGeometry::GridView >::GridIndex > extractNodalDofs_(const typename GridGeometry::GridView::template Codim< 0 >::Entity &element, const GridGeometry &gridGeometry)
Creates a container with the nodal dofs within an element.
Definition couplingmapperbase.hh:210
void update_(const BulkFVG &bulkFvGridGeometry, const LowDimFVG &lowDimFvGridGeometry, std::shared_ptr< const Embeddings > embeddings, AddCouplingEntryPolicy &&addCouplingEntryPolicy)
Update coupling maps.
Definition couplingmapperbase.hh:175
static constexpr auto bulkGridId
Export grid ids.
Definition couplingmapperbase.hh:129
BulkCouplingMap & couplingMap_(GridIdType< bulkId >, GridIdType< lowDimId >)
returns non-const coupling data for bulk -> lowDim
Definition couplingmapperbase.hh:225
Implementation for the coupling mapper that sets up and stores the coupling maps between two domains ...