3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * 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 2 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_CCMPFA_FACETCOUPLING_MAPPER_HH
25#define DUMUX_CCMPFA_FACETCOUPLING_MAPPER_HH
26
27#include <dune/common/indices.hh>
28#include <dune/common/float_cmp.hh>
29
35
36namespace Dumux {
37
50template<class BulkFVG, class LowDimFVG, std::size_t bulkId, std::size_t lowDimId>
51class FacetCouplingMapper<BulkFVG, LowDimFVG, bulkId, lowDimId, DiscretizationMethods::CCMpfa>
52: public virtual FacetCouplingMapperBase<BulkFVG, LowDimFVG, bulkId, lowDimId>
53{
55
56 using BulkFVElementGeometry = typename BulkFVG::LocalView;
57 using BulkSubControlVolumeFace = typename BulkFVG::SubControlVolumeFace;
58
59 using BulkIndexType = typename IndexTraits<typename BulkFVG::GridView>::GridIndex;
60 using LowDimIndexType = typename IndexTraits<typename LowDimFVG::GridView>::GridIndex;
61
62 using LowDimElement = typename LowDimFVG::GridView::template Codim<0>::Entity;
63 using GlobalPosition = typename LowDimElement::Geometry::GlobalCoordinate;
64 static constexpr int lowDimDim = LowDimFVG::GridView::dimension;
65
66public:
68 using ParentType::bulkGridId;
69 using ParentType::facetGridId;
70
79 template< class Embeddings >
80 void update(const BulkFVG& bulkFvGridGeometry,
81 const LowDimFVG& lowDimFvGridGeometry,
82 std::shared_ptr<const Embeddings> embeddings)
83 {
84 // forward to update function with instantiated vertex adapter
86 update(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, GridAdapter(embeddings));
87 }
88
98 template< class Embeddings, class CodimOneGridAdapter >
99 void update(const BulkFVG& bulkFvGridGeometry,
100 const LowDimFVG& lowDimFvGridGeometry,
101 std::shared_ptr<const Embeddings> embeddings,
102 const CodimOneGridAdapter& codimOneGridAdapter)
103 {
104 // define the policy how to add map entries for given lowdim element and adjoined entity indices
105 auto addCouplingEntryPolicy = [&] (auto&& adjoinedEntityIndices,
106 const LowDimElement& lowDimElement,
107 const LowDimFVG& lowDimFvGridGeometry,
108 const BulkFVG& bulkFvGridGeometry)
109 {
110 const auto lowDimElemIdx = lowDimFvGridGeometry.elementMapper().index(lowDimElement);
111 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
112
113 // determine corner indices (in bulk grid indices)
114 const auto& lowDimGeometry = lowDimElement.geometry();
115 const auto numElementCorners = lowDimElement.subEntities(lowDimDim);
116 std::vector<BulkIndexType> elemCornerIndices(numElementCorners);
117 for (int i = 0; i < numElementCorners; ++i)
118 elemCornerIndices[i] = codimOneGridAdapter.bulkGridVertexIndex(lowDimElement.template subEntity<lowDimDim>(i));
119
120 // Find the scvfs in the adjoined entities coinciding with the low dim element
121 auto fvGeometry = localView(bulkFvGridGeometry);
122 for (auto bulkElemIdx : adjoinedEntityIndices)
123 {
124 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
125 fvGeometry.bind(bulkElement);
126
127 std::vector<BulkIndexType> embeddedScvfIndices;
128 for (const auto& scvf : scvfs(fvGeometry))
129 {
130 // for non-boundary faces, it suffices to check if one
131 // of the outside scv indices is in the set of embedments
132 if (!scvf.boundary())
133 {
134 if ( std::count(adjoinedEntityIndices.begin(), adjoinedEntityIndices.end(), scvf.outsideScvIdx()) )
135 embeddedScvfIndices.push_back(scvf.index());
136 }
137
138 // otherwise, do float comparison of element and scvf facet corner
139 else
140 {
141 const auto eps = lowDimGeometry.volume()*1e-8;
142 const auto diffVec = lowDimGeometry.center()-scvf.facetCorner();
143
144 if ( Dune::FloatCmp::eq<GlobalPosition, Dune::FloatCmp::CmpStyle::absolute>(diffVec, GlobalPosition(0.0), eps) )
145 embeddedScvfIndices.push_back(scvf.index());
146 }
147 }
148
149 // Error tracking. The boundary scvf detection might has to be improved for very fine grids!?
150 if ( embeddedScvfIndices.size() != numElementCorners )
151 DUNE_THROW(Dune::InvalidStateException, "Could not find all coupling scvfs in embedment");
152
153 // add each dof in the low dim element to coupling stencil of the bulk element and vice versa
154 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
155 const auto lowDimElementDofs = LowDimFVG::discMethod == DiscretizationMethods::box
156 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
157 : std::vector<LowDimIndexType>( {lowDimElemIdx} );
158
159 // add coupling info to all elements/scvfs in interaction volumes
160 for (auto dofIdx : lowDimElementDofs)
161 for (auto scvfIdx : embeddedScvfIndices)
162 this->addCouplingsFromIV_(bulkFvGridGeometry, fvGeometry.scvf(scvfIdx), fvGeometry, lowDimElemIdx, dofIdx);
163
164 // sort the scvfs according to the corners of the low dim element if box is used
165 // that allows identifying which scvf flux enters which low dim scv later
166 if (LowDimFVG::discMethod == DiscretizationMethods::box)
167 {
168 const auto copy = embeddedScvfIndices;
169
170 for (unsigned int i = 0; i < numElementCorners; ++i)
171 {
172 const auto& scvf = fvGeometry.scvf(copy[i]);
173 auto it = std::find(elemCornerIndices.begin(), elemCornerIndices.end(), scvf.vertexIndex());
174 assert(it != elemCornerIndices.end());
175 embeddedScvfIndices[ std::distance(elemCornerIndices.begin(), it) ] = copy[i];
176 }
177 }
178
179 // add info on which scvfs coincide with which low dim element
180 auto& elemToScvfMap = bulkData.elementToScvfMap[lowDimElemIdx];
181 elemToScvfMap.insert(elemToScvfMap.end(), embeddedScvfIndices.begin(), embeddedScvfIndices.end());
182
183 // add embedment
184 lowDimData.embedments.emplace_back( bulkElemIdx, std::move(embeddedScvfIndices) );
185 }
186 };
187
188 // let the parent do the update subject to the execution policy defined above
189 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
190
191 // lambda to make a container unique
192 auto makeUnique = [] (auto& c)
193 {
194 std::sort(c.begin(), c.end());
195 c.erase( std::unique(c.begin(), c.end()), c.end() );
196 };
197
198 // lambda to make bulk coupling map entry unique
199 auto makeBulkMapEntryUnique = [&makeUnique] (auto& data)
200 {
201 makeUnique(data.second.couplingStencil);
202 makeUnique(data.second.couplingElementStencil);
203 std::for_each(data.second.dofToCouplingScvfMap.begin(),
204 data.second.dofToCouplingScvfMap.end(),
205 [&makeUnique] (auto& pair) { makeUnique(pair.second); });
206 };
207
208 // make bulk map unique
209 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
210 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeBulkMapEntryUnique);
211
212 // coupling stencils might not be unique in lowdim domain
213 auto& lowDimCouplingData = this->couplingMap_(facetGridId, bulkGridId);
214 std::for_each(lowDimCouplingData.begin(),
215 lowDimCouplingData.end(),
216 [&makeUnique] (auto& pair) { makeUnique(pair.second.couplingStencil); });
217 }
218
219private:
220 void addCouplingsFromIV_(const BulkFVG& bulkFvGridGeometry,
221 const BulkSubControlVolumeFace& scvf,
222 const BulkFVElementGeometry& fvGeometry,
223 LowDimIndexType lowDimElemIdx,
224 LowDimIndexType lowDimDofIdx)
225 {
226 const auto& gridIvIndexSets = bulkFvGridGeometry.gridInteractionVolumeIndexSets();
227 if (bulkFvGridGeometry.vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
228 addCouplingsFromIVIndexSet_(gridIvIndexSets.secondaryIndexSet(scvf), fvGeometry, lowDimElemIdx, lowDimDofIdx);
229 else
230 addCouplingsFromIVIndexSet_(gridIvIndexSets.primaryIndexSet(scvf), fvGeometry, lowDimElemIdx, lowDimDofIdx);
231 }
232
233 template< class IVIndexSet >
234 void addCouplingsFromIVIndexSet_(const IVIndexSet& indexSet,
235 const BulkFVElementGeometry& fvGeometry,
236 LowDimIndexType lowDimElemIdx,
237 LowDimIndexType lowDimDofIdx)
238 {
239 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
240 for (auto bulkElemIdx : indexSet.gridScvIndices())
241 {
242 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
243
244 lowDimData.couplingStencil.push_back( bulkElemIdx );
245 bulkData.couplingStencil.push_back( lowDimDofIdx );
246 bulkData.couplingElementStencil.push_back( lowDimElemIdx );
247 }
248
249 // insert scvf couplings stemming from this interaction volume
250 for (auto scvfIdx : indexSet.gridScvfIndices())
251 {
252 const auto& scvf = fvGeometry.scvf(scvfIdx);
253 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[scvf.insideScvIdx()];
254 auto& couplingScvfs = bulkData.dofToCouplingScvfMap[lowDimDofIdx];
255 couplingScvfs.insert(couplingScvfs.end(), scvfIdx);
256 }
257 }
258};
259
260} // end namespace Dumux
261
262#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
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. This is the standard interface required by any mapper implementation.
Definition: facet/cellcentered/mpfa/couplingmapper.hh:99
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:80
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