3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
boundary/freeflowporenetwork/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 *****************************************************************************/
25#ifndef DUMUX_MULTIDOMAIN_FREEFLOW_POROUSMEDIUM_COUPLINGMAPPER_HH
26#define DUMUX_MULTIDOMAIN_FREEFLOW_POROUSMEDIUM_COUPLINGMAPPER_HH
27
28#include <algorithm>
29#include <iostream>
30#include <map>
31#include <unordered_map>
32#include <tuple>
33#include <vector>
34
35#include <dune/common/timer.hh>
36#include <dune/common/exceptions.hh>
37#include <dune/common/indices.hh>
38#include <dune/common/reservedvector.hh>
39#include <dune/geometry/axisalignedcubegeometry.hh>
40
45
46namespace Dumux {
47
52// template<class MDTraits, class CouplingManager>
54{
55 using MapType = std::unordered_map<std::size_t, std::vector<std::size_t>>;
56
57public:
61 template<class FreeFlowMomentumGridGeometry, class FreeFlowMassGridGeometry, class PoreNetworkGridGeometry>
62 void update(const FreeFlowMomentumGridGeometry& ffMomentumGridGeometry,
63 const FreeFlowMassGridGeometry& ffMassGridGeometry,
64 const PoreNetworkGridGeometry& pnmGridGeometry)
65 {
66 clear_();
67 resize_(ffMomentumGridGeometry, pnmGridGeometry);
68 Dune::Timer watch;
69 std::cout << "Initializing the coupling map..." << std::endl;
70
71 auto ffFvGeometry = localView(ffMomentumGridGeometry);
72 auto pnmFvGeometry = localView(pnmGridGeometry);
73
74 using GlobalPosition = typename FreeFlowMomentumGridGeometry::GridView::template Codim<0>::Entity::Geometry::GlobalCoordinate;
75
76 for (const auto& pnmElement : elements(pnmGridGeometry.gridView()))
77 {
78 const auto pnmElementIdx = pnmGridGeometry.elementMapper().index(pnmElement);
79 pnmFvGeometry.bindElement(pnmElement);
80 for (const auto& pnmScv : scvs(pnmFvGeometry))
81 {
82 // skip the dof if it is not on the boundary
83 if (!pnmGridGeometry.dofOnBoundary(pnmScv.dofIndex()))
84 continue;
85
86 // get the intersection bulk element
87 const auto pnmPos = pnmScv.dofPosition();
88 const auto pnmDofIdx = pnmScv.dofIndex();
89 const auto& otherPNMScv = pnmFvGeometry.scv(1 - pnmScv.indexInElement());
90 const auto otherPNMScvDofIdx = otherPNMScv.dofIndex();
91
92 // check for intersections, skip if no intersection was found
93 const auto directlyCoupledFreeFlowElements = intersectingEntities(pnmPos, ffMomentumGridGeometry.boundingBoxTree());
94 if (directlyCoupledFreeFlowElements.empty())
95 continue;
96 else
97 isCoupledPNMDof_[pnmDofIdx] = true;
98
99 // determine the normal direction of the local coupling interface heuristically:
100 // find all element intersections touching the pore and take the normal which
101 // occurs most frequently
102 const std::size_t couplingNormalDirectionIndex = [&]
103 {
104 using Key = std::pair<std::size_t, bool>;
105 std::map<Key, std::size_t> result;
106 for (const auto eIdx : directlyCoupledFreeFlowElements)
107 {
108 for (const auto& intersection : intersections(ffMomentumGridGeometry.gridView(), ffMomentumGridGeometry.element(eIdx)))
109 {
110 if (intersectsPointGeometry(pnmPos, intersection.geometry()))
111 {
112 const auto& normal = intersection.centerUnitOuterNormal();
113 const auto normalAxis = Dumux::normalAxis(normal);
114 const Key key = std::make_pair(normalAxis, std::signbit(normal[normalAxis]));
115 ++result[key];
116 }
117 }
118 }
119
120 // TODO how to properly handle this corner (literally) case
121 if (directlyCoupledFreeFlowElements.size() == 1 && result.size() > 1)
122 DUNE_THROW(Dune::InvalidStateException, "Pore may not intersect with faces of different orientation when coupled to only one element");
123
124 return std::max_element(result.begin(), result.end(), [](const auto& x, const auto& y) { return x.second < y.second;})->first.first;
125 }();
126
127 using Scalar = typename FreeFlowMomentumGridGeometry::GridView::ctype;
128 const Scalar couplingPoreRadius = pnmGridGeometry.poreInscribedRadius(pnmDofIdx);
129 GlobalPosition lowerLeft = pnmPos - GlobalPosition(couplingPoreRadius);
130 lowerLeft[couplingNormalDirectionIndex] = pnmPos[couplingNormalDirectionIndex];
131 GlobalPosition upperRight = pnmPos + GlobalPosition(couplingPoreRadius);
132 upperRight[couplingNormalDirectionIndex] = pnmPos[couplingNormalDirectionIndex];
133
134 auto axes = std::move(std::bitset<FreeFlowMomentumGridGeometry::Grid::dimensionworld>{}.set());
135 axes.set(couplingNormalDirectionIndex, false);
136
137 using PoreIntersectionGeometryType = Dune::AxisAlignedCubeGeometry<Scalar,
138 FreeFlowMomentumGridGeometry::GridView::dimension-1,
139 FreeFlowMomentumGridGeometry::GridView::dimensionworld>;
140
141 PoreIntersectionGeometryType poreIntersectionGeometry(lowerLeft, upperRight, axes);
142 const auto allCoupledFreeFlowElements = intersectingEntities(std::move(poreIntersectionGeometry), ffMomentumGridGeometry.boundingBoxTree());
143
144
145 for (const auto& ffElementInfo : allCoupledFreeFlowElements)
146 {
147 const auto freeFlowElementIndex = ffElementInfo.second();
148 pnmElementToFreeFlowElementsMap_[pnmElementIdx].push_back(freeFlowElementIndex);
149 freeFlowElementToPNMElementMap_[freeFlowElementIndex] = pnmElementIdx;
150
151 pnmToFreeFlowMassStencils_[pnmElementIdx].push_back(freeFlowElementIndex);
152 freeFlowMassToPNMStencils_[freeFlowElementIndex].push_back(pnmDofIdx);
153
154 ffFvGeometry.bindElement(ffMomentumGridGeometry.element(freeFlowElementIndex));
155 const auto coupledFreeFlowMomentumDofIndices = coupledFFMomentumDofs_(ffFvGeometry, ffMassGridGeometry, pnmPos, couplingPoreRadius, couplingNormalDirectionIndex);
156
157 pnmToFreeFlowMomentumStencils_[pnmElementIdx].push_back(coupledFreeFlowMomentumDofIndices.coupledFrontalDof);
158 freeFlowMomentumToPNMStencils_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof].push_back(pnmDofIdx);
159 freeFlowMomentumToPNMStencils_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof].push_back(otherPNMScvDofIdx);
160
161 isCoupledFreeFlowMomentumDof_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof] = true;
162 isCoupledFreeFlowMomentumDofOnInterface_[coupledFreeFlowMomentumDofIndices.coupledFrontalDof] = true;
163
164 // treat the coupled ff dofs not directly on interface
165 for (const auto ffDofIdx : coupledFreeFlowMomentumDofIndices.coupledLateralDofs)
166 {
167 freeFlowMomentumToPNMStencils_[ffDofIdx].push_back(pnmDofIdx);
168 freeFlowMomentumToPNMStencils_[ffDofIdx].push_back(otherPNMScvDofIdx);
169 isCoupledFreeFlowMomentumDof_[ffDofIdx] = true;
170 }
171 }
172 }
173 }
174
175 std::cout << "took " << watch.elapsed() << " seconds." << std::endl;
176 }
177
183 const std::vector<std::size_t>& poreNetworkToFreeFlowMomentumCouplingStencil(const std::size_t eIdxI) const
184 {
186 return pnmToFreeFlowMomentumStencils_.at(eIdxI);
187 else
188 return emptyStencil_;
189 }
190
196 const std::vector<std::size_t>& poreNetworkToFreeFlowMassCouplingStencil(const std::size_t eIdxI) const
197 {
199 return pnmToFreeFlowMassStencils_.at(eIdxI);
200 else
201 return emptyStencil_;
202 }
203
209 const std::vector<std::size_t>& freeFlowMassToPoreNetworkCouplingStencil(const std::size_t eIdxI) const
210 {
211 if (isCoupledFreeFlowElement(eIdxI))
212 return freeFlowMassToPNMStencils_.at(eIdxI);
213 else
214 return emptyStencil_;
215 }
216
222 const std::vector<std::size_t>& freeFlowMomentumToPoreNetworkCouplingStencil(const std::size_t dofIndex) const
223 {
224 if (isCoupledFreeFlowMomentumDof(dofIndex))
225 return freeFlowMomentumToPNMStencils_.at(dofIndex);
226 else
227 return emptyStencil_;
228 }
229
233 bool isCoupledFreeFlowElement(std::size_t eIdx) const
234 {
235 return static_cast<bool>(freeFlowElementToPNMElementMap_.count(eIdx));
236 }
237
241 bool isCoupledFreeFlowMomentumDof(std::size_t dofIdx) const
242 {
243 return isCoupledFreeFlowMomentumDof_[dofIdx];
244 }
245
249 bool isCoupledPoreNetworkElement(std::size_t eIdx) const
250 {
251 return static_cast<bool>(pnmElementToFreeFlowElementsMap_.count(eIdx));
252 }
253
257 bool isCoupledPoreNetworkDof(std::size_t dofIdx) const
258 {
259 return isCoupledPNMDof_[dofIdx];
260 }
261
262 bool isCoupledFreeFlowMomentumScvf(std::size_t scvfIdx) const
263 {
264 return isCoupledFrontalFreeFlowMomentumScvf_.count(scvfIdx);
265 }
266
267 bool isCoupledFreeFlowMomentumLateralScvf(std::size_t scvfIdx) const
268 {
269 return isCoupledLateralFreeFlowMomentumScvf_.count(scvfIdx);
270 }
271
272 bool isCoupledFreeFlowMassScvf(std::size_t scvfIdx) const
273 {
274 return isCoupledFreeFlowMassScvf_.count(scvfIdx);
275 }
276
278 { return pnmElementToFreeFlowElementsMap_;}
279
281 { return freeFlowElementToPNMElementMap_; }
282
283private:
284
285 void clear_()
286 {
287 pnmElementToFreeFlowElementsMap_.clear();
288 freeFlowElementToPNMElementMap_.clear();
289 isCoupledPNMDof_.clear();
290 isCoupledFreeFlowMomentumDof_.clear();
291 isCoupledFreeFlowMomentumDofOnInterface_.clear();
292 pnmToFreeFlowMassStencils_.clear();
293 pnmToFreeFlowMomentumStencils_.clear();
294 freeFlowMassToPNMStencils_.clear();
295 freeFlowMomentumToPNMStencils_.clear();
296 }
297
298 template<class FreeFlowMomentumGridGeometry, class PoreNetworkGridGeometry>
299 void resize_(const FreeFlowMomentumGridGeometry& ffMomentumGridGeometry,
300 const PoreNetworkGridGeometry& pnmGridGeometry)
301 {
302 const auto numPNMDofs = pnmGridGeometry.numDofs();
303 const auto numFreeFlowMomentumDofs = ffMomentumGridGeometry.numDofs();
304 isCoupledPNMDof_.resize(numPNMDofs, false);
305 isCoupledFreeFlowMomentumDof_.resize(numFreeFlowMomentumDofs, false);
306 isCoupledFreeFlowMomentumDofOnInterface_.resize(numFreeFlowMomentumDofs, false);
307 }
308
310 template<class FVElementGeometry, class FreeFlowMassGridGeometry, class GlobalPosition, class Scalar>
311 auto coupledFFMomentumDofs_(const FVElementGeometry& fvGeometry,
312 const FreeFlowMassGridGeometry& ffMassGridGeometry,
313 const GlobalPosition& pnmPos,
314 const Scalar couplingPoreRadius,
315 const int couplingInterfaceDirectionIdx)
316 {
317
318 struct Result
319 {
320 Dune::ReservedVector<std::size_t, FVElementGeometry::maxNumElementScvs> coupledLateralDofs;
321 std::size_t coupledFrontalDof;
322 } result;
323
324
325 using std::abs;
326 for (const auto& scv : scvs(fvGeometry))
327 {
328 const Scalar eps = diameter(scv.geometry())*1e-6; // TODO
329 assert(eps < couplingPoreRadius);
330
331 if (scv.dofAxis() == couplingInterfaceDirectionIdx) // the free flow dofs that lie within the coupling interface
332 {
333 if (abs(scv.dofPosition()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
334 {
335 result.coupledFrontalDof = scv.dofIndex();
336
337 // treat scvfs
338 for (const auto& scvf : scvfs(fvGeometry, scv))
339 {
340 // add lateral faces "standing" on coupling interface
341 if (scvf.isLateral() && !scvf.boundary())
342 isCoupledLateralFreeFlowMomentumScvf_[scvf.index()] = true;
343 else if (scvf.isFrontal() && scvf.boundary()) // add face lying on interface
344 {
345 isCoupledFrontalFreeFlowMomentumScvf_[scvf.index()] = true;
346
347 const auto& element = ffMassGridGeometry.element(fvGeometry.elementIndex()); // this local variable is needed to prevent a memory error
348 const auto ffMassFVGeometry = localView(ffMassGridGeometry).bindElement(element);
349 for (const auto& ffMassScvf : scvfs(ffMassFVGeometry))
350 {
351 if (abs(ffMassScvf.center()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
352 isCoupledFreeFlowMassScvf_[ffMassScvf.index()] = true;
353 }
354 }
355 }
356 }
357 }
358 else // the free flow dofs perpendicular to the coupling interface
359 {
360 bool isCoupledDof = false;
361
362 for (int dimIdx = 0; dimIdx < GlobalPosition::dimension; ++dimIdx)
363 {
364 if (dimIdx == couplingInterfaceDirectionIdx || scv.boundary())
365 continue;
366
367 isCoupledDof = abs(scv.dofPosition()[dimIdx] - pnmPos[dimIdx]) < couplingPoreRadius + eps;
368 }
369
370 if (isCoupledDof)
371 {
372 result.coupledLateralDofs.push_back(scv.dofIndex());
373
374 // treat scvfs
375 for (const auto& scvf : scvfs(fvGeometry, scv))
376 {
377 if (scvf.isLateral() && scvf.boundary())
378 {
379 // add lateral scvfs lying on interface
380 if (abs(scvf.ipGlobal()[couplingInterfaceDirectionIdx] - pnmPos[couplingInterfaceDirectionIdx]) < eps)
381 isCoupledLateralFreeFlowMomentumScvf_[scvf.index()] = true;
382 }
383 }
384 }
385 }
386 }
387
388 return result;
389 }
390
391 std::vector<std::size_t> emptyStencil_;
392
393 std::unordered_map<std::size_t, std::vector<std::size_t>> pnmElementToFreeFlowElementsMap_;
394 std::unordered_map<std::size_t, std::size_t> freeFlowElementToPNMElementMap_;
395
396 std::vector<bool> isCoupledPNMDof_;
397 std::vector<bool> isCoupledFreeFlowMomentumDof_;
398 std::vector<bool> isCoupledFreeFlowMomentumDofOnInterface_;
399
400 std::unordered_map<std::size_t, bool> isCoupledLateralFreeFlowMomentumScvf_;
401 std::unordered_map<std::size_t, bool> isCoupledFrontalFreeFlowMomentumScvf_;
402 std::unordered_map<std::size_t, bool> isCoupledFreeFlowMassScvf_;
403
404 MapType pnmToFreeFlowMassStencils_;
405 MapType pnmToFreeFlowMomentumStencils_;
406 MapType freeFlowMassToPNMStencils_;
407 MapType freeFlowMomentumToPNMStencils_;
408};
409
410} // end namespace Dumux
411
412#endif
Returns the normal axis index of a unit vector (0 = x, 1 = y, 2 = z)
The available discretization methods in Dumux.
A function to compute a geometry's diameter, i.e. the longest distance between points of a geometry.
Algorithms that finds which geometric entites intersect.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
std::vector< std::size_t > intersectingEntities(const Dune::FieldVector< ctype, dimworld > &point, const BoundingBoxTree< EntitySet > &tree, bool isCartesianGrid=false)
Compute all intersections between entities and a point.
Definition: intersectingentities.hh:112
Geometry::ctype diameter(const Geometry &geo)
Computes the longest distance between points of a geometry.
Definition: diameter.hh:36
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: adapt.hh:29
static std::size_t normalAxis(const Vector &v)
Definition: normalaxis.hh:39
Coupling mapper for staggered free-flow and pore-network models.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:54
bool isCoupledFreeFlowMomentumLateralScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:267
const std::vector< std::size_t > & freeFlowMassToPoreNetworkCouplingStencil(const std::size_t eIdxI) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:209
void update(const FreeFlowMomentumGridGeometry &ffMomentumGridGeometry, const FreeFlowMassGridGeometry &ffMassGridGeometry, const PoreNetworkGridGeometry &pnmGridGeometry)
Main update routine.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:62
const std::vector< std::size_t > & freeFlowMomentumToPoreNetworkCouplingStencil(const std::size_t dofIndex) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:222
bool isCoupledFreeFlowMomentumScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:262
bool isCoupledPoreNetworkDof(std::size_t dofIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:257
bool isCoupledFreeFlowMassScvf(std::size_t scvfIdx) const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:272
const std::vector< std::size_t > & poreNetworkToFreeFlowMomentumCouplingStencil(const std::size_t eIdxI) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:183
const auto & freeFlowElementToPNMElementMap() const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:280
const auto & pnmElementToFreeFlowElementsMap() const
Definition: boundary/freeflowporenetwork/couplingmapper.hh:277
bool isCoupledFreeFlowMomentumDof(std::size_t dofIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:241
const std::vector< std::size_t > & poreNetworkToFreeFlowMassCouplingStencil(const std::size_t eIdxI) const
returns an iteratable container of all indices of degrees of freedom of domain j that couple with / i...
Definition: boundary/freeflowporenetwork/couplingmapper.hh:196
bool isCoupledFreeFlowElement(std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:233
bool isCoupledPoreNetworkElement(std::size_t eIdx) const
Return if an element residual with index eIdx of domain i is coupled to domain j.
Definition: boundary/freeflowporenetwork/couplingmapper.hh:249