Line | Branch | Exec | Source |
---|---|---|---|
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-FileCopyrightText: Copyright © DuMux Project contributors, see AUTHORS.md in root folder | ||
5 | // SPDX-License-Identifier: GPL-3.0-or-later | ||
6 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup Geometry | ||
10 | * \brief A class representing the intersection entities of two geometric entity sets | ||
11 | */ | ||
12 | |||
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 | |||
32 | #include <dumux/geometry/boundingboxtree.hh> | ||
33 | #include <dumux/geometry/intersectingentities.hh> | ||
34 | |||
35 | namespace Dumux::Detail::Intersections { | ||
36 | |||
37 | /*! | ||
38 | * \brief A class representing an intersection entity | ||
39 | */ | ||
40 | template<class DomainEntitySet, class TargetEntitySet> | ||
41 |
5/8✗ Branch 0 not taken.
✗ Branch 1 not taken.
✓ Branch 2 taken 681360 times.
✓ Branch 3 taken 168 times.
✓ Branch 4 taken 2 times.
✓ Branch 5 taken 321312 times.
✓ Branch 6 taken 310 times.
✗ Branch 7 not taken.
|
1685383 | class IntersectionEntity |
42 | { | ||
43 | static constexpr int dimDomain = DomainEntitySet::Entity::Geometry::mydimension; | ||
44 | static constexpr int dimTarget = TargetEntitySet::Entity::Geometry::mydimension; | ||
45 | using ctype = typename Dune::PromotionTraits<typename DomainEntitySet::ctype, typename TargetEntitySet::ctype>::PromotedType; | ||
46 | |||
47 | static constexpr int dimWorld = DomainEntitySet::dimensionworld; | ||
48 | static_assert(dimWorld == int(TargetEntitySet::dimensionworld), "Entity sets must have the same world dimension"); | ||
49 | |||
50 | static constexpr int dimIs = std::min(dimDomain, dimTarget); | ||
51 | |||
52 | using GlobalPosition = Dune::FieldVector<ctype, dimWorld>; | ||
53 | using Geometry = Dune::AffineGeometry<ctype, dimIs, dimWorld>; // geometries are always simplices | ||
54 | |||
55 | // we can only have multiple neighbors in the mixeddimensional case and then only for the side with the largest dimension | ||
56 | using IndexStorage = std::pair<std::conditional_t<dimDomain <= dimTarget, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>, | ||
57 | std::conditional_t<dimTarget <= dimDomain, Dune::ReservedVector<std::size_t, 1>, std::vector<std::size_t>>>; | ||
58 | |||
59 | static constexpr auto domainIdx = Dune::index_constant<0>{}; | ||
60 | static constexpr auto targetIdx = Dune::index_constant<1>{}; | ||
61 | |||
62 | public: | ||
63 | 397553 | IntersectionEntity(const DomainEntitySet& domainEntitySet, const TargetEntitySet& targetEntitySet) | |
64 | 395267 | : domainEntitySet_(domainEntitySet) | |
65 | 395267 | , targetEntitySet_(targetEntitySet) | |
66 | {} | ||
67 | |||
68 | //! set the intersection geometry corners | ||
69 | 402336 | void setCorners(const std::vector<GlobalPosition>& corners) | |
70 | { | ||
71 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 397553 times.
|
402336 | corners_ = corners; |
72 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 397553 times.
|
402336 | assert(corners.size() == dimIs + 1); // Only simplex intersections are allowed |
73 | 402336 | } | |
74 | |||
75 | //! add a pair of neighbor elements | ||
76 | 403333 | void addNeighbors(std::size_t domain, std::size_t target) | |
77 | { | ||
78 |
3/4✓ Branch 0 taken 397553 times.
✓ Branch 1 taken 997 times.
✓ Branch 2 taken 397553 times.
✗ Branch 3 not taken.
|
403333 | if (numDomainNeighbors() == 0 && numTargetNeighbors() == 0) |
79 | { | ||
80 | 402336 | std::get<domainIdx>(neighbors_).push_back(domain); | |
81 | 402336 | std::get<targetIdx>(neighbors_).push_back(target); | |
82 | } | ||
83 | else if (dimDomain > dimTarget) | ||
84 | 997 | std::get<domainIdx>(neighbors_).push_back(domain); | |
85 | |||
86 | else if (dimTarget > dimDomain) | ||
87 | ✗ | std::get<targetIdx>(neighbors_).push_back(target); | |
88 | |||
89 | else | ||
90 |
0/20✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✗ Branch 16 not taken.
✗ Branch 17 not taken.
✗ Branch 19 not taken.
✗ Branch 20 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 25 not taken.
✗ Branch 26 not taken.
✗ Branch 28 not taken.
✗ Branch 29 not taken.
|
42661 | DUNE_THROW(Dune::InvalidStateException, "Cannot add more than one neighbor per side for equidimensional intersection!"); |
91 | 403333 | } | |
92 | |||
93 | //! get the intersection geometry | ||
94 | 414626 | Geometry geometry() const | |
95 |
4/7✓ Branch 1 taken 395273 times.
✗ Branch 2 not taken.
✓ Branch 5 taken 19353 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 8 taken 2873 times.
✓ Branch 4 taken 45264 times.
|
829252 | { return Geometry(Dune::GeometryTypes::simplex(dimIs), corners_); } |
96 | |||
97 | //! get the number of domain neighbors of this intersection | ||
98 | 1833702 | std::size_t numDomainNeighbors() const | |
99 |
12/12✓ Branch 0 taken 394957 times.
✓ Branch 1 taken 997 times.
✓ Branch 2 taken 324547 times.
✓ Branch 3 taken 347403 times.
✓ Branch 4 taken 9559 times.
✓ Branch 5 taken 9241 times.
✓ Branch 6 taken 25786 times.
✓ Branch 7 taken 318 times.
✓ Branch 9 taken 25381 times.
✓ Branch 10 taken 3118 times.
✓ Branch 8 taken 25786 times.
✓ Branch 11 taken 2773 times.
|
1810954 | { return std::get<domainIdx>(neighbors_).size(); } |
100 | |||
101 | //! get the number of target neighbors of this intersection | ||
102 | 401807 | std::size_t numTargetNeighbors() const | |
103 |
5/10✓ Branch 0 taken 394957 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 2278 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 2286 times.
✓ Branch 5 taken 1968 times.
✗ Branch 6 not taken.
✓ Branch 7 taken 318 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
|
401807 | { return std::get<targetIdx>(neighbors_).size(); } |
104 | |||
105 | //! get the nth domain neighbor entity | ||
106 | 1069876 | typename DomainEntitySet::Entity domainEntity(unsigned int n = 0) const | |
107 |
4/8✓ Branch 1 taken 395822 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 638774 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 7056 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 28224 times.
✗ Branch 11 not taken.
|
1069876 | { return domainEntitySet_.entity(std::get<domainIdx>(neighbors_)[n]); } |
108 | |||
109 | //! get the nth target neighbor entity | ||
110 | 692136 | typename TargetEntitySet::Entity targetEntity(unsigned int n = 0) const | |
111 |
5/9✓ Branch 1 taken 314256 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 314256 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 7056 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 7056 times.
✗ Branch 11 not taken.
✓ Branch 3 taken 1968 times.
|
692136 | { return targetEntitySet_.entity(std::get<targetIdx>(neighbors_)[n]); } |
112 | |||
113 | private: | ||
114 | IndexStorage neighbors_; | ||
115 | std::vector<GlobalPosition> corners_; | ||
116 | |||
117 | const DomainEntitySet& domainEntitySet_; | ||
118 | const TargetEntitySet& targetEntitySet_; | ||
119 | }; | ||
120 | |||
121 | } // end namespace Dumux::Detail::Intersections | ||
122 | |||
123 | namespace Dumux { | ||
124 | |||
125 | /*! | ||
126 | * \ingroup Geometry | ||
127 | * \brief A class representing the intersection entities two geometric entity sets | ||
128 | */ | ||
129 | template<class DomainEntitySet, class TargetEntitySet> | ||
130 | class IntersectionEntitySet | ||
131 | { | ||
132 | using DomainTree = BoundingBoxTree<DomainEntitySet>; | ||
133 | using TargetTree = BoundingBoxTree<TargetEntitySet>; | ||
134 | |||
135 | using ctype = typename Dune::PromotionTraits<typename DomainEntitySet::ctype, typename TargetEntitySet::ctype>::PromotedType; | ||
136 | |||
137 | static constexpr int dimWorld = DomainEntitySet::dimensionworld; | ||
138 | static_assert(dimWorld == int(TargetEntitySet::dimensionworld), "Entity sets must have the same world dimension"); | ||
139 | |||
140 | using GlobalPosition = Dune::FieldVector<ctype, dimWorld>; | ||
141 | |||
142 | static constexpr int dimDomain = DomainEntitySet::Entity::Geometry::mydimension; | ||
143 | static constexpr int dimTarget = TargetEntitySet::Entity::Geometry::mydimension; | ||
144 | static constexpr bool isMixedDimensional = dimDomain != dimTarget; | ||
145 | |||
146 | using IntersectionEntity = Dumux::Detail::Intersections::IntersectionEntity<DomainEntitySet, TargetEntitySet>; | ||
147 | using Intersections = std::vector<IntersectionEntity>; | ||
148 | |||
149 | using DomainGeometry = typename DomainEntitySet::Entity::Geometry; | ||
150 | using TargetGeometry = typename TargetEntitySet::Entity::Geometry; | ||
151 | |||
152 | public: | ||
153 | //! make intersection entity type available | ||
154 | using Entity = IntersectionEntity; | ||
155 | //! make entity iterator type available | ||
156 | using EntityIterator = typename Intersections::const_iterator; | ||
157 | |||
158 | /*! | ||
159 | * \brief Default constructor | ||
160 | */ | ||
161 | 28 | IntersectionEntitySet() = default; | |
162 | |||
163 | /*! | ||
164 | * \brief Build intersections | ||
165 | */ | ||
166 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 37 times.
|
114 | void build(std::shared_ptr<const DomainEntitySet> domainSet, std::shared_ptr<const TargetEntitySet> targetSet) |
167 | { | ||
168 | // specialized algorithm if the target geometry is only a single entity | ||
169 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 37 times.
|
114 | if (targetSet->size() == 1) |
170 | { | ||
171 |
2/2✓ Branch 0 taken 35 times.
✓ Branch 1 taken 4 times.
|
75 | domainTree_ = std::make_shared<DomainTree>(domainSet); |
172 | 40 | targetEntitySet_ = targetSet; | |
173 |
0/2✗ Branch 0 not taken.
✗ Branch 1 not taken.
|
40 | build_(*domainTree_, targetEntitySet_->entity(0).geometry()); |
174 | } | ||
175 | // specialized algorithm if the domain geometry is only a single entity | ||
176 |
2/3✗ Branch 0 not taken.
✓ Branch 1 taken 36 times.
✓ Branch 2 taken 1 times.
|
74 | else if (domainSet->size() == 1) |
177 | { | ||
178 | ✗ | domainEntitySet_ = domainSet; | |
179 | ✗ | targetTree_ = std::make_shared<TargetTree>(targetSet); | |
180 | ✗ | build_(domainEntitySet_->entity(0).geometry(), *targetTree_); | |
181 | } | ||
182 | else | ||
183 | { | ||
184 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 37 times.
|
74 | domainTree_ = std::make_shared<DomainTree>(domainSet); |
185 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 37 times.
|
74 | targetTree_ = std::make_shared<TargetTree>(targetSet); |
186 | 74 | build(*domainTree_, *targetTree_); | |
187 | } | ||
188 | 114 | } | |
189 | |||
190 | /*! | ||
191 | * \brief Build intersections | ||
192 | */ | ||
193 | void build(std::shared_ptr<const DomainTree> domainTree, std::shared_ptr<const TargetTree> targetTree) | ||
194 | { | ||
195 | // make sure the tree don't get out of scope | ||
196 | domainTree_ = domainTree; | ||
197 | targetTree_ = targetTree; | ||
198 | build(*domainTree_, *targetTree_); | ||
199 | } | ||
200 | |||
201 | /*! | ||
202 | * \brief Build intersections | ||
203 | * \note If you call this, make sure the bounding box tree stays alive for the life-time of this object | ||
204 | */ | ||
205 | 118 | void build(const DomainTree& domainTree, const TargetTree& targetTree) | |
206 | { | ||
207 | 118 | Dune::Timer timer; | |
208 | |||
209 |
1/2✓ Branch 1 taken 79 times.
✗ Branch 2 not taken.
|
118 | const auto rawIntersections = intersectingEntities(domainTree, targetTree); |
210 |
1/2✓ Branch 1 taken 79 times.
✗ Branch 2 not taken.
|
118 | build_(domainTree.entitySet(), targetTree.entitySet(), rawIntersections); |
211 | |||
212 |
3/6✓ Branch 1 taken 79 times.
✗ Branch 2 not taken.
✓ Branch 5 taken 79 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 79 times.
✗ Branch 9 not taken.
|
118 | std::cout << "Computed " << size() << " intersection entities in " << timer.elapsed() << std::endl; |
213 | 118 | } | |
214 | |||
215 | //! return begin iterator to intersection container | ||
216 | 122 | typename Intersections::const_iterator ibegin() const | |
217 | 122 | { return intersections_.begin(); } | |
218 | |||
219 | //! return end iterator to intersection container | ||
220 | 122 | typename Intersections::const_iterator iend() const | |
221 | 122 | { return intersections_.end(); } | |
222 | |||
223 | //! the number of intersections | ||
224 | 181 | std::size_t size() const | |
225 |
5/14✓ Branch 1 taken 77 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 2 times.
✗ Branch 5 not taken.
✓ Branch 7 taken 1 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
✓ Branch 10 taken 1 times.
✗ Branch 12 not taken.
✗ Branch 13 not taken.
✗ Branch 14 not taken.
✓ Branch 15 taken 36 times.
✗ Branch 17 not taken.
✗ Branch 18 not taken.
|
165 | { return intersections_.size(); } |
226 | |||
227 | /*! | ||
228 | * \brief Range generator to iterate with range-based for loops over all intersections | ||
229 | * as follows: for (const auto& is : intersections(glue)) { ... } | ||
230 | * \note free function | ||
231 | */ | ||
232 | 122 | friend Dune::IteratorRange<EntityIterator> intersections(const IntersectionEntitySet& set) | |
233 | 122 | { return {set.ibegin(), set.iend()}; } | |
234 | |||
235 | private: | ||
236 | 40 | void build_(const DomainTree& domainTree, const TargetGeometry& targetGeometry) | |
237 | { | ||
238 | 40 | Dune::Timer timer; | |
239 | |||
240 | 40 | const auto rawIntersections = intersectingEntities(targetGeometry, domainTree); | |
241 | // because in intersectingEntities geometry always comes first (domain set) | ||
242 | // but here we want the geometry to be the second (target set), we need to flip the index sets | ||
243 | // of the resulting intersections | ||
244 | 40 | const bool flipNeighborIndices = true; | |
245 |
1/2✓ Branch 1 taken 39 times.
✗ Branch 2 not taken.
|
40 | build_(domainTree.entitySet(), *targetEntitySet_, rawIntersections, flipNeighborIndices); |
246 | |||
247 |
4/8✓ Branch 1 taken 39 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 39 times.
✗ Branch 5 not taken.
✓ Branch 8 taken 39 times.
✗ Branch 9 not taken.
✓ Branch 11 taken 39 times.
✗ Branch 12 not taken.
|
40 | std::cout << "Computed " << size() << " intersection entities in " << timer.elapsed() << std::endl; |
248 | 40 | } | |
249 | |||
250 | ✗ | void build_(const DomainGeometry& domainGeometry, const TargetTree& targetTree) | |
251 | { | ||
252 | ✗ | Dune::Timer timer; | |
253 | |||
254 | ✗ | const auto rawIntersections = intersectingEntities(domainGeometry, targetTree); | |
255 | ✗ | build_(*domainEntitySet_, targetTree.entitySet(), rawIntersections); | |
256 | |||
257 | ✗ | std::cout << "Computed " << size() << " intersection entities in " << timer.elapsed() << std::endl; | |
258 | ✗ | } | |
259 | |||
260 | template<class RawIntersections> | ||
261 | 158 | void build_(const DomainEntitySet& domainEntitySet, const TargetEntitySet& targetEntitySet, | |
262 | const RawIntersections& rawIntersections, | ||
263 | bool flipNeighborIndices = false) | ||
264 | { | ||
265 | // create a map to check whether intersection geometries were already inserted | ||
266 | // Note that this can only occur if the grids have different dimensionality. | ||
267 | // If this is the case, we keep track of the intersections using the indices of the lower- | ||
268 | // dimensional entities which is identical for all geometrically identical intersections. | ||
269 | 158 | std::vector<std::vector<std::vector<GlobalPosition>>> intersectionMap; | |
270 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
158 | std::vector<std::vector<std::size_t>> intersectionIndex; |
271 | if constexpr (isMixedDimensional) | ||
272 | { | ||
273 | 33 | const auto numLowDimEntities = dimTarget < dimDomain ? targetEntitySet.size() : domainEntitySet.size(); | |
274 |
1/2✓ Branch 1 taken 30 times.
✗ Branch 2 not taken.
|
33 | intersectionMap.resize(numLowDimEntities); |
275 |
1/2✓ Branch 1 taken 30 times.
✗ Branch 2 not taken.
|
33 | intersectionIndex.resize(numLowDimEntities); |
276 | } | ||
277 | |||
278 | // reserve memory for storing the intersections. In case of grids of | ||
279 | // different dimensionality this might be an overestimate. We get rid | ||
280 | // of the overhead memory at the end of this function though. | ||
281 |
2/2✓ Branch 0 taken 31 times.
✓ Branch 1 taken 87 times.
|
158 | intersections_.clear(); |
282 |
1/2✓ Branch 1 taken 118 times.
✗ Branch 2 not taken.
|
158 | intersections_.reserve(rawIntersections.size()); |
283 | |||
284 |
2/2✓ Branch 0 taken 398550 times.
✓ Branch 1 taken 118 times.
|
403491 | for (const auto& rawIntersection : rawIntersections) |
285 | { | ||
286 | 403333 | bool add = true; | |
287 | |||
288 | // Check if intersection was already inserted. | ||
289 | // In this case we only add new neighbor information as the geometry is identical. | ||
290 | if constexpr (isMixedDimensional) | ||
291 | { | ||
292 | 85322 | const auto lowDimNeighborIdx = getLowDimNeighborIndex_(rawIntersection, flipNeighborIndices); | |
293 |
2/2✓ Branch 0 taken 490493 times.
✓ Branch 1 taken 37193 times.
|
532869 | for (int i = 0; i < intersectionMap[lowDimNeighborIdx].size(); ++i) |
294 | { | ||
295 |
2/2✓ Branch 0 taken 997 times.
✓ Branch 1 taken 489496 times.
|
491205 | if (rawIntersection.cornersMatch(intersectionMap[lowDimNeighborIdx][i])) |
296 | { | ||
297 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 997 times.
|
997 | add = false; |
298 | // only add the pair of neighbors using the insertionIndex | ||
299 | 997 | auto idx = intersectionIndex[lowDimNeighborIdx][i]; | |
300 |
1/2✓ Branch 1 taken 997 times.
✗ Branch 2 not taken.
|
997 | intersections_[idx].addNeighbors(domainNeighborIndex_(rawIntersection, flipNeighborIndices), targetNeighborIndex_(rawIntersection, flipNeighborIndices)); |
301 | break; | ||
302 | } | ||
303 | } | ||
304 | } | ||
305 | |||
306 | if(add) | ||
307 | { | ||
308 | // maybe add to the map | ||
309 | if constexpr (isMixedDimensional) | ||
310 | { | ||
311 | 41664 | const auto lowDimNeighborIdx = getLowDimNeighborIndex_(rawIntersection, flipNeighborIndices); | |
312 | 41664 | intersectionMap[lowDimNeighborIdx].push_back(rawIntersection.corners()); | |
313 |
1/2✓ Branch 1 taken 37193 times.
✗ Branch 2 not taken.
|
41664 | intersectionIndex[lowDimNeighborIdx].push_back(intersections_.size()); |
314 | } | ||
315 | |||
316 | // add new intersection and add the neighbors | ||
317 |
1/2✓ Branch 1 taken 397553 times.
✗ Branch 2 not taken.
|
402336 | intersections_.emplace_back(domainEntitySet, targetEntitySet); |
318 |
1/2✓ Branch 1 taken 397553 times.
✗ Branch 2 not taken.
|
402336 | intersections_.back().setCorners(rawIntersection.corners()); |
319 |
3/4✓ Branch 0 taken 45794 times.
✓ Branch 1 taken 351759 times.
✓ Branch 3 taken 397553 times.
✗ Branch 4 not taken.
|
804672 | intersections_.back().addNeighbors(domainNeighborIndex_(rawIntersection, flipNeighborIndices), targetNeighborIndex_(rawIntersection, flipNeighborIndices)); |
320 | } | ||
321 | } | ||
322 | |||
323 |
2/2✓ Branch 0 taken 32 times.
✓ Branch 1 taken 86 times.
|
158 | intersections_.shrink_to_fit(); |
324 | 158 | } | |
325 | |||
326 | template<class RawIntersection, | ||
327 | bool enable = isMixedDimensional, std::enable_if_t<enable, int> = 0> | ||
328 | 38190 | auto getLowDimNeighborIndex_(const RawIntersection& is, bool flipNeighborIndices) | |
329 | { | ||
330 | if constexpr (dimTarget < dimDomain) | ||
331 | 36222 | return targetNeighborIndex_(is, flipNeighborIndices); | |
332 | else | ||
333 | 1968 | return domainNeighborIndex_(is, flipNeighborIndices); | |
334 | } | ||
335 | |||
336 | template<class RawIntersection> | ||
337 | 400518 | auto domainNeighborIndex_(const RawIntersection& is, bool flipNeighborIndices) | |
338 |
6/11✓ Branch 0 taken 45794 times.
✓ Branch 1 taken 349163 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 2278 times.
✓ Branch 5 taken 318 times.
✗ Branch 6 not taken.
✓ Branch 8 taken 1968 times.
✗ Branch 9 not taken.
✗ Branch 10 not taken.
✓ Branch 11 taken 1968 times.
✗ Branch 4 not taken.
|
402486 | { return flipNeighborIndices ? is.second() : is.first(); } |
339 | |||
340 | template<class RawIntersection> | ||
341 | 433775 | auto targetNeighborIndex_(const RawIntersection& is, bool flipNeighborIndices) | |
342 |
7/13✓ Branch 0 taken 45794 times.
✓ Branch 1 taken 350160 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 1307 times.
✗ Branch 4 not taken.
✓ Branch 5 taken 35225 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
✓ Branch 9 taken 318 times.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✓ Branch 12 taken 2286 times.
✓ Branch 8 taken 34907 times.
|
469997 | { return flipNeighborIndices ? is.first() : is.second(); } |
343 | |||
344 | Intersections intersections_; | ||
345 | |||
346 | std::shared_ptr<const DomainTree> domainTree_; | ||
347 | std::shared_ptr<const TargetTree> targetTree_; | ||
348 | std::shared_ptr<const DomainEntitySet> domainEntitySet_; | ||
349 | std::shared_ptr<const TargetEntitySet> targetEntitySet_; | ||
350 | }; | ||
351 | |||
352 | } // end namespace Dumux | ||
353 | |||
354 | #endif | ||
355 |