13#ifndef DUMUX_GEOMETRY_INTERSECTION_ENTITY_SET_HH
14#define DUMUX_GEOMETRY_INTERSECTION_ENTITY_SET_HH
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>
41template<
class DomainEntitySet,
class TargetEntitySet>
47 using ctype =
typename Dune::PromotionTraits<typename DomainEntitySet::ctype, typename TargetEntitySet::ctype>::PromotedType;
49 static constexpr int dimWorld = DomainEntitySet::dimensionworld;
50 static_assert(dimWorld == int(TargetEntitySet::dimensionworld),
"Entity sets must have the same world dimension");
52 using GlobalPosition = Dune::FieldVector<ctype, dimWorld>;
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;
61 class IntersectionEntity
63 static constexpr int dimIs = std::min(dimDomain, dimTarget);
64 using Geometry = Dune::AffineGeometry<ctype, dimIs, dimWorld>;
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>>>;
70 static constexpr auto domainIdx = Dune::index_constant<0>{};
71 static constexpr auto targetIdx = Dune::index_constant<1>{};
75 : domainTree_(domainTree)
76 , targetTree_(targetTree)
80 void setCorners(
const std::vector<GlobalPosition>& corners)
83 assert(corners.size() == dimIs + 1);
87 void addNeighbors(std::size_t domain, std::size_t target)
89 if (numDomainNeighbors() == 0 && numTargetNeighbors() == 0)
91 std::get<domainIdx>(neighbors_).push_back(domain);
92 std::get<targetIdx>(neighbors_).push_back(target);
94 else if (dimDomain > dimTarget)
95 std::get<domainIdx>(neighbors_).push_back(domain);
97 else if (dimTarget > dimDomain)
98 std::get<targetIdx>(neighbors_).push_back(target);
101 DUNE_THROW(Dune::InvalidStateException,
"Cannot add more than one neighbor per side for equidimensional intersection!");
105 Geometry geometry()
const
106 {
return Geometry(Dune::GeometryTypes::simplex(dimIs), corners_); }
109 std::size_t numDomainNeighbors()
const
110 {
return std::get<domainIdx>(neighbors_).size(); }
113 std::size_t numTargetNeighbors()
const
114 {
return std::get<targetIdx>(neighbors_).size(); }
117 typename DomainEntitySet::Entity domainEntity(
unsigned int n = 0)
const
118 {
return domainTree_.entitySet().entity(std::get<domainIdx>(neighbors_)[n]); }
121 typename TargetEntitySet::Entity targetEntity(
unsigned int n = 0)
const
122 {
return targetTree_.entitySet().entity(std::get<targetIdx>(neighbors_)[n]); }
125 IndexStorage neighbors_;
126 std::vector<GlobalPosition> corners_;
132 using Intersections = std::vector<IntersectionEntity>;
148 void build(std::shared_ptr<const DomainEntitySet> domainSet, std::shared_ptr<const TargetEntitySet> targetSet)
150 domainTree_ = std::make_shared<DomainTree>(domainSet);
151 targetTree_ = std::make_shared<TargetTree>(targetSet);
152 build(*domainTree_, *targetTree_);
158 void build(std::shared_ptr<const DomainTree> domainTree, std::shared_ptr<const TargetTree> targetTree)
161 domainTree_ = domainTree;
162 targetTree_ = targetTree;
163 build(*domainTree_, *targetTree_);
181 std::vector<std::vector<std::vector<GlobalPosition>>> intersectionMap;
182 std::vector<std::vector<std::size_t>> intersectionIndex;
183 if constexpr (isMixedDimensional)
185 const auto numLowDimEntities = dimTarget < dimDomain ? targetTree.
entitySet().size()
187 intersectionMap.resize(numLowDimEntities);
188 intersectionIndex.resize(numLowDimEntities);
194 intersections_.clear();
195 intersections_.reserve(rawIntersections.size());
197 for (
const auto& rawIntersection : rawIntersections)
203 if constexpr (isMixedDimensional)
205 const auto lowDimNeighborIdx = getLowDimNeighborIdx_(rawIntersection);
206 for (
int i = 0; i < intersectionMap[lowDimNeighborIdx].size(); ++i)
208 if (rawIntersection.cornersMatch(intersectionMap[lowDimNeighborIdx][i]))
212 auto idx = intersectionIndex[lowDimNeighborIdx][i];
213 intersections_[idx].addNeighbors(rawIntersection.first(), rawIntersection.second());
222 if constexpr (isMixedDimensional)
224 intersectionMap[getLowDimNeighborIdx_(rawIntersection)].push_back(rawIntersection.corners());
225 intersectionIndex[getLowDimNeighborIdx_(rawIntersection)].push_back(intersections_.size());
229 intersections_.emplace_back(domainTree, targetTree);
230 intersections_.back().setCorners(rawIntersection.corners());
231 intersections_.back().addNeighbors(rawIntersection.first(), rawIntersection.second());
235 intersections_.shrink_to_fit();
236 std::cout <<
"Computed " <<
size() <<
" intersection entities in " << timer.elapsed() << std::endl;
240 typename Intersections::const_iterator
ibegin()
const
241 {
return intersections_.begin(); }
244 typename Intersections::const_iterator
iend()
const
245 {
return intersections_.end(); }
249 {
return intersections_.size(); }
260 template<
class RawIntersection,
261 bool enable = isMixedDimensional, std::enable_if_t<enable, int> = 0>
262 auto getLowDimNeighborIdx_(
const RawIntersection& is)
264 if constexpr (dimTarget < dimDomain)
270 Intersections intersections_;
272 std::shared_ptr<const DomainTree> domainTree_;
273 std::shared_ptr<const TargetTree> targetTree_;
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.