3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
facet/cellcentered/tpfa/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_CCTPFA_FACETCOUPLING_MAPPER_HH
25#define DUMUX_CCTPFA_FACETCOUPLING_MAPPER_HH
26
27#include <dune/common/indices.hh>
28#include <dune/common/float_cmp.hh>
29
34
35namespace Dumux {
36
49template<class BulkFVG, class LowDimFVG, std::size_t bulkId, std::size_t lowDimId>
50class FacetCouplingMapper<BulkFVG, LowDimFVG, bulkId, lowDimId, DiscretizationMethods::CCTpfa>
51: public virtual FacetCouplingMapperBase<BulkFVG, LowDimFVG, bulkId, lowDimId>
52{
54 using LowDimElement = typename LowDimFVG::GridView::template Codim<0>::Entity;
55 using GlobalPosition = typename LowDimElement::Geometry::GlobalCoordinate;
56
57public:
59 using ParentType::bulkGridId;
60 using ParentType::facetGridId;
61
70 template< class Embeddings >
71 void update(const BulkFVG& bulkFvGridGeometry,
72 const LowDimFVG& lowDimFvGridGeometry,
73 std::shared_ptr<const Embeddings> embeddings)
74 {
75 // define the policy how to add map entries for given lowdim element and adjoined entity indices
76 auto addCouplingEntryPolicy = [&] (auto&& adjoinedEntityIndices,
77 const LowDimElement& lowDimElement,
78 const LowDimFVG& lowDimFvGridGeometry,
79 const BulkFVG& bulkFvGridGeometry)
80 {
81 using LowDimIndexType = typename IndexTraits<typename LowDimFVG::GridView>::GridIndex;
82 using BulkIndexType = typename IndexTraits<typename BulkFVG::GridView>::GridIndex;
83
84 const auto lowDimGeometry = lowDimElement.geometry();
85 const auto lowDimElemIdx = lowDimFvGridGeometry.elementMapper().index(lowDimElement);
86 auto& lowDimData = this->couplingMap_(facetGridId, bulkGridId)[lowDimElemIdx];
87
88 auto fvGeometry = localView(bulkFvGridGeometry);
89 // find the scvfs in the adjoined entities coinciding with the low dim element
90 // since the bulk domain uses tpfa, there is always only going to be one scvf
91 for (auto bulkElemIdx : adjoinedEntityIndices)
92 {
93 const auto bulkElement = bulkFvGridGeometry.element(bulkElemIdx);
94 fvGeometry.bindElement(bulkElement);
95
96 bool found = false;
97 BulkIndexType embeddedScvfIdx;
98 for (const auto& scvf : scvfs(fvGeometry))
99 {
100 // for non-boundary faces, it suffices to check if one
101 // of the outside scv indices is in the set of embedments
102 if (!scvf.boundary())
103 {
104 if ( std::find(adjoinedEntityIndices.begin(),
105 adjoinedEntityIndices.end(),
106 scvf.outsideScvIdx()) != adjoinedEntityIndices.end() )
107 {
108 embeddedScvfIdx = scvf.index();
109 found = true; break;
110 }
111 }
112
113 // otherwise, do float comparison of element and scvf center
114 else
115 {
116 const auto eps = lowDimGeometry.volume()*1e-8;
117 const auto diffVec = lowDimGeometry.center()-scvf.center();
118
119 if ( Dune::FloatCmp::eq<GlobalPosition, Dune::FloatCmp::CmpStyle::absolute>(diffVec, GlobalPosition(0.0), eps) )
120 {
121 embeddedScvfIdx = scvf.index();
122 found = true; break;
123 }
124 }
125 }
126
127 // Error tracking. The boundary scvf detection might has to be improved for very fine grids!?
128 if (!found)
129 DUNE_THROW(Dune::InvalidStateException, "Could not find coupling scvf in embedment");
130
131 // add each dof in the low dim element to coupling stencil of the bulk element
132 auto& bulkData = this->couplingMap_(bulkGridId, facetGridId)[bulkElemIdx];
133 const auto lowDimElementDofs = LowDimFVG::discMethod == DiscretizationMethods::box
134 ? this->extractNodalDofs_(lowDimElement, lowDimFvGridGeometry)
135 : std::vector<LowDimIndexType>( {lowDimElemIdx} );
136
137 for (auto dofIdx : lowDimElementDofs)
138 {
139 bulkData.couplingStencil.push_back( dofIdx );
140 bulkData.dofToCouplingScvfMap[dofIdx].push_back( embeddedScvfIdx );
141 }
142
143 // add info on which scvfs coincide with which low dim element
144 bulkData.couplingElementStencil.push_back(lowDimElemIdx);
145 bulkData.elementToScvfMap[lowDimElemIdx].push_back( embeddedScvfIdx );
146
147 // add embedment (coupling stencil will be done below)
148 lowDimData.embedments.emplace_back( bulkElemIdx, std::vector<BulkIndexType>({embeddedScvfIdx}) );
149 }
150
151 // adjoint entity indices = coupling stencil for tpfa
152 lowDimData.couplingStencil = std::move(adjoinedEntityIndices);
153 };
154
155 // let the parent do the update subject to the execution policy defined above
156 ParentType::update_(bulkFvGridGeometry, lowDimFvGridGeometry, embeddings, addCouplingEntryPolicy);
157
158 // coupling stencils might not be unique if box is used in lowdim domain
159 if (LowDimFVG::discMethod == DiscretizationMethods::box)
160 {
161 auto makeStencilUnique = [] (auto& data)
162 {
163 auto& cs = data.second.couplingStencil;
164 std::sort(cs.begin(), cs.end());
165 cs.erase( std::unique(cs.begin(), cs.end()), cs.end() );
166 };
167
168 auto& bulkCouplingData = this->couplingMap_(bulkGridId, facetGridId);
169 std::for_each(bulkCouplingData.begin(), bulkCouplingData.end(), makeStencilUnique);
170 }
171 }
172};
173
174} // end namespace Dumux
175
176#endif
Defines the index types used for grid and local indices.
The available discretization methods in Dumux.
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)
Update coupling maps. This is the standard interface required by any mapper implementation.
Definition: facet/cellcentered/tpfa/couplingmapper.hh:71
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