version 3.10-dev
facet/cellcentered/mpfa/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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_CCMPFA_FACETCOUPLING_MAPPER_HH
13#define DUMUX_CCMPFA_FACETCOUPLING_MAPPER_HH
14
15#include <dune/common/indices.hh>
16#include <dune/common/float_cmp.hh>
17
23
24namespace Dumux {
25
38template<class BulkFVG, class LowDimFVG, std::size_t bulkId, std::size_t lowDimId>
39class FacetCouplingMapper<BulkFVG, LowDimFVG, bulkId, lowDimId, DiscretizationMethods::CCMpfa>
40: public virtual FacetCouplingMapperBase<BulkFVG, LowDimFVG, bulkId, lowDimId>
41{
43
44 using BulkFVElementGeometry = typename BulkFVG::LocalView;
45 using BulkSubControlVolumeFace = typename BulkFVG::SubControlVolumeFace;
46
47 using BulkIndexType = typename IndexTraits<typename BulkFVG::GridView>::GridIndex;
48 using LowDimIndexType = typename IndexTraits<typename LowDimFVG::GridView>::GridIndex;
49
50 using LowDimElement = typename LowDimFVG::GridView::template Codim<0>::Entity;
51 using GlobalPosition = typename LowDimElement::Geometry::GlobalCoordinate;
52 static constexpr int lowDimDim = LowDimFVG::GridView::dimension;
53
54public:
56 using ParentType::bulkGridId;
57 using ParentType::facetGridId;
58
67 template< class Embeddings >
68 void update(const BulkFVG& bulkFvGridGeometry,
69 const LowDimFVG& lowDimFvGridGeometry,
70 std::shared_ptr<const Embeddings> embeddings)
71 {
72 // forward to update function with instantiated vertex adapter
74 update(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, GridAdapter(embeddings));
75 }
76
86 template< class Embeddings, class CodimOneGridAdapter >
87 void update(const BulkFVG& bulkFvGridGeometry,
88 const LowDimFVG& lowDimFvGridGeometry,
89 std::shared_ptr<const Embeddings> embeddings,
90 const CodimOneGridAdapter& codimOneGridAdapter)
91 {
92 // define the policy how to add map entries for given lowdim element and adjoined entity indices
93 auto addCouplingEntryPolicy = [&] (auto&& adjoinedEntityIndices,
94 const LowDimElement& lowDimElement,
95 const LowDimFVG& lowDimFvGridGeometry,
96 const BulkFVG& bulkFvGridGeometry)
97 {
98 const auto lowDimElemIdx = lowDimFvGridGeometry.elementMapper().index(lowDimElement);
99 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
100
101 // determine corner indices (in bulk grid indices)
102 const auto& lowDimGeometry = lowDimElement.geometry();
103 const auto numElementCorners = lowDimElement.subEntities(lowDimDim);
104 std::vector<BulkIndexType> elemCornerIndices(numElementCorners);
105 for (int i = 0; i < numElementCorners; ++i)
106 elemCornerIndices[i] = codimOneGridAdapter.bulkGridVertexIndex(lowDimElement.template subEntity<lowDimDim>(i));
107
108 // Find the scvfs in the adjoined entities coinciding with the low dim element
109 auto fvGeometry = localView(bulkFvGridGeometry);
110 for (auto bulkElemIdx : adjoinedEntityIndices)
111 {
112 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
113 fvGeometry.bind(bulkElement);
114
115 std::vector<BulkIndexType> embeddedScvfIndices;
116 for (const auto& scvf : scvfs(fvGeometry))
117 {
118 // for non-boundary faces, it suffices to check if one
119 // of the outside scv indices is in the set of embedments
120 if (!scvf.boundary())
121 {
122 if ( std::count(adjoinedEntityIndices.begin(), adjoinedEntityIndices.end(), scvf.outsideScvIdx()) )
123 embeddedScvfIndices.push_back(scvf.index());
124 }
125
126 // otherwise, do float comparison of element and scvf facet corner
127 else
128 {
129 const auto eps = lowDimGeometry.volume()*1e-8;
130 const auto diffVec = lowDimGeometry.center()-fvGeometry.facetCorner(scvf);
131
132 if ( Dune::FloatCmp::eq<GlobalPosition, Dune::FloatCmp::CmpStyle::absolute>(diffVec, GlobalPosition(0.0), eps) )
133 embeddedScvfIndices.push_back(scvf.index());
134 }
135 }
136
137 // Error tracking. The boundary scvf detection might has to be improved for very fine grids!?
138 if ( embeddedScvfIndices.size() != numElementCorners )
139 DUNE_THROW(Dune::InvalidStateException, "Could not find all coupling scvfs in embedment");
140
141 // add each dof in the low dim element to coupling stencil of the bulk element and vice versa
142 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
143 const auto lowDimElementDofs = LowDimFVG::discMethod == DiscretizationMethods::box
144 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
145 : std::vector<LowDimIndexType>( {lowDimElemIdx} );
146
147 // add coupling info to all elements/scvfs in interaction volumes
148 for (auto dofIdx : lowDimElementDofs)
149 for (auto scvfIdx : embeddedScvfIndices)
150 this->addCouplingsFromIV_(bulkFvGridGeometry, fvGeometry.scvf(scvfIdx), fvGeometry, lowDimElemIdx, dofIdx);
151
152 // sort the scvfs according to the corners of the low dim element if box is used
153 // that allows identifying which scvf flux enters which low dim scv later
154 if (LowDimFVG::discMethod == DiscretizationMethods::box)
155 {
156 const auto copy = embeddedScvfIndices;
157
158 for (unsigned int i = 0; i < numElementCorners; ++i)
159 {
160 const auto& scvf = fvGeometry.scvf(copy[i]);
161 auto it = std::find(elemCornerIndices.begin(), elemCornerIndices.end(), scvf.vertexIndex());
162 assert(it != elemCornerIndices.end());
163 embeddedScvfIndices[ std::distance(elemCornerIndices.begin(), it) ] = copy[i];
164 }
165 }
166
167 // add info on which scvfs coincide with which low dim element
168 auto& elemToScvfMap = bulkData.elementToScvfMap[lowDimElemIdx];
169 elemToScvfMap.insert(elemToScvfMap.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
170
171 // add embedment
172 lowDimData.embedments.emplace_back( bulkElemIdx, std::move(embeddedScvfIndices) );
173 }
174 };
175
176 // let the parent do the update subject to the execution policy defined above
177 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
178
179 // lambda to make a container unique
180 auto makeUnique = [] (auto& c)
181 {
182 std::sort(c.begin(), c.end());
183 c.erase( std::unique(c.begin(), c.end()), c.end() );
184 };
185
186 // lambda to make bulk coupling map entry unique
187 auto makeBulkMapEntryUnique = [&makeUnique] (auto& data)
188 {
189 makeUnique(data.second.couplingStencil);
190 makeUnique(data.second.couplingElementStencil);
191 std::for_each(data.second.dofToCouplingScvfMap.begin(),
192 data.second.dofToCouplingScvfMap.end(),
193 [&makeUnique] (auto& pair) { makeUnique(pair.second); });
194 };
195
196 // make bulk map unique
197 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
198 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeBulkMapEntryUnique);
199
200 // coupling stencils might not be unique in lowdim domain
201 auto& lowDimCouplingData = this->couplingMap_(facetGridId, bulkGridId);
202 std::for_each(lowDimCouplingData.begin(),
203 lowDimCouplingData.end(),
204 [&makeUnique] (auto& pair) { makeUnique(pair.second.couplingStencil); });
205 }
206
207private:
208 void addCouplingsFromIV_(const BulkFVG& bulkFvGridGeometry,
209 const BulkSubControlVolumeFace& scvf,
210 const BulkFVElementGeometry& fvGeometry,
211 LowDimIndexType lowDimElemIdx,
212 LowDimIndexType lowDimDofIdx)
213 {
214 const auto& gridIvIndexSets = bulkFvGridGeometry.gridInteractionVolumeIndexSets();
215 if (bulkFvGridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
216 addCouplingsFromIVIndexSet_(gridIvIndexSets.secondaryIndexSet(scvf), fvGeometry, lowDimElemIdx, lowDimDofIdx);
217 else
218 addCouplingsFromIVIndexSet_(gridIvIndexSets.primaryIndexSet(scvf), fvGeometry, lowDimElemIdx, lowDimDofIdx);
219 }
220
221 template< class IVIndexSet >
222 void addCouplingsFromIVIndexSet_(const IVIndexSet& indexSet,
223 const BulkFVElementGeometry& fvGeometry,
224 LowDimIndexType lowDimElemIdx,
225 LowDimIndexType lowDimDofIdx)
226 {
227 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
228 for (auto bulkElemIdx : indexSet.gridScvIndices())
229 {
230 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
231
232 lowDimData.couplingStencil.push_back( bulkElemIdx );
233 bulkData.couplingStencil.push_back( lowDimDofIdx );
234 bulkData.couplingElementStencil.push_back( lowDimElemIdx );
235 }
236
237 // insert scvf couplings stemming from this interaction volume
238 for (auto scvfIdx : indexSet.gridScvfIndices())
239 {
240 const auto& scvf = fvGeometry.scvf(scvfIdx);
241 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[scvf.insideScvIdx()];
242 auto& couplingScvfs = bulkData.dofToCouplingScvfMap[lowDimDofIdx];
243 couplingScvfs.insert(couplingScvfs.end(), scvfIdx);
244 }
245 }
246};
247
248} // end namespace Dumux
249
250#endif
Adapter that allows retrieving information on a d-dimensional grid for entities of a (d-1)-dimensiona...
Definition: codimonegridadapter.hh:40
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:133
void update(const BulkFVG &bulkFvGridGeometry, const LowDimFVG &lowDimFvGridGeometry, std::shared_ptr< const Embeddings > embeddings, const CodimOneGridAdapter &codimOneGridAdapter)
Update coupling maps. This is the standard interface required by any mapper implementation.
Definition: facet/cellcentered/mpfa/couplingmapper.hh:87
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/cellcentered/mpfa/couplingmapper.hh:68
Base class for the coupling mapper that sets up and stores the coupling maps between two domains of d...
Definition: couplingmapperbase.hh:41
Implementation for the coupling mapper that sets up and stores the coupling maps between two domains ...
Definition: facet/couplingmapper.hh:42
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
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
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
typename GridView::IndexSet::IndexType GridIndex
Definition: indextraits.hh:27