3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
facet/box/couplingmapper.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_BOX_FACETCOUPLING_MAPPER_HH
25#define DUMUX_BOX_FACETCOUPLING_MAPPER_HH
26
27#include <dune/common/indices.hh>
28
34
35namespace Dumux {
36
48template<class BulkFVG, class LowDimFVG, std::size_t bulkId, std::size_t lowDimId>
49class FacetCouplingMapper<BulkFVG, LowDimFVG, bulkId, lowDimId, DiscretizationMethods::Box>
50: public virtual FacetCouplingMapperBase<BulkFVG, LowDimFVG, bulkId, lowDimId>
51{
53 using LowDimElement = typename LowDimFVG::GridView::template Codim<0>::Entity;
54
55 // dimensions of the two grids
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;
60
61public:
63 using ParentType::bulkGridId;
64 using ParentType::facetGridId;
65
74 template< class Embeddings >
75 void update(const BulkFVG& bulkFvGridGeometry,
76 const LowDimFVG& lowDimFvGridGeometry,
77 std::shared_ptr<const Embeddings> embeddings)
78 {
79 // forward to update function with instantiated vertex adapter
81 update(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, GridAdapter(embeddings));
82 }
83
92 template< class Embeddings, class CodimOneGridAdapter >
93 void update(const BulkFVG& bulkFvGridGeometry,
94 const LowDimFVG& lowDimFvGridGeometry,
95 std::shared_ptr<const Embeddings> embeddings,
96 const CodimOneGridAdapter& codimOneGridAdapter)
97 {
98 // define the policy how to add map entries for given lowdim element and adjoined entity indices
99 auto addCouplingEntryPolicy = [&] (auto&& adjoinedEntityIndices,
100 const LowDimElement& lowDimElement,
101 const LowDimFVG& lowDimFvGridGeometry,
102 const BulkFVG& bulkFvGridGeometry)
103 {
104 using LowDimIndexType = typename IndexTraits<LowDimGridView>::GridIndex;
105 using BulkIndexType = typename IndexTraits<BulkGridView>::GridIndex;
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
121 auto fvGeometry = localView(bulkFvGridGeometry);
122 for (auto bulkElemIdx : adjoinedEntityIndices)
123 {
124 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
125 const auto bulkRefElem = referenceElement(bulkElement);
126
127 // find the bulk element facet that lies on this low dim element (assumes conformity!)
128 const auto coupledFacetIndex = [&]
129 {
130 bool found = false;
131 unsigned int coupledFacetIndex = 0;
132 std::vector<unsigned int> handledFacets;
133 for (const auto& is : intersections(bulkFvGridGeometry.gridView(), bulkElement))
134 {
135 // skip already handled facets (necessary for e.g. Dune::FoamGrid)
136 if (std::count(handledFacets.begin(), handledFacets.end(), is.indexInInside()))
137 continue;
138
139 handledFacets.push_back(is.indexInInside());
140
141 // determine if it lies on low dim element by comparing corner indices
142 const auto numCorners = is.geometry().corners();
143 std::vector<BulkIndexType> facetIndices(numCorners);
144 for (int i = 0; i < numCorners; ++i)
145 {
146 const auto vIdxLocal = bulkRefElem.subEntity(is.indexInInside(), 1, i, bulkDim);
147 facetIndices[i] = bulkFvGridGeometry.vertexMapper().vertexIndex(bulkElement.template subEntity<bulkDim>(vIdxLocal));
148 }
149
150 std::sort(facetIndices.begin(), facetIndices.end());
151 if ( std::equal(facetIndices.begin(), facetIndices.end(), elementCorners.begin(), elementCorners.end()) )
152 {
153 coupledFacetIndex = is.indexInInside();
154 found = true; break;
155 }
156 }
157
158 // ensure that we found the facet!
159 if (!found)
160 DUNE_THROW(Dune::InvalidStateException, "Could not find the bulk element coupling facet!");
161
162 return coupledFacetIndex;
163 }();
164
165 // we should always find numElementCorners coupling scvfs
166 fvGeometry.bindElement(bulkElement);
167
168 unsigned int foundCounter = 0;
169 std::vector<BulkIndexType> embeddedScvfIndices(numElementCorners);
170 for (const auto& scvf : scvfs(fvGeometry))
171 {
172 if (scvf.interiorBoundary() && scvf.facetIndexInElement() == coupledFacetIndex)
173 {
174 // we want to order the set of scvfs lying on the lower-dimensional element such that the i-th scvf
175 // coincides with the i-th low dim element corner. This will help later to identify which scvf fluxes
176 // enter which scv of the low dim element if the lower-dimensional domain uses the box scheme
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();
183 foundCounter++;
184 }
185 }
186
187 // ensure we found all scvfs
188 if (foundCounter != numElementCorners)
189 DUNE_THROW(Dune::InvalidStateException, "Found " << foundCounter << " instead of " << numElementCorners << " coupling scvfs in the bulk element");
190
191 // add each dof in the low dim element to coupling stencil of the bulk element
192 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
193 const auto lowDimElementDofs = LowDimFVG::discMethod == DiscretizationMethods::box
194 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
195 : std::vector<LowDimIndexType>({lowDimElemIdx});
196
197 for (auto lowDimDofIdx : lowDimElementDofs)
198 {
199 bulkData.couplingStencil.push_back(lowDimDofIdx);
200 auto& couplingScvfs = bulkData.dofToCouplingScvfMap[lowDimDofIdx];
201 couplingScvfs.insert(couplingScvfs.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
202 }
203
204 // add info on which scvfs coincide with which low dim element
205 bulkData.couplingElementStencil.push_back(lowDimElemIdx);
206 auto& elemToScvfMap = bulkData.elementToScvfMap[lowDimElemIdx];
207 elemToScvfMap.insert(elemToScvfMap.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
208
209 // add embedment and the bulk cell dofs to coupling stencil of low dim element
210 lowDimData.embedments.emplace_back(bulkElemIdx, std::move(embeddedScvfIndices));
211
212 const auto bulkElementDofs = this->extractNodalDofs_(bulkElement, bulkFvGridGeometry);
213 for (auto bulkDofIdx : bulkElementDofs)
214 lowDimData.couplingStencil.push_back(bulkDofIdx);
215 }
216 };
217
218 // let the parent do the update subject to the execution policy defined above
219 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
220
221 // coupling stencils might not be unique with the policy above
222 auto makeStencilUnique = [] (auto& data)
223 {
224 auto& cs = data.second.couplingStencil;
225 std::sort(cs.begin(), cs.end());
226 cs.erase( std::unique(cs.begin(), cs.end()), cs.end() );
227 };
228
229 auto& lowDimCouplingData = this->couplingMap_(facetGridId, bulkGridId);
230 std::for_each(lowDimCouplingData.begin(), lowDimCouplingData.end(), makeStencilUnique);
231
232 // bulk coupling stencil is only non-unique if box is used
233 if (LowDimFVG::discMethod == DiscretizationMethods::box)
234 {
235 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
236 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeStencilUnique);
237 }
238 }
239};
240
241} // end namespace Dumux
242
243#endif
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