59 using ctype =
typename Dune::PromotionTraits<typename DomainEntitySet::ctype, typename TargetEntitySet::ctype>::PromotedType;
61 static constexpr int dimWorld = DomainEntitySet::dimensionworld;
62 static_assert(dimWorld == int(TargetEntitySet::dimensionworld),
"Entity sets must have the same world dimension");
64 using GlobalPosition = Dune::FieldVector<ctype, dimWorld>;
66 static constexpr int dimDomain = DomainEntitySet::Entity::Geometry::mydimension;
67 static constexpr int dimTarget = TargetEntitySet::Entity::Geometry::mydimension;
68 static constexpr bool isMixedDimensional = dimDomain != dimTarget;
73 class IntersectionEntity
75 static constexpr int dimIs = std::min(dimDomain, dimTarget);
76 using Geometry = Dune::AffineGeometry<ctype, dimIs, dimWorld>;
79 using IndexStorage = std::pair<std::conditional_t<dimDomain <= dimTarget, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>,
80 std::conditional_t<dimTarget <= dimDomain, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>>;
82 static constexpr auto domainIdx = Dune::index_constant<0>{};
83 static constexpr auto targetIdx = Dune::index_constant<1>{};
86 IntersectionEntity(
const DomainTree& domainTree,
const TargetTree& targetTree)
87 : domainTree_(domainTree)
88 , targetTree_(targetTree)
92 void setCorners(
const std::vector<GlobalPosition>& corners)
95 assert(corners.size() == dimIs + 1);
99 void addNeighbors(std::size_t domain, std::size_t target)
101 if (numDomainNeighbors() == 0 && numTargetNeighbors() == 0)
103 std::get<domainIdx>(neighbors_).push_back(domain);
104 std::get<targetIdx>(neighbors_).push_back(target);
106 else if (dimDomain > dimTarget)
107 std::get<domainIdx>(neighbors_).push_back(domain);
109 else if (dimTarget > dimDomain)
110 std::get<targetIdx>(neighbors_).push_back(target);
113 DUNE_THROW(Dune::InvalidStateException,
"Cannot add more than one neighbor per side for equidimensional intersection!");
117 Geometry geometry()
const
118 {
return Geometry(Dune::GeometryTypes::simplex(dimIs), corners_); }
121 std::size_t numDomainNeighbors()
const
122 {
return std::get<domainIdx>(neighbors_).size(); }
125 std::size_t numTargetNeighbors()
const
126 {
return std::get<targetIdx>(neighbors_).size(); }
129 typename DomainEntitySet::Entity domainEntity(
unsigned int n = 0)
const
130 {
return domainTree_.entitySet().entity(std::get<domainIdx>(neighbors_)[n]); }
133 typename TargetEntitySet::Entity targetEntity(
unsigned int n = 0)
const
134 {
return targetTree_.entitySet().entity(std::get<targetIdx>(neighbors_)[n]); }
137 IndexStorage neighbors_;
138 std::vector<GlobalPosition> corners_;
140 const DomainTree& domainTree_;
141 const TargetTree& targetTree_;
144 using Intersections = std::vector<IntersectionEntity>;
160 void build(std::shared_ptr<const DomainEntitySet> domainSet, std::shared_ptr<const TargetEntitySet> targetSet)
162 domainTree_ = std::make_shared<DomainTree>(domainSet);
163 targetTree_ = std::make_shared<TargetTree>(targetSet);
164 build(*domainTree_, *targetTree_);
170 void build(std::shared_ptr<const DomainTree> domainTree, std::shared_ptr<const TargetTree> targetTree)
173 domainTree_ = domainTree;
174 targetTree_ = targetTree_;
175 build(*domainTree_, *targetTree_);
182 void build(
const DomainTree& domainTree,
const TargetTree& targetTree)
193 std::vector<std::vector<std::vector<GlobalPosition>>> intersectionMap;
194 std::vector<std::vector<std::size_t>> intersectionIndex;
195 if constexpr (isMixedDimensional)
197 const auto numLowDimEntities = dimTarget < dimDomain ? targetTree.
entitySet().size()
199 intersectionMap.resize(numLowDimEntities);
200 intersectionIndex.resize(numLowDimEntities);
206 intersections_.clear();
207 intersections_.reserve(rawIntersections.size());
209 for (
const auto& rawIntersection : rawIntersections)
215 if constexpr (isMixedDimensional)
217 const auto lowDimNeighborIdx = getLowDimNeighborIdx_(rawIntersection);
218 for (
int i = 0; i < intersectionMap[lowDimNeighborIdx].size(); ++i)
220 if (rawIntersection.cornersMatch(intersectionMap[lowDimNeighborIdx][i]))
224 auto idx = intersectionIndex[lowDimNeighborIdx][i];
225 intersections_[idx].addNeighbors(rawIntersection.first(), rawIntersection.second());
234 if constexpr (isMixedDimensional)
236 intersectionMap[getLowDimNeighborIdx_(rawIntersection)].push_back(rawIntersection.corners());
237 intersectionIndex[getLowDimNeighborIdx_(rawIntersection)].push_back(intersections_.size());
241 intersections_.emplace_back(domainTree, targetTree);
242 intersections_.back().setCorners(rawIntersection.corners());
243 intersections_.back().addNeighbors(rawIntersection.first(), rawIntersection.second());
247 intersections_.shrink_to_fit();
248 std::cout <<
"Computed " <<
size() <<
" intersection entities in " << timer.elapsed() << std::endl;
252 typename Intersections::const_iterator
ibegin()
const
253 {
return intersections_.begin(); }
256 typename Intersections::const_iterator
iend()
const
257 {
return intersections_.end(); }
261 {
return intersections_.size(); }
272 template<
class RawIntersection,
273 bool enable = isMixedDimensional, std::enable_if_t<enable, int> = 0>
274 auto getLowDimNeighborIdx_(
const RawIntersection& is)
276 if constexpr (dimTarget < dimDomain)
282 Intersections intersections_;
284 std::shared_ptr<const DomainTree> domainTree_;
285 std::shared_ptr<const TargetTree> targetTree_;
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:99