version 3.8
intersectionentityset.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//
13#ifndef DUMUX_GEOMETRY_INTERSECTION_ENTITY_SET_HH
14#define DUMUX_GEOMETRY_INTERSECTION_ENTITY_SET_HH
15
16#include <algorithm>
17#include <cassert>
18#include <iostream>
19#include <memory>
20#include <type_traits>
21#include <utility>
22#include <vector>
23
24#include <dune/common/indices.hh>
25#include <dune/common/timer.hh>
26#include <dune/common/iteratorrange.hh>
27#include <dune/common/promotiontraits.hh>
28#include <dune/common/reservedvector.hh>
29#include <dune/geometry/affinegeometry.hh>
30#include <dune/geometry/type.hh>
31
34
35namespace Dumux {
36
41template<class DomainEntitySet, class TargetEntitySet>
43{
46
47 using ctype = typename Dune::PromotionTraits<typename DomainEntitySet::ctype, typename TargetEntitySet::ctype>::PromotedType;
48
49 static constexpr int dimWorld = DomainEntitySet::dimensionworld;
50 static_assert(dimWorld == int(TargetEntitySet::dimensionworld), "Entity sets must have the same world dimension");
51
52 using GlobalPosition = Dune::FieldVector<ctype, dimWorld>;
53
54 static constexpr int dimDomain = DomainEntitySet::Entity::Geometry::mydimension;
55 static constexpr int dimTarget = TargetEntitySet::Entity::Geometry::mydimension;
56 static constexpr bool isMixedDimensional = dimDomain != dimTarget;
57
61 class IntersectionEntity
62 {
63 static constexpr int dimIs = std::min(dimDomain, dimTarget);
64 using Geometry = Dune::AffineGeometry<ctype, dimIs, dimWorld>; // geometries are always simplices
65
66 // we can only have multiple neighbors in the mixeddimensional case and then only for the side with the largest dimension
67 using IndexStorage = std::pair<std::conditional_t<dimDomain <= dimTarget, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>,
68 std::conditional_t<dimTarget <= dimDomain, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>>;
69
70 static constexpr auto domainIdx = Dune::index_constant<0>{};
71 static constexpr auto targetIdx = Dune::index_constant<1>{};
72
73 public:
74 IntersectionEntity(const DomainTree& domainTree, const TargetTree& targetTree)
75 : domainTree_(domainTree)
76 , targetTree_(targetTree)
77 {}
78
80 void setCorners(const std::vector<GlobalPosition>& corners)
81 {
82 corners_ = corners;
83 assert(corners.size() == dimIs + 1); // Only simplex intersections are allowed
84 }
85
87 void addNeighbors(std::size_t domain, std::size_t target)
88 {
89 if (numDomainNeighbors() == 0 && numTargetNeighbors() == 0)
90 {
91 std::get<domainIdx>(neighbors_).push_back(domain);
92 std::get<targetIdx>(neighbors_).push_back(target);
93 }
94 else if (dimDomain > dimTarget)
95 std::get<domainIdx>(neighbors_).push_back(domain);
96
97 else if (dimTarget > dimDomain)
98 std::get<targetIdx>(neighbors_).push_back(target);
99
100 else
101 DUNE_THROW(Dune::InvalidStateException, "Cannot add more than one neighbor per side for equidimensional intersection!");
102 }
103
105 Geometry geometry() const
106 { return Geometry(Dune::GeometryTypes::simplex(dimIs), corners_); }
107
109 std::size_t numDomainNeighbors() const
110 { return std::get<domainIdx>(neighbors_).size(); }
111
113 std::size_t numTargetNeighbors() const
114 { return std::get<targetIdx>(neighbors_).size(); }
115
117 typename DomainEntitySet::Entity domainEntity(unsigned int n = 0) const
118 { return domainTree_.entitySet().entity(std::get<domainIdx>(neighbors_)[n]); }
119
121 typename TargetEntitySet::Entity targetEntity(unsigned int n = 0) const
122 { return targetTree_.entitySet().entity(std::get<targetIdx>(neighbors_)[n]); }
123
124 private:
125 IndexStorage neighbors_;
126 std::vector<GlobalPosition> corners_;
127
128 const DomainTree& domainTree_;
129 const TargetTree& targetTree_;
130 };
131
132 using Intersections = std::vector<IntersectionEntity>;
133
134public:
136 using Entity = IntersectionEntity;
138 using EntityIterator = typename Intersections::const_iterator;
139
144
148 void build(std::shared_ptr<const DomainEntitySet> domainSet, std::shared_ptr<const TargetEntitySet> targetSet)
149 {
150 domainTree_ = std::make_shared<DomainTree>(domainSet);
151 targetTree_ = std::make_shared<TargetTree>(targetSet);
152 build(*domainTree_, *targetTree_);
153 }
154
158 void build(std::shared_ptr<const DomainTree> domainTree, std::shared_ptr<const TargetTree> targetTree)
159 {
160 // make sure the tree don't get out of scope
161 domainTree_ = domainTree;
162 targetTree_ = targetTree;
163 build(*domainTree_, *targetTree_);
164 }
165
170 void build(const DomainTree& domainTree, const TargetTree& targetTree)
171 {
172 Dune::Timer timer;
173
174 // compute raw intersections
175 const auto rawIntersections = intersectingEntities(domainTree, targetTree);
176
177 // create a map to check whether intersection geometries were already inserted
178 // Note that this can only occur if the grids have different dimensionality.
179 // If this is the case, we keep track of the intersections using the indices of the lower-
180 // dimensional entities which is identical for all geometrically identical intersections.
181 std::vector<std::vector<std::vector<GlobalPosition>>> intersectionMap;
182 std::vector<std::vector<std::size_t>> intersectionIndex;
183 if constexpr (isMixedDimensional)
184 {
185 const auto numLowDimEntities = dimTarget < dimDomain ? targetTree.entitySet().size()
186 : domainTree.entitySet().size();
187 intersectionMap.resize(numLowDimEntities);
188 intersectionIndex.resize(numLowDimEntities);
189 }
190
191 // reserve memory for storing the intersections. In case of grids of
192 // different dimensionality this might be an overestimate. We get rid
193 // of the overhead memory at the end of this function though.
194 intersections_.clear();
195 intersections_.reserve(rawIntersections.size());
196
197 for (const auto& rawIntersection : rawIntersections)
198 {
199 bool add = true;
200
201 // Check if intersection was already inserted.
202 // In this case we only add new neighbor information as the geometry is identical.
203 if constexpr (isMixedDimensional)
204 {
205 const auto lowDimNeighborIdx = getLowDimNeighborIdx_(rawIntersection);
206 for (int i = 0; i < intersectionMap[lowDimNeighborIdx].size(); ++i)
207 {
208 if (rawIntersection.cornersMatch(intersectionMap[lowDimNeighborIdx][i]))
209 {
210 add = false;
211 // only add the pair of neighbors using the insertionIndex
212 auto idx = intersectionIndex[lowDimNeighborIdx][i];
213 intersections_[idx].addNeighbors(rawIntersection.first(), rawIntersection.second());
214 break;
215 }
216 }
217 }
218
219 if(add)
220 {
221 // maybe add to the map
222 if constexpr (isMixedDimensional)
223 {
224 intersectionMap[getLowDimNeighborIdx_(rawIntersection)].push_back(rawIntersection.corners());
225 intersectionIndex[getLowDimNeighborIdx_(rawIntersection)].push_back(intersections_.size());
226 }
227
228 // add new intersection and add the neighbors
229 intersections_.emplace_back(domainTree, targetTree);
230 intersections_.back().setCorners(rawIntersection.corners());
231 intersections_.back().addNeighbors(rawIntersection.first(), rawIntersection.second());
232 }
233 }
234
235 intersections_.shrink_to_fit();
236 std::cout << "Computed " << size() << " intersection entities in " << timer.elapsed() << std::endl;
237 }
238
240 typename Intersections::const_iterator ibegin() const
241 { return intersections_.begin(); }
242
244 typename Intersections::const_iterator iend() const
245 { return intersections_.end(); }
246
248 std::size_t size() const
249 { return intersections_.size(); }
250
256 friend Dune::IteratorRange<EntityIterator> intersections(const IntersectionEntitySet& set)
257 { return {set.ibegin(), set.iend()}; }
258
259private:
260 template<class RawIntersection,
261 bool enable = isMixedDimensional, std::enable_if_t<enable, int> = 0>
262 auto getLowDimNeighborIdx_(const RawIntersection& is)
263 {
264 if constexpr (dimTarget < dimDomain)
265 return is.second();
266 else
267 return is.first();
268 }
269
270 Intersections intersections_;
271
272 std::shared_ptr<const DomainTree> domainTree_;
273 std::shared_ptr<const TargetTree> targetTree_;
274};
275
276} // end namespace Dumux
277
278#endif
An axis-aligned bounding box volume hierarchy for dune grids.
An axis-aligned bounding box volume tree implementation.
Definition: boundingboxtree.hh:56
const EntitySet & entitySet() const
the entity set this tree was built with
Definition: boundingboxtree.hh:123
A class representing the intersection entities two geometric entity sets.
Definition: intersectionentityset.hh:43
IntersectionEntity Entity
make intersection entity type available
Definition: intersectionentityset.hh:136
IntersectionEntitySet()=default
Default constructor.
void build(std::shared_ptr< const DomainEntitySet > domainSet, std::shared_ptr< const TargetEntitySet > targetSet)
Build intersections.
Definition: intersectionentityset.hh:148
std::size_t size() const
the number of intersections
Definition: intersectionentityset.hh:248
void build(const DomainTree &domainTree, const TargetTree &targetTree)
Build intersections.
Definition: intersectionentityset.hh:170
void build(std::shared_ptr< const DomainTree > domainTree, std::shared_ptr< const TargetTree > targetTree)
Build intersections.
Definition: intersectionentityset.hh:158
typename Intersections::const_iterator EntityIterator
make entity iterator type available
Definition: intersectionentityset.hh:138
friend Dune::IteratorRange< EntityIterator > intersections(const IntersectionEntitySet &set)
Range generator to iterate with range-based for loops over all intersections as follows: for (const a...
Definition: intersectionentityset.hh:256
Intersections::const_iterator ibegin() const
return begin iterator to intersection container
Definition: intersectionentityset.hh:240
Intersections::const_iterator iend() const
return end iterator to intersection container
Definition: intersectionentityset.hh:244
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:102
Algorithms that finds which geometric entities intersect.
Definition: adapt.hh:17