CDT  1.4.1
C++ library for constrained Delaunay triangulation
Loading...
Searching...
No Matches
Triangulation.hpp
Go to the documentation of this file.
1/* This Source Code Form is subject to the terms of the Mozilla Public
2 * License, v. 2.0. If a copy of the MPL was not distributed with this
3 * file, You can obtain one at https://mozilla.org/MPL/2.0/. */
4
10#include "Triangulation.h"
11#include "portable_nth_element.hpp"
12
13#include <algorithm>
14#include <cassert>
15#include <cmath>
16#include <deque>
17#include <stdexcept>
18
19namespace CDT
20{
21
22typedef std::deque<TriInd> TriDeque;
23
24namespace detail
25{
26
28template <typename T>
29array<T, 3> arr3(const T& v0, const T& v1, const T& v2)
30{
31 const array<T, 3> out = {v0, v1, v2};
32 return out;
33}
34
35namespace defaults
36{
37
38const std::size_t nTargetVerts = 0;
40const VertexInsertionOrder::Enum vertexInsertionOrder =
42const IntersectingConstraintEdges::Enum intersectingEdgesStrategy =
44const float minDistToConstraintEdge(0);
45
46} // namespace defaults
47
48} // namespace detail
49
50template <typename T, typename TNearPointLocator>
52 : m_nTargetVerts(detail::defaults::nTargetVerts)
53 , m_superGeomType(detail::defaults::superGeomType)
54 , m_vertexInsertionOrder(detail::defaults::vertexInsertionOrder)
55 , m_intersectingEdgesStrategy(detail::defaults::intersectingEdgesStrategy)
56 , m_minDistToConstraintEdge(detail::defaults::minDistToConstraintEdge)
57{}
58
59template <typename T, typename TNearPointLocator>
61 const VertexInsertionOrder::Enum vertexInsertionOrder)
62 : m_nTargetVerts(detail::defaults::nTargetVerts)
63 , m_superGeomType(detail::defaults::superGeomType)
64 , m_vertexInsertionOrder(vertexInsertionOrder)
65 , m_intersectingEdgesStrategy(detail::defaults::intersectingEdgesStrategy)
66 , m_minDistToConstraintEdge(detail::defaults::minDistToConstraintEdge)
67{}
68
69template <typename T, typename TNearPointLocator>
71 const VertexInsertionOrder::Enum vertexInsertionOrder,
72 const IntersectingConstraintEdges::Enum intersectingEdgesStrategy,
73 const T minDistToConstraintEdge)
74 : m_nTargetVerts(detail::defaults::nTargetVerts)
75 , m_superGeomType(detail::defaults::superGeomType)
76 , m_vertexInsertionOrder(vertexInsertionOrder)
77 , m_intersectingEdgesStrategy(intersectingEdgesStrategy)
78 , m_minDistToConstraintEdge(minDistToConstraintEdge)
79{}
80
81template <typename T, typename TNearPointLocator>
83 const VertexInsertionOrder::Enum vertexInsertionOrder,
84 const TNearPointLocator& nearPtLocator,
85 const IntersectingConstraintEdges::Enum intersectingEdgesStrategy,
86 const T minDistToConstraintEdge)
87 : m_nearPtLocator(nearPtLocator)
88 , m_nTargetVerts(detail::defaults::nTargetVerts)
89 , m_superGeomType(detail::defaults::superGeomType)
90 , m_vertexInsertionOrder(vertexInsertionOrder)
91 , m_intersectingEdgesStrategy(intersectingEdgesStrategy)
92 , m_minDistToConstraintEdge(minDistToConstraintEdge)
93{}
94
95template <typename T, typename TNearPointLocator>
97{
98 if(m_dummyTris.empty())
99 return;
100 const TriIndUSet dummySet(m_dummyTris.begin(), m_dummyTris.end());
101 TriIndUMap triIndMap;
102 triIndMap[noNeighbor] = noNeighbor;
103 for(TriInd iT(0), iTnew(0); iT < TriInd(triangles.size()); ++iT)
104 {
105 if(dummySet.count(iT))
106 continue;
107 triIndMap[iT] = iTnew;
108 triangles[iTnew] = triangles[iT];
109 iTnew++;
110 }
111 triangles.erase(triangles.end() - dummySet.size(), triangles.end());
112
113 // remap adjacent triangle indices for vertices
114 for(TriIndVec::iterator iT = m_vertTris.begin(); iT != m_vertTris.end();
115 ++iT)
116 {
117 *iT = triIndMap[*iT];
118 }
119 // remap neighbor indices for triangles
120 for(TriangleVec::iterator t = triangles.begin(); t != triangles.end(); ++t)
121 {
122 NeighborsArr3& nn = t->neighbors;
123 for(NeighborsArr3::iterator iN = nn.begin(); iN != nn.end(); ++iN)
124 *iN = triIndMap[*iN];
125 }
126 // clear dummy triangles
127 m_dummyTris = std::vector<TriInd>();
128}
129
130template <typename T, typename TNearPointLocator>
132{
133 if(m_superGeomType != SuperGeometryType::SuperTriangle)
134 return;
135 // find triangles adjacent to super-triangle's vertices
136 TriIndUSet toErase;
137 for(TriInd iT(0); iT < TriInd(triangles.size()); ++iT)
138 {
139 if(touchesSuperTriangle(triangles[iT]))
140 toErase.insert(iT);
141 }
142 finalizeTriangulation(toErase);
143}
144
145template <typename T, typename TNearPointLocator>
147{
148 // make dummy triangles adjacent to super-triangle's vertices
149 assert(m_vertTris[0] != noNeighbor);
150 const std::stack<TriInd> seed(std::deque<TriInd>(1, m_vertTris[0]));
151 const TriIndUSet toErase = growToBoundary(seed);
152 finalizeTriangulation(toErase);
153}
154
155template <typename T, typename TNearPointLocator>
157{
158 const std::vector<LayerDepth> triDepths = calculateTriangleDepths();
159 TriIndUSet toErase;
160 toErase.reserve(triangles.size());
161 for(std::size_t iT = 0; iT != triangles.size(); ++iT)
162 {
163 if(triDepths[iT] % 2 == 0)
164 toErase.insert(static_cast<TriInd>(iT));
165 }
166 finalizeTriangulation(toErase);
167}
168
171{
172 return Edge(VertInd(e.v1() - 3), VertInd(e.v2() - 3));
173}
174
175template <typename T, typename TNearPointLocator>
177 const TriIndUSet& removedTriangles)
178{
179 if(removedTriangles.empty())
180 return;
181 // remove triangles and calculate triangle index mapping
182 TriIndUMap triIndMap;
183 for(TriInd iT(0), iTnew(0); iT < TriInd(triangles.size()); ++iT)
184 {
185 if(removedTriangles.count(iT))
186 continue;
187 triIndMap[iT] = iTnew;
188 triangles[iTnew] = triangles[iT];
189 iTnew++;
190 }
191 triangles.erase(triangles.end() - removedTriangles.size(), triangles.end());
192 // adjust triangles' neighbors
193 for(TriInd iT(0); iT < triangles.size(); ++iT)
194 {
195 Triangle& t = triangles[iT];
196 // update neighbors to account for removed triangles
197 NeighborsArr3& nn = t.neighbors;
198 for(NeighborsArr3::iterator n = nn.begin(); n != nn.end(); ++n)
199 {
200 if(removedTriangles.count(*n))
201 {
202 *n = noNeighbor;
203 }
204 else if(*n != noNeighbor)
205 {
206 *n = triIndMap[*n];
207 }
208 }
209 }
210}
211
212template <typename T, typename TNearPointLocator>
217
218template <typename T, typename TNearPointLocator>
220{
221 return m_vertTris;
222}
223
224template <typename T, typename TNearPointLocator>
226 const TriIndUSet& removedTriangles)
227{
228 eraseDummies();
229 m_vertTris = TriIndVec();
230 // remove super-triangle
231 if(m_superGeomType == SuperGeometryType::SuperTriangle)
232 {
233 vertices.erase(vertices.begin(), vertices.begin() + 3);
234 // Edge re-mapping
235 { // fixed edges
236 EdgeUSet updatedFixedEdges;
237 typedef CDT::EdgeUSet::const_iterator It;
238 for(It e = fixedEdges.begin(); e != fixedEdges.end(); ++e)
239 {
240 updatedFixedEdges.insert(RemapNoSuperTriangle(*e));
241 }
242 fixedEdges = updatedFixedEdges;
243 }
244 { // overlap count
245 unordered_map<Edge, BoundaryOverlapCount> updatedOverlapCount;
246 typedef unordered_map<Edge, BoundaryOverlapCount>::const_iterator
247 It;
248 for(It it = overlapCount.begin(); it != overlapCount.end(); ++it)
249 {
250 updatedOverlapCount.insert(std::make_pair(
251 RemapNoSuperTriangle(it->first), it->second));
252 }
253 overlapCount = updatedOverlapCount;
254 }
255 { // split edges mapping
256 unordered_map<Edge, EdgeVec> updatedPieceToOriginals;
257 typedef unordered_map<Edge, EdgeVec>::const_iterator It;
258 for(It it = pieceToOriginals.begin(); it != pieceToOriginals.end();
259 ++it)
260 {
261 EdgeVec ee = it->second;
262 for(EdgeVec::iterator eeIt = ee.begin(); eeIt != ee.end();
263 ++eeIt)
264 {
265 *eeIt = RemapNoSuperTriangle(*eeIt);
266 }
267 updatedPieceToOriginals.insert(
268 std::make_pair(RemapNoSuperTriangle(it->first), ee));
269 }
270 pieceToOriginals = updatedPieceToOriginals;
271 }
272 }
273 // remove other triangles
274 removeTriangles(removedTriangles);
275 // adjust triangle vertices: account for removed super-triangle
276 if(m_superGeomType == SuperGeometryType::SuperTriangle)
277 {
278 for(TriangleVec::iterator t = triangles.begin(); t != triangles.end();
279 ++t)
280 {
281 VerticesArr3& vv = t->vertices;
282 for(VerticesArr3::iterator v = vv.begin(); v != vv.end(); ++v)
283 {
284 *v -= 3;
285 }
286 }
287 }
288}
289
290template <typename T, typename TNearPointLocator>
292{
293 m_nearPtLocator.initialize(vertices);
294 m_nTargetVerts = static_cast<IndexSizeType>(vertices.size());
295 m_superGeomType = SuperGeometryType::Custom;
296}
297
298template <typename T, typename TNearPointLocator>
300 std::stack<TriInd> seeds) const
301{
302 TriIndUSet traversed;
303 while(!seeds.empty())
304 {
305 const TriInd iT = seeds.top();
306 seeds.pop();
307 traversed.insert(iT);
308 const Triangle& t = triangles[iT];
309 for(Index i(0); i < Index(3); ++i)
310 {
311 const Edge opEdge(t.vertices[ccw(i)], t.vertices[cw(i)]);
312 if(fixedEdges.count(opEdge))
313 continue;
314 const TriInd iN = t.neighbors[opoNbr(i)];
315 if(iN != noNeighbor && traversed.count(iN) == 0)
316 seeds.push(iN);
317 }
318 }
319 return traversed;
320}
321
322template <typename T, typename TNearPointLocator>
323void Triangulation<T, TNearPointLocator>::makeDummy(const TriInd iT)
324{
325 m_dummyTris.push_back(iT);
326}
327
328template <typename T, typename TNearPointLocator>
329TriInd Triangulation<T, TNearPointLocator>::addTriangle(const Triangle& t)
330{
331 if(m_dummyTris.empty())
332 {
333 triangles.push_back(t);
334 return TriInd(triangles.size() - 1);
335 }
336 const TriInd nxtDummy = m_dummyTris.back();
337 m_dummyTris.pop_back();
338 triangles[nxtDummy] = t;
339 return nxtDummy;
340}
341
342template <typename T, typename TNearPointLocator>
343TriInd Triangulation<T, TNearPointLocator>::addTriangle()
344{
345 if(m_dummyTris.empty())
346 {
347 const Triangle dummy = {
348 {noVertex, noVertex, noVertex},
349 {noNeighbor, noNeighbor, noNeighbor}};
350 triangles.push_back(dummy);
351 return TriInd(triangles.size() - 1);
352 }
353 const TriInd nxtDummy = m_dummyTris.back();
354 m_dummyTris.pop_back();
355 return nxtDummy;
356}
357
358template <typename T, typename TNearPointLocator>
360 const std::vector<Edge>& edges)
361{
362 insertEdges(edges.begin(), edges.end(), edge_get_v1, edge_get_v2);
363}
364
365template <typename T, typename TNearPointLocator>
367 const std::vector<Edge>& edges)
368{
369 conformToEdges(edges.begin(), edges.end(), edge_get_v1, edge_get_v2);
370}
371
372template <typename T, typename TNearPointLocator>
374{
375 if(!fixedEdges.insert(edge).second)
376 {
377 ++overlapCount[edge]; // if edge is already fixed increment the counter
378 }
379}
380
381namespace detail
382{
383
384// add element to 'to' if not already in 'to'
385template <typename T, typename Allocator1>
386void insert_unique(std::vector<T, Allocator1>& to, const T& elem)
387{
388 if(std::find(to.begin(), to.end(), elem) == to.end())
389 {
390 to.push_back(elem);
391 }
392}
393
394// add elements of 'from' that are not present in 'to' to 'to'
395template <typename T, typename Allocator1, typename Allocator2>
396void insert_unique(
397 std::vector<T, Allocator1>& to,
398 const std::vector<T, Allocator2>& from)
399{
400 typedef typename std::vector<T, Allocator2>::const_iterator Cit;
401 to.reserve(to.size() + from.size());
402 for(Cit cit = from.begin(); cit != from.end(); ++cit)
403 {
404 insert_unique(to, *cit);
405 }
406}
407
408} // namespace detail
409
410template <typename T, typename TNearPointLocator>
411void Triangulation<T, TNearPointLocator>::splitFixedEdge(
412 const Edge& edge,
413 const VertInd iSplitVert)
414{
415 // split constraint (fixed) edge that already exists in triangulation
416 const Edge half1(edge.v1(), iSplitVert);
417 const Edge half2(iSplitVert, edge.v2());
418 // remove the edge that and add its halves
419 fixedEdges.erase(edge);
420 fixEdge(half1);
421 fixEdge(half2);
422 // maintain overlaps
423 typedef unordered_map<Edge, BoundaryOverlapCount>::const_iterator It;
424 const It overlapIt = overlapCount.find(edge);
425 if(overlapIt != overlapCount.end())
426 {
427 overlapCount[half1] += overlapIt->second;
428 overlapCount[half2] += overlapIt->second;
429 overlapCount.erase(overlapIt);
430 }
431 // maintain piece-to-original mapping
432 EdgeVec newOriginals(1, edge);
433 const unordered_map<Edge, EdgeVec>::const_iterator originalsIt =
434 pieceToOriginals.find(edge);
435 if(originalsIt != pieceToOriginals.end())
436 { // edge being split was split before: pass-through originals
437 newOriginals = originalsIt->second;
438 pieceToOriginals.erase(originalsIt);
439 }
440 detail::insert_unique(pieceToOriginals[half1], newOriginals);
441 detail::insert_unique(pieceToOriginals[half2], newOriginals);
442}
443
444template <typename T, typename TNearPointLocator>
445VertInd Triangulation<T, TNearPointLocator>::addSplitEdgeVertex(
446 const V2d<T>& splitVert,
447 const TriInd iT,
448 const TriInd iTopo)
449{
450 // add a new point on the edge that splits an edge in two
451 const VertInd iSplitVert = static_cast<VertInd>(vertices.size());
452 addNewVertex(splitVert, noNeighbor);
453 std::stack<TriInd> triStack = insertVertexOnEdge(iSplitVert, iT, iTopo);
454 tryAddVertexToLocator(iSplitVert);
455 ensureDelaunayByEdgeFlips(iSplitVert, triStack);
456 return iSplitVert;
457}
458
459template <typename T, typename TNearPointLocator>
460VertInd Triangulation<T, TNearPointLocator>::splitFixedEdgeAt(
461 const Edge& edge,
462 const V2d<T>& splitVert,
463 const TriInd iT,
464 const TriInd iTopo)
465{
466 const VertInd iSplitVert = addSplitEdgeVertex(splitVert, iT, iTopo);
467 splitFixedEdge(edge, iSplitVert);
468 return iSplitVert;
469}
470
471template <typename T, typename TNearPointLocator>
472void Triangulation<T, TNearPointLocator>::fixEdge(
473 const Edge& edge,
474 const Edge& originalEdge)
475{
476 fixEdge(edge);
477 if(edge != originalEdge)
478 detail::insert_unique(pieceToOriginals[edge], originalEdge);
479}
480
481namespace detail
482{
483
484template <typename T>
485T lerp(const T& a, const T& b, const T t)
486{
487 return (T(1) - t) * a + t * b;
488}
489
490// Precondition: ab and cd intersect normally
491template <typename T>
492V2d<T> intersectionPosition(
493 const V2d<T>& a,
494 const V2d<T>& b,
495 const V2d<T>& c,
496 const V2d<T>& d)
497{
498 using namespace predicates::adaptive;
499
500 // note: for better accuracy we interpolate x and y separately
501 // on a segment with the shortest x/y-projection correspondingly
502 const T a_cd = orient2d(c.x, c.y, d.x, d.y, a.x, a.y);
503 const T b_cd = orient2d(c.x, c.y, d.x, d.y, b.x, b.y);
504 const T t_ab = a_cd / (a_cd - b_cd);
505
506 const T c_ab = orient2d(a.x, a.y, b.x, b.y, c.x, c.y);
507 const T d_ab = orient2d(a.x, a.y, b.x, b.y, d.x, d.y);
508 const T t_cd = c_ab / (c_ab - d_ab);
509
510 return V2d<T>::make(
511 std::fabs(a.x - b.x) < std::fabs(c.x - d.x) ? lerp(a.x, b.x, t_ab)
512 : lerp(c.x, d.x, t_cd),
513 std::fabs(a.y - b.y) < std::fabs(c.y - d.y) ? lerp(a.y, b.y, t_ab)
514 : lerp(c.y, d.y, t_cd));
515}
516
517} // namespace detail
518
519template <typename T, typename TNearPointLocator>
520void Triangulation<T, TNearPointLocator>::insertEdgeIteration(
521 const Edge edge,
522 const Edge originalEdge,
523 EdgeVec& remaining,
524 std::vector<TriangulatePseudoPolygonTask>& tppIterations)
525{
526 const VertInd iA = edge.v1();
527 VertInd iB = edge.v2();
528 if(iA == iB) // edge connects a vertex to itself
529 return;
530
531 if(hasEdge(iA, iB))
532 {
533 fixEdge(edge, originalEdge);
534 return;
535 }
536
537 const V2d<T>& a = vertices[iA];
538 const V2d<T>& b = vertices[iB];
539 const T distanceTolerance =
540 m_minDistToConstraintEdge == T(0)
541 ? T(0)
542 : m_minDistToConstraintEdge * distance(a, b);
543
544 TriInd iT;
545 // Note: 'L' is left and 'R' is right of the inserted constraint edge
546 VertInd iVL, iVR;
547 tie(iT, iVL, iVR) = intersectedTriangle(iA, a, b, distanceTolerance);
548 // if one of the triangle vertices is on the edge, move edge start
549 if(iT == noNeighbor)
550 {
551 const Edge edgePart(iA, iVL);
552 fixEdge(edgePart, originalEdge);
553 remaining.push_back(Edge(iVL, iB));
554 return;
555 }
556 Triangle t = triangles[iT];
557 std::vector<TriInd> intersected(1, iT);
558 std::vector<VertInd> polyL, polyR;
559 polyL.reserve(2);
560 polyL.push_back(iA);
561 polyL.push_back(iVL);
562 polyR.reserve(2);
563 polyR.push_back(iA);
564 polyR.push_back(iVR);
565 unordered_map<Edge, TriInd> outerTris;
566 outerTris[Edge(iA, iVL)] = edgeNeighbor(t, iA, iVL);
567 outerTris[Edge(iA, iVR)] = edgeNeighbor(t, iA, iVR);
568 VertInd iV = iA;
569
570 while(!t.containsVertex(iB))
571 {
572 const TriInd iTopo = opposedTriangle(t, iV);
573 const Triangle& tOpo = triangles[iTopo];
574 const VertInd iVopo = opposedVertex(tOpo, iT);
575
576 switch(m_intersectingEdgesStrategy)
577 {
579 if(fixedEdges.count(Edge(iVL, iVR)))
580 {
581 // make sure to report original input edges in the exception
582 Edge e1 = originalEdge;
583 Edge e2 = Edge(iVL, iVR);
584 e2 = pieceToOriginals.count(e2)
585 ? pieceToOriginals.at(e2).front()
586 : e2;
587 // don't count super-triangle vertices
588 e1 = Edge(e1.v1() - m_nTargetVerts, e1.v2() - m_nTargetVerts);
589 e2 = Edge(e2.v1() - m_nTargetVerts, e2.v2() - m_nTargetVerts);
590 throw IntersectingConstraintsError(
591 e1,
592 pieceToOriginals.count(e2) ? pieceToOriginals.at(e2).front()
593 : e2,
595 }
596 break;
598 if(!fixedEdges.count(Edge(iVL, iVR)))
599 break;
600 // split edge at the intersection of two constraint edges
601 const V2d<T> newV = detail::intersectionPosition(
602 vertices[iA], vertices[iB], vertices[iVL], vertices[iVR]);
603 const VertInd iNewVert =
604 splitFixedEdgeAt(Edge(iVL, iVR), newV, iT, iTopo);
605 // TODO: is it's possible to re-use pseudo-polygons
606 // for inserting [iA, iNewVert] edge half?
607 remaining.push_back(Edge(iA, iNewVert));
608 remaining.push_back(Edge(iNewVert, iB));
609 return;
610 }
612 assert(!fixedEdges.count(Edge(iVL, iVR)));
613 break;
614 }
615
616 const PtLineLocation::Enum loc =
617 locatePointLine(vertices[iVopo], a, b, distanceTolerance);
618 if(loc == PtLineLocation::Left)
619 {
620 const Edge e(polyL.back(), iVopo);
621 const TriInd outer = edgeNeighbor(tOpo, e.v1(), e.v2());
622 if(!outerTris.insert(std::make_pair(e, outer)).second)
623 outerTris.at(e) = noNeighbor; // hanging edge detected
624 polyL.push_back(iVopo);
625 iV = iVL;
626 iVL = iVopo;
627 }
628 else if(loc == PtLineLocation::Right)
629 {
630 const Edge e(polyR.back(), iVopo);
631 const TriInd outer = edgeNeighbor(tOpo, e.v1(), e.v2());
632 if(!outerTris.insert(std::make_pair(e, outer)).second)
633 outerTris.at(e) = noNeighbor; // hanging edge detected
634 polyR.push_back(iVopo);
635 iV = iVR;
636 iVR = iVopo;
637 }
638 else // encountered point on the edge
639 iB = iVopo;
640
641 intersected.push_back(iTopo);
642 iT = iTopo;
643 t = triangles[iT];
644 }
645 outerTris[Edge(polyL.back(), iB)] = edgeNeighbor(t, polyL.back(), iB);
646 outerTris[Edge(polyR.back(), iB)] = edgeNeighbor(t, polyR.back(), iB);
647 polyL.push_back(iB);
648 polyR.push_back(iB);
649
650 assert(!intersected.empty());
651 // make sure start/end vertices have a valid adjacent triangle
652 // that is not intersected by an edge
653 if(m_vertTris[iA] == intersected.front())
654 pivotVertexTriangleCW(iA);
655 if(m_vertTris[iB] == intersected.back())
656 pivotVertexTriangleCW(iB);
657 // Remove intersected triangles
658 typedef std::vector<TriInd>::const_iterator TriIndCit;
659 for(TriIndCit it = intersected.begin(); it != intersected.end(); ++it)
660 makeDummy(*it);
661 { // Triangulate pseudo-polygons on both sides
662 std::reverse(polyR.begin(), polyR.end());
663 const TriInd iTL = addTriangle();
664 const TriInd iTR = addTriangle();
665 triangulatePseudoPolygon(polyL, outerTris, iTL, iTR, tppIterations);
666 triangulatePseudoPolygon(polyR, outerTris, iTR, iTL, tppIterations);
667 }
668
669 if(iB != edge.v2()) // encountered point on the edge
670 {
671 // fix edge part
672 const Edge edgePart(iA, iB);
673 fixEdge(edgePart, originalEdge);
674 remaining.push_back(Edge(iB, edge.v2()));
675 return;
676 }
677 else
678 {
679 fixEdge(edge, originalEdge);
680 }
681}
682
683template <typename T, typename TNearPointLocator>
684void Triangulation<T, TNearPointLocator>::insertEdge(
685 Edge edge,
686 const Edge originalEdge,
687 EdgeVec& remaining,
688 std::vector<TriangulatePseudoPolygonTask>& tppIterations)
689{
690 // use iteration over recursion to avoid stack overflows
691 remaining.clear();
692 remaining.push_back(edge);
693 while(!remaining.empty())
694 {
695 edge = remaining.back();
696 remaining.pop_back();
697 insertEdgeIteration(edge, originalEdge, remaining, tppIterations);
698 }
699}
700
701template <typename T, typename TNearPointLocator>
702void Triangulation<T, TNearPointLocator>::conformToEdgeIteration(
703 Edge edge,
704 const EdgeVec& originals,
705 BoundaryOverlapCount overlaps,
706 std::vector<ConformToEdgeTask>& remaining)
707{
708 const VertInd iA = edge.v1();
709 VertInd iB = edge.v2();
710 if(iA == iB) // edge connects a vertex to itself
711 return;
712
713 if(hasEdge(iA, iB))
714 {
715 fixEdge(edge);
716 if(overlaps > 0)
717 overlapCount[edge] = overlaps;
718 // avoid marking edge as a part of itself
719 if(!originals.empty() && edge != originals.front())
720 {
721 detail::insert_unique(pieceToOriginals[edge], originals);
722 }
723 return;
724 }
725
726 const V2d<T>& a = vertices[iA];
727 const V2d<T>& b = vertices[iB];
728 const T distanceTolerance =
729 m_minDistToConstraintEdge == T(0)
730 ? T(0)
731 : m_minDistToConstraintEdge * distance(a, b);
732 TriInd iT;
733 VertInd iVleft, iVright;
734 tie(iT, iVleft, iVright) = intersectedTriangle(iA, a, b, distanceTolerance);
735 // if one of the triangle vertices is on the edge, move edge start
736 if(iT == noNeighbor)
737 {
738 const Edge edgePart(iA, iVleft);
739 fixEdge(edgePart);
740 if(overlaps > 0)
741 overlapCount[edgePart] = overlaps;
742 detail::insert_unique(pieceToOriginals[edgePart], originals);
743#ifdef CDT_CXX11_IS_SUPPORTED
744 remaining.emplace_back(Edge(iVleft, iB), originals, overlaps);
745#else
746 remaining.push_back(make_tuple(Edge(iVleft, iB), originals, overlaps));
747#endif
748 return;
749 }
750
751 VertInd iV = iA;
752 Triangle t = triangles[iT];
753 while(std::find(t.vertices.begin(), t.vertices.end(), iB) ==
754 t.vertices.end())
755 {
756 const TriInd iTopo = opposedTriangle(t, iV);
757 const Triangle& tOpo = triangles[iTopo];
758 const VertInd iVopo = opposedVertex(tOpo, iT);
759 const V2d<T> vOpo = vertices[iVopo];
760
761 switch(m_intersectingEdgesStrategy)
762 {
764 if(fixedEdges.count(Edge(iVleft, iVright)))
765 {
766 // make sure to report original input edges in the exception
767 Edge e1 = pieceToOriginals.count(edge)
768 ? pieceToOriginals.at(edge).front()
769 : edge;
770 Edge e2(iVleft, iVright);
771 e2 = pieceToOriginals.count(e2)
772 ? pieceToOriginals.at(e2).front()
773 : e2;
774 // don't count super-triangle vertices
775 e1 = Edge(e1.v1() - m_nTargetVerts, e1.v2() - m_nTargetVerts);
776 e2 = Edge(e2.v1() - m_nTargetVerts, e2.v2() - m_nTargetVerts);
777 throw IntersectingConstraintsError(e1, e2, CDT_SOURCE_LOCATION);
778 }
779 break;
781 if(!fixedEdges.count(Edge(iVleft, iVright)))
782 break;
783 // split edge at the intersection of two constraint edges
784 const V2d<T> newV = detail::intersectionPosition(
785 vertices[iA],
786 vertices[iB],
787 vertices[iVleft],
788 vertices[iVright]);
789 const VertInd iNewVert =
790 splitFixedEdgeAt(Edge(iVleft, iVright), newV, iT, iTopo);
791#ifdef CDT_CXX11_IS_SUPPORTED
792 remaining.emplace_back(Edge(iNewVert, iB), originals, overlaps);
793 remaining.emplace_back(Edge(iA, iNewVert), originals, overlaps);
794#else
795 remaining.push_back(
796 make_tuple(Edge(iNewVert, iB), originals, overlaps));
797 remaining.push_back(
798 make_tuple(Edge(iA, iNewVert), originals, overlaps));
799#endif
800 return;
801 }
803 assert(!fixedEdges.count(Edge(iVleft, iVright)));
804 break;
805 }
806
807 iT = iTopo;
808 t = triangles[iT];
809
810 const PtLineLocation::Enum loc =
811 locatePointLine(vOpo, a, b, distanceTolerance);
812 if(loc == PtLineLocation::Left)
813 {
814 iV = iVleft;
815 iVleft = iVopo;
816 }
817 else if(loc == PtLineLocation::Right)
818 {
819 iV = iVright;
820 iVright = iVopo;
821 }
822 else // encountered point on the edge
823 iB = iVopo;
824 }
825
826 // encountered one or more points on the edge: add remaining edge part
827 if(iB != edge.v2())
828 {
829#ifdef CDT_CXX11_IS_SUPPORTED
830 remaining.emplace_back(Edge(iB, edge.v2()), originals, overlaps);
831#else
832 remaining.push_back(
833 make_tuple(Edge(iB, edge.v2()), originals, overlaps));
834#endif
835 }
836
837 // add mid-point to triangulation
838 const VertInd iMid = static_cast<VertInd>(vertices.size());
839 const V2d<T>& start = vertices[iA];
840 const V2d<T>& end = vertices[iB];
841 addNewVertex(
842 V2d<T>::make((start.x + end.x) / T(2), (start.y + end.y) / T(2)),
843 noNeighbor);
844 const std::vector<Edge> flippedFixedEdges =
845 insertVertex_FlipFixedEdges(iMid);
846
847#ifdef CDT_CXX11_IS_SUPPORTED
848 remaining.emplace_back(Edge(iMid, iB), originals, overlaps);
849 remaining.emplace_back(Edge(iA, iMid), originals, overlaps);
850#else
851 remaining.push_back(make_tuple(Edge(iMid, iB), originals, overlaps));
852 remaining.push_back(make_tuple(Edge(iA, iMid), originals, overlaps));
853#endif
854
855 // re-introduce fixed edges that were flipped
856 // and make sure overlap count is preserved
857 for(std::vector<Edge>::const_iterator it = flippedFixedEdges.begin();
858 it != flippedFixedEdges.end();
859 ++it)
860 {
861 const Edge& flippedFixedEdge = *it;
862 fixedEdges.erase(flippedFixedEdge);
863
864 BoundaryOverlapCount prevOverlaps = 0;
865 const unordered_map<Edge, BoundaryOverlapCount>::const_iterator
866 overlapsIt = overlapCount.find(flippedFixedEdge);
867 if(overlapsIt != overlapCount.end())
868 {
869 prevOverlaps = overlapsIt->second;
870 overlapCount.erase(overlapsIt);
871 }
872 // override overlapping boundaries count when re-inserting an edge
873 EdgeVec prevOriginals(1, flippedFixedEdge);
874 const unordered_map<Edge, EdgeVec>::const_iterator originalsIt =
875 pieceToOriginals.find(flippedFixedEdge);
876 if(originalsIt != pieceToOriginals.end())
877 {
878 prevOriginals = originalsIt->second;
879 }
880#ifdef CDT_CXX11_IS_SUPPORTED
881 remaining.emplace_back(flippedFixedEdge, prevOriginals, prevOverlaps);
882#else
883 remaining.push_back(
884 make_tuple(flippedFixedEdge, prevOriginals, prevOverlaps));
885#endif
886 }
887}
888
889template <typename T, typename TNearPointLocator>
890void Triangulation<T, TNearPointLocator>::conformToEdge(
891 Edge edge,
892 EdgeVec originals,
893 BoundaryOverlapCount overlaps,
894 std::vector<ConformToEdgeTask>& remaining)
895{
896 // use iteration over recursion to avoid stack overflows
897 remaining.clear();
898#ifdef CDT_CXX11_IS_SUPPORTED
899 remaining.emplace_back(edge, originals, overlaps);
900#else
901 remaining.push_back(make_tuple(edge, originals, overlaps));
902#endif
903 while(!remaining.empty())
904 {
905 tie(edge, originals, overlaps) = remaining.back();
906 remaining.pop_back();
907 conformToEdgeIteration(edge, originals, overlaps, remaining);
908 }
909}
910
921template <typename T, typename TNearPointLocator>
922tuple<TriInd, VertInd, VertInd>
923Triangulation<T, TNearPointLocator>::intersectedTriangle(
924 const VertInd iA,
925 const V2d<T>& a,
926 const V2d<T>& b,
927 const T orientationTolerance) const
928{
929 const TriInd startTri = m_vertTris[iA];
930 TriInd iT = startTri;
931 do
932 {
933 const Triangle t = triangles[iT];
934 const Index i = vertexInd(t.vertices, iA);
935 const VertInd iP2 = t.vertices[ccw(i)];
936 const T orientP2 = orient2D(vertices[iP2], a, b);
937 const PtLineLocation::Enum locP2 = classifyOrientation(orientP2);
938 if(locP2 == PtLineLocation::Right)
939 {
940 const VertInd iP1 = t.vertices[cw(i)];
941 const T orientP1 = orient2D(vertices[iP1], a, b);
942 const PtLineLocation::Enum locP1 = classifyOrientation(orientP1);
943 if(locP1 == PtLineLocation::OnLine)
944 {
945 return make_tuple(noNeighbor, iP1, iP1);
946 }
947 if(locP1 == PtLineLocation::Left)
948 {
949 if(orientationTolerance)
950 {
951 T closestOrient;
952 VertInd iClosestP;
953 if(std::abs(orientP1) <= std::abs(orientP2))
954 {
955 closestOrient = orientP1;
956 iClosestP = iP1;
957 }
958 else
959 {
960 closestOrient = orientP2;
961 iClosestP = iP2;
962 }
964 closestOrient, orientationTolerance) ==
965 PtLineLocation::OnLine)
966 {
967 return make_tuple(noNeighbor, iClosestP, iClosestP);
968 }
969 }
970 return make_tuple(iT, iP1, iP2);
971 }
972 }
973 iT = t.next(iA).first;
974 } while(iT != startTri);
975
976 throw Error(
977 "Could not find vertex triangle intersected by an edge.",
979}
980
981template <typename T, typename TNearPointLocator>
982void Triangulation<T, TNearPointLocator>::addSuperTriangle(const Box2d<T>& box)
983{
984 m_nTargetVerts = 3;
985 m_superGeomType = SuperGeometryType::SuperTriangle;
986
987 const V2d<T> center = {
988 (box.min.x + box.max.x) / T(2), (box.min.y + box.max.y) / T(2)};
989 const T w = box.max.x - box.min.x;
990 const T h = box.max.y - box.min.y;
991 T r = std::max(w, h); // incircle radius upper bound
992
993 // Note: make sure radius is big enough. Constants chosen experimentally.
994 // - for tiny bounding boxes: use 1.0 as the smallest radius
995 // - multiply radius by 2.0 for extra safety margin
996 r = std::max(T(2) * r, T(1));
997
998 // Note: for very large floating point numbers rounding can lead to wrong
999 // super-triangle coordinates. This is a very rare corner-case so the
1000 // handling is very primitive.
1001 { // note: '<=' means '==' but avoids the warning
1002 while(center.y <= center.y - r)
1003 r *= T(2);
1004 }
1005
1006 const T R = T(2) * r; // excircle radius
1007 const T cos_30_deg = 0.8660254037844386; // note: (std::sqrt(3.0) / 2.0)
1008 const T shiftX = R * cos_30_deg;
1009 const V2d<T> posV1 = {center.x - shiftX, center.y - r};
1010 const V2d<T> posV2 = {center.x + shiftX, center.y - r};
1011 const V2d<T> posV3 = {center.x, center.y + R};
1012 addNewVertex(posV1, TriInd(0));
1013 addNewVertex(posV2, TriInd(0));
1014 addNewVertex(posV3, TriInd(0));
1015 const Triangle superTri = {
1016 {VertInd(0), VertInd(1), VertInd(2)},
1017 {noNeighbor, noNeighbor, noNeighbor}};
1018 addTriangle(superTri);
1019 if(m_vertexInsertionOrder != VertexInsertionOrder::Auto)
1020 {
1021 m_nearPtLocator.initialize(vertices);
1022 }
1023}
1024
1025template <typename T, typename TNearPointLocator>
1026void Triangulation<T, TNearPointLocator>::addNewVertex(
1027 const V2d<T>& pos,
1028 const TriInd iT)
1029{
1030 vertices.push_back(pos);
1031 m_vertTris.push_back(iT);
1032}
1033
1034template <typename T, typename TNearPointLocator>
1035std::vector<Edge>
1036Triangulation<T, TNearPointLocator>::insertVertex_FlipFixedEdges(
1037 const VertInd iV1)
1038{
1039 std::vector<Edge> flippedFixedEdges;
1040
1041 const V2d<T>& v1 = vertices[iV1];
1042 const VertInd startVertex = m_nearPtLocator.nearPoint(v1, vertices);
1043 array<TriInd, 2> trisAt = walkingSearchTrianglesAt(iV1, startVertex);
1044 std::stack<TriInd> triStack =
1045 trisAt[1] == noNeighbor ? insertVertexInsideTriangle(iV1, trisAt[0])
1046 : insertVertexOnEdge(iV1, trisAt[0], trisAt[1]);
1047
1048 TriInd iTopo, n1, n2, n3, n4;
1049 VertInd iV2, iV3, iV4;
1050 while(!triStack.empty())
1051 {
1052 const TriInd iT = triStack.top();
1053 triStack.pop();
1054
1055 edgeFlipInfo(iT, iV1, iTopo, iV2, iV3, iV4, n1, n2, n3, n4);
1056 if(iTopo != noNeighbor && isFlipNeeded(iV1, iV2, iV3, iV4))
1057 {
1058 // if flipped edge is fixed, remember it
1059 const Edge flippedEdge(iV2, iV4);
1060 if(!fixedEdges.empty() &&
1061 fixedEdges.find(flippedEdge) != fixedEdges.end())
1062 {
1063 flippedFixedEdges.push_back(flippedEdge);
1064 }
1065
1066 flipEdge(iT, iTopo, iV1, iV2, iV3, iV4, n1, n2, n3, n4);
1067 triStack.push(iT);
1068 triStack.push(iTopo);
1069 }
1070 }
1071
1072 tryAddVertexToLocator(iV1);
1073 return flippedFixedEdges;
1074}
1075
1076template <typename T, typename TNearPointLocator>
1077void Triangulation<T, TNearPointLocator>::insertVertex(
1078 const VertInd iVert,
1079 const VertInd walkStart)
1080{
1081 const array<TriInd, 2> trisAt = walkingSearchTrianglesAt(iVert, walkStart);
1082 std::stack<TriInd> triStack =
1083 trisAt[1] == noNeighbor
1084 ? insertVertexInsideTriangle(iVert, trisAt[0])
1085 : insertVertexOnEdge(iVert, trisAt[0], trisAt[1]);
1086 ensureDelaunayByEdgeFlips(iVert, triStack);
1087}
1088
1089template <typename T, typename TNearPointLocator>
1090void Triangulation<T, TNearPointLocator>::insertVertex(const VertInd iVert)
1091{
1092 const V2d<T>& v = vertices[iVert];
1093 const VertInd walkStart = m_nearPtLocator.nearPoint(v, vertices);
1094 insertVertex(iVert, walkStart);
1095 tryAddVertexToLocator(iVert);
1096}
1097
1098template <typename T, typename TNearPointLocator>
1099void Triangulation<T, TNearPointLocator>::ensureDelaunayByEdgeFlips(
1100 const VertInd iV1,
1101 std::stack<TriInd>& triStack)
1102{
1103 TriInd iTopo, n1, n2, n3, n4;
1104 VertInd iV2, iV3, iV4;
1105 while(!triStack.empty())
1106 {
1107 const TriInd iT = triStack.top();
1108 triStack.pop();
1109
1110 edgeFlipInfo(iT, iV1, iTopo, iV2, iV3, iV4, n1, n2, n3, n4);
1111 if(iTopo != noNeighbor && isFlipNeeded(iV1, iV2, iV3, iV4))
1112 {
1113 flipEdge(iT, iTopo, iV1, iV2, iV3, iV4, n1, n2, n3, n4);
1114 triStack.push(iT);
1115 triStack.push(iTopo);
1116 }
1117 }
1118}
1119
1120/*
1121 * v4 original edge: (v1, v3)
1122 * /|\ flip-candidate edge: (v, v2)
1123 * / | \
1124 * n3 / | \ n4
1125 * / | \
1126 * new vertex--> v1 T | Topo v3
1127 * \ | /
1128 * n1 \ | / n2
1129 * \ | /
1130 * \|/
1131 * v2
1132 */
1133template <typename T, typename TNearPointLocator>
1134void Triangulation<T, TNearPointLocator>::edgeFlipInfo(
1135 const TriInd iT,
1136 const VertInd iV1,
1137 TriInd& iTopo,
1138 VertInd& iV2,
1139 VertInd& iV3,
1140 VertInd& iV4,
1141 TriInd& n1,
1142 TriInd& n2,
1143 TriInd& n3,
1144 TriInd& n4)
1145{
1146 /* v[2]
1147 / \
1148 n[2]/ \n[1]
1149 /_____\
1150 v[0] n[0] v[1] */
1151 const Triangle& t = triangles[iT];
1152 if(t.vertices[0] == iV1)
1153 {
1154 iV2 = t.vertices[1];
1155 iV4 = t.vertices[2];
1156 n1 = t.neighbors[0];
1157 n3 = t.neighbors[2];
1158 iTopo = t.neighbors[1];
1159 }
1160 else if(t.vertices[1] == iV1)
1161 {
1162 iV2 = t.vertices[2];
1163 iV4 = t.vertices[0];
1164 n1 = t.neighbors[1];
1165 n3 = t.neighbors[0];
1166 iTopo = t.neighbors[2];
1167 }
1168 else
1169 {
1170 iV2 = t.vertices[0];
1171 iV4 = t.vertices[1];
1172 n1 = t.neighbors[2];
1173 n3 = t.neighbors[1];
1174 iTopo = t.neighbors[0];
1175 }
1176 if(iTopo == noNeighbor)
1177 return;
1178 const Triangle& tOpo = triangles[iTopo];
1179 if(tOpo.neighbors[0] == iT)
1180 {
1181 iV3 = tOpo.vertices[2];
1182 n2 = tOpo.neighbors[1];
1183 n4 = tOpo.neighbors[2];
1184 }
1185 else if(tOpo.neighbors[1] == iT)
1186 {
1187 iV3 = tOpo.vertices[0];
1188 n2 = tOpo.neighbors[2];
1189 n4 = tOpo.neighbors[0];
1190 }
1191 else
1192 {
1193 iV3 = tOpo.vertices[1];
1194 n2 = tOpo.neighbors[0];
1195 n4 = tOpo.neighbors[1];
1196 }
1197}
1198
1209/*
1210 * v4 original edge: (v2, v4)
1211 * /|\ flip-candidate edge: (v1, v3)
1212 * / | \
1213 * / | \
1214 * / | \
1215 * new vertex--> v1 | v3
1216 * \ | /
1217 * \ | /
1218 * \ | /
1219 * \|/
1220 * v2
1221 */
1222template <typename T, typename TNearPointLocator>
1223bool Triangulation<T, TNearPointLocator>::isFlipNeeded(
1224 const VertInd iV1,
1225 const VertInd iV2,
1226 const VertInd iV3,
1227 const VertInd iV4) const
1228{
1229 if(fixedEdges.count(Edge(iV2, iV4)))
1230 return false; // flip not needed if the original edge is fixed
1231 const V2d<T>& v1 = vertices[iV1];
1232 const V2d<T>& v2 = vertices[iV2];
1233 const V2d<T>& v3 = vertices[iV3];
1234 const V2d<T>& v4 = vertices[iV4];
1235 if(m_superGeomType == SuperGeometryType::SuperTriangle)
1236 {
1237 // If flip-candidate edge touches super-triangle in-circumference
1238 // test has to be replaced with orient2d test against the line
1239 // formed by two non-artificial vertices (that don't belong to
1240 // super-triangle)
1241 if(iV1 < 3) // flip-candidate edge touches super-triangle
1242 {
1243 // does original edge also touch super-triangle?
1244 if(iV2 < 3)
1245 return locatePointLine(v2, v3, v4) ==
1246 locatePointLine(v1, v3, v4);
1247 if(iV4 < 3)
1248 return locatePointLine(v4, v2, v3) ==
1249 locatePointLine(v1, v2, v3);
1250 return false; // original edge does not touch super-triangle
1251 }
1252 if(iV3 < 3) // flip-candidate edge touches super-triangle
1253 {
1254 // does original edge also touch super-triangle?
1255 if(iV2 < 3)
1256 {
1257 return locatePointLine(v2, v1, v4) ==
1258 locatePointLine(v3, v1, v4);
1259 }
1260 if(iV4 < 3)
1261 {
1262 return locatePointLine(v4, v2, v1) ==
1263 locatePointLine(v3, v2, v1);
1264 }
1265 return false; // original edge does not touch super-triangle
1266 }
1267 // flip-candidate edge does not touch super-triangle
1268 if(iV2 < 3)
1269 return locatePointLine(v2, v3, v4) == locatePointLine(v1, v3, v4);
1270 if(iV4 < 3)
1271 return locatePointLine(v4, v2, v3) == locatePointLine(v1, v2, v3);
1272 }
1273 return isInCircumcircle(v1, v2, v3, v4);
1274}
1275
1276/* Flip edge between T and Topo:
1277 *
1278 * v4 | - old edge
1279 * /|\ ~ - new edge
1280 * / | \
1281 * n3 / T' \ n4
1282 * / | \
1283 * / | \
1284 * T -> v1 ~~~~~~~~ v3 <- Topo
1285 * \ | /
1286 * \ | /
1287 * n1 \Topo'/ n2
1288 * \ | /
1289 * \|/
1290 * v2
1291 */
1292template <typename T, typename TNearPointLocator>
1294 const TriInd iT,
1295 const TriInd iTopo,
1296 const VertInd v1,
1297 const VertInd v2,
1298 const VertInd v3,
1299 const VertInd v4,
1300 const TriInd n1,
1301 const TriInd n2,
1302 const TriInd n3,
1303 const TriInd n4)
1304{
1305 // change vertices and neighbors
1306 using detail::arr3;
1307 triangles[iT] = Triangle::make(arr3(v4, v1, v3), arr3(n3, iTopo, n4));
1308 triangles[iTopo] = Triangle::make(arr3(v2, v3, v1), arr3(n2, iT, n1));
1309 // adjust neighboring triangles and vertices
1310 changeNeighbor(n1, iT, iTopo);
1311 changeNeighbor(n4, iTopo, iT);
1312 // only adjust adjacent triangles if triangulation is not finalized:
1313 // can happen when called from outside on an already finalized
1314 // triangulation
1315 if(!isFinalized())
1316 {
1317 setAdjacentTriangle(v4, iT);
1318 setAdjacentTriangle(v2, iTopo);
1319 }
1320}
1321
1322/* Insert point into triangle: split into 3 triangles:
1323 * - create 2 new triangles
1324 * - re-use old triangle for the 3rd
1325 * v3
1326 * / | \
1327 * / | \ <-- original triangle (t)
1328 * / | \
1329 * n3 / | \ n2
1330 * /newT2|newT1\
1331 * / v \
1332 * / __/ \__ \
1333 * / __/ \__ \
1334 * / _/ t' \_ \
1335 * v1 ___________________ v2
1336 * n1
1337 */
1338template <typename T, typename TNearPointLocator>
1339std::stack<TriInd>
1340Triangulation<T, TNearPointLocator>::insertVertexInsideTriangle(
1341 VertInd v,
1342 TriInd iT)
1343{
1344 const TriInd iNewT1 = addTriangle();
1345 const TriInd iNewT2 = addTriangle();
1346
1347 Triangle& t = triangles[iT];
1348 const array<VertInd, 3> vv = t.vertices;
1349 const array<TriInd, 3> nn = t.neighbors;
1350 const VertInd v1 = vv[0], v2 = vv[1], v3 = vv[2];
1351 const TriInd n1 = nn[0], n2 = nn[1], n3 = nn[2];
1352 // make two new triangles and convert current triangle to 3rd new
1353 // triangle
1354 using detail::arr3;
1355 triangles[iNewT1] = Triangle::make(arr3(v2, v3, v), arr3(n2, iNewT2, iT));
1356 triangles[iNewT2] = Triangle::make(arr3(v3, v1, v), arr3(n3, iT, iNewT1));
1357 t = Triangle::make(arr3(v1, v2, v), arr3(n1, iNewT1, iNewT2));
1358 // adjust adjacent triangles
1359 setAdjacentTriangle(v, iT);
1360 setAdjacentTriangle(v3, iNewT1);
1361 // change triangle neighbor's neighbors to new triangles
1362 changeNeighbor(n2, iT, iNewT1);
1363 changeNeighbor(n3, iT, iNewT2);
1364 // return newly added triangles
1365 std::stack<TriInd> newTriangles;
1366 newTriangles.push(iT);
1367 newTriangles.push(iNewT1);
1368 newTriangles.push(iNewT2);
1369 return newTriangles;
1370}
1371
1372/* Inserting a point on the edge between two triangles
1373 * T1 (top) v1
1374 * /|\
1375 * n1 / | \ n4
1376 * / | \
1377 * / T1' | Tnew1\
1378 * v2-------v-------v4
1379 * \ T2' | Tnew2/
1380 * \ | /
1381 * n2 \ | / n3
1382 * \|/
1383 * T2 (bottom) v3
1384 */
1385template <typename T, typename TNearPointLocator>
1386std::stack<TriInd> Triangulation<T, TNearPointLocator>::insertVertexOnEdge(
1387 VertInd v,
1388 TriInd iT1,
1389 TriInd iT2)
1390{
1391 const TriInd iTnew1 = addTriangle();
1392 const TriInd iTnew2 = addTriangle();
1393
1394 Triangle& t1 = triangles[iT1];
1395 Triangle& t2 = triangles[iT2];
1396 Index i = opposedVertexInd(t1.neighbors, iT2);
1397 const VertInd v1 = t1.vertices[i];
1398 const VertInd v2 = t1.vertices[ccw(i)];
1399 const TriInd n1 = t1.neighbors[i];
1400 const TriInd n4 = t1.neighbors[cw(i)];
1401 i = opposedVertexInd(t2.neighbors, iT1);
1402 const VertInd v3 = t2.vertices[i];
1403 const VertInd v4 = t2.vertices[ccw(i)];
1404 const TriInd n3 = t2.neighbors[i];
1405 const TriInd n2 = t2.neighbors[cw(i)];
1406 // add new triangles and change existing ones
1407 using detail::arr3;
1408 t1 = Triangle::make(arr3(v, v1, v2), arr3(iTnew1, n1, iT2));
1409 t2 = Triangle::make(arr3(v, v2, v3), arr3(iT1, n2, iTnew2));
1410 triangles[iTnew1] = Triangle::make(arr3(v, v4, v1), arr3(iTnew2, n4, iT1));
1411 triangles[iTnew2] = Triangle::make(arr3(v, v3, v4), arr3(iT2, n3, iTnew1));
1412 // adjust adjacent triangles
1413 setAdjacentTriangle(v, iT1);
1414 setAdjacentTriangle(v4, iTnew1);
1415 // adjust neighboring triangles and vertices
1416 changeNeighbor(n4, iT1, iTnew1);
1417 changeNeighbor(n3, iT2, iTnew2);
1418 // return newly added triangles
1419 std::stack<TriInd> newTriangles;
1420 newTriangles.push(iT1);
1421 newTriangles.push(iTnew2);
1422 newTriangles.push(iT2);
1423 newTriangles.push(iTnew1);
1424 return newTriangles;
1425}
1426
1427template <typename T, typename TNearPointLocator>
1428array<TriInd, 2>
1429Triangulation<T, TNearPointLocator>::trianglesAt(const V2d<T>& pos) const
1430{
1431 array<TriInd, 2> out = {noNeighbor, noNeighbor};
1432 for(TriInd i = TriInd(0); i < TriInd(triangles.size()); ++i)
1433 {
1434 const Triangle& t = triangles[i];
1435 const V2d<T>& v1 = vertices[t.vertices[0]];
1436 const V2d<T>& v2 = vertices[t.vertices[1]];
1437 const V2d<T>& v3 = vertices[t.vertices[2]];
1438 const PtTriLocation::Enum loc = locatePointTriangle(pos, v1, v2, v3);
1439 if(loc == PtTriLocation::Outside)
1440 continue;
1441 out[0] = i;
1442 if(isOnEdge(loc))
1443 out[1] = t.neighbors[edgeNeighbor(loc)];
1444 return out;
1445 }
1446 throw Error("No triangle was found at position", CDT_SOURCE_LOCATION);
1447}
1448
1449template <typename T, typename TNearPointLocator>
1450TriInd Triangulation<T, TNearPointLocator>::walkTriangles(
1451 const VertInd startVertex,
1452 const V2d<T>& pos) const
1453{
1454 // begin walk in search of triangle at pos
1455 TriInd currTri = m_vertTris[startVertex];
1456 bool found = false;
1457 detail::SplitMix64RandGen prng;
1458 while(!found)
1459 {
1460 const Triangle& t = triangles[currTri];
1461 found = true;
1462 // stochastic offset to randomize which edge we check first
1463 const Index offset(prng() % 3);
1464 for(Index i_(0); i_ < Index(3); ++i_)
1465 {
1466 const Index i((i_ + offset) % 3);
1467 const V2d<T>& vStart = vertices[t.vertices[i]];
1468 const V2d<T>& vEnd = vertices[t.vertices[ccw(i)]];
1469 const PtLineLocation::Enum edgeCheck =
1470 locatePointLine(pos, vStart, vEnd);
1471 const TriInd iN = t.neighbors[i];
1472 if(edgeCheck == PtLineLocation::Right && iN != noNeighbor)
1473 {
1474 found = false;
1475 currTri = iN;
1476 break;
1477 }
1478 }
1479 }
1480 return currTri;
1481}
1482
1483template <typename T, typename TNearPointLocator>
1484array<TriInd, 2> Triangulation<T, TNearPointLocator>::walkingSearchTrianglesAt(
1485 const VertInd iV,
1486 const VertInd startVertex) const
1487{
1488 const V2d<T> v = vertices[iV];
1489 array<TriInd, 2> out = {noNeighbor, noNeighbor};
1490 const TriInd iT = walkTriangles(startVertex, v);
1491 // Finished walk, locate point in current triangle
1492 const Triangle& t = triangles[iT];
1493 const V2d<T>& v1 = vertices[t.vertices[0]];
1494 const V2d<T>& v2 = vertices[t.vertices[1]];
1495 const V2d<T>& v3 = vertices[t.vertices[2]];
1496 const PtTriLocation::Enum loc = locatePointTriangle(v, v1, v2, v3);
1497
1498 if(loc == PtTriLocation::Outside)
1499 throw Error("No triangle was found at position", CDT_SOURCE_LOCATION);
1500 if(loc == PtTriLocation::OnVertex)
1501 {
1502 const VertInd iDupe = v1 == v ? t.vertices[0]
1503 : v2 == v ? t.vertices[1]
1504 : t.vertices[2];
1505 throw DuplicateVertexError(
1506 iV - m_nTargetVerts, iDupe - m_nTargetVerts, CDT_SOURCE_LOCATION);
1507 }
1508
1509 out[0] = iT;
1510 if(isOnEdge(loc))
1511 out[1] = t.neighbors[edgeNeighbor(loc)];
1512 return out;
1513}
1514
1515/* Flip edge between T and Topo:
1516 *
1517 * v4 | - old edge
1518 * /|\ ~ - new edge
1519 * / | \
1520 * n3 / T' \ n4
1521 * / | \
1522 * / | \
1523 * T -> v1~~~~~~~~~v3 <- Topo
1524 * \ | /
1525 * \ | /
1526 * n1 \Topo'/ n2
1527 * \ | /
1528 * \|/
1529 * v2
1530 */
1531template <typename T, typename TNearPointLocator>
1533 const TriInd iT,
1534 const TriInd iTopo)
1535{
1536 Triangle& t = triangles[iT];
1537 Triangle& tOpo = triangles[iTopo];
1538 const array<TriInd, 3>& triNs = t.neighbors;
1539 const array<TriInd, 3>& triOpoNs = tOpo.neighbors;
1540 const array<VertInd, 3>& triVs = t.vertices;
1541 const array<VertInd, 3>& triOpoVs = tOpo.vertices;
1542 // find vertices and neighbors
1543 Index i = opposedVertexInd(t.neighbors, iTopo);
1544 const VertInd v1 = triVs[i];
1545 const VertInd v2 = triVs[ccw(i)];
1546 const TriInd n1 = triNs[i];
1547 const TriInd n3 = triNs[cw(i)];
1548 i = opposedVertexInd(tOpo.neighbors, iT);
1549 const VertInd v3 = triOpoVs[i];
1550 const VertInd v4 = triOpoVs[ccw(i)];
1551 const TriInd n4 = triOpoNs[i];
1552 const TriInd n2 = triOpoNs[cw(i)];
1553 // change vertices and neighbors
1554 using detail::arr3;
1555 t = Triangle::make(arr3(v4, v1, v3), arr3(n3, iTopo, n4));
1556 tOpo = Triangle::make(arr3(v2, v3, v1), arr3(n2, iT, n1));
1557 // adjust neighboring triangles and vertices
1558 changeNeighbor(n1, iT, iTopo);
1559 changeNeighbor(n4, iTopo, iT);
1560 // only adjust adjacent triangles if triangulation is not finalized:
1561 // can happen when called from outside on an already finalized
1562 // triangulation
1563 if(!isFinalized())
1564 {
1565 setAdjacentTriangle(v4, iT);
1566 setAdjacentTriangle(v2, iTopo);
1567 }
1568}
1569
1570template <typename T, typename TNearPointLocator>
1572 const TriInd iT,
1573 const TriInd oldNeighbor,
1574 const TriInd newNeighbor)
1575{
1576 if(iT == noNeighbor)
1577 return;
1578 NeighborsArr3& nn = triangles[iT].neighbors;
1579 assert(
1580 nn[0] == oldNeighbor || nn[1] == oldNeighbor || nn[2] == oldNeighbor);
1581 if(nn[0] == oldNeighbor)
1582 nn[0] = newNeighbor;
1583 else if(nn[1] == oldNeighbor)
1584 nn[1] = newNeighbor;
1585 else
1586 nn[2] = newNeighbor;
1587}
1588
1589template <typename T, typename TNearPointLocator>
1590void Triangulation<T, TNearPointLocator>::changeNeighbor(
1591 const TriInd iT,
1592 const VertInd iVedge1,
1593 const VertInd iVedge2,
1594 const TriInd newNeighbor)
1595{
1596 assert(iT != noNeighbor);
1597 Triangle& t = triangles[iT];
1598 t.neighbors[edgeNeighborInd(t.vertices, iVedge1, iVedge2)] = newNeighbor;
1599}
1600
1601template <typename T, typename TNearPointLocator>
1602void Triangulation<T, TNearPointLocator>::triangulatePseudoPolygon(
1603 const std::vector<VertInd>& poly,
1604 unordered_map<Edge, TriInd>& outerTris,
1605 TriInd iT,
1606 TriInd iN,
1607 std::vector<TriangulatePseudoPolygonTask>& iterations)
1608{
1609 assert(poly.size() > 2);
1610 // note: uses iteration instead of recursion to avoid stack overflows
1611 iterations.clear();
1612 iterations.push_back(make_tuple(
1613 IndexSizeType(0),
1614 static_cast<IndexSizeType>(poly.size() - 1),
1615 iT,
1616 iN,
1617 Index(0)));
1618 while(!iterations.empty())
1619 {
1620 triangulatePseudoPolygonIteration(poly, outerTris, iterations);
1621 }
1622}
1623
1624template <typename T, typename TNearPointLocator>
1625void Triangulation<T, TNearPointLocator>::triangulatePseudoPolygonIteration(
1626 const std::vector<VertInd>& poly,
1627 unordered_map<Edge, TriInd>& outerTris,
1628 std::vector<TriangulatePseudoPolygonTask>& iterations)
1629{
1630 IndexSizeType iA, iB;
1631 TriInd iT, iParent;
1632 Index iInParent;
1633 assert(!iterations.empty());
1634 tie(iA, iB, iT, iParent, iInParent) = iterations.back();
1635 iterations.pop_back();
1636 assert(iB - iA > 1 && iT != noNeighbor && iParent != noNeighbor);
1637 Triangle& t = triangles[iT];
1638 // find Delaunay point
1639 const IndexSizeType iC = findDelaunayPoint(poly, iA, iB);
1640
1641 const VertInd a = poly[iA];
1642 const VertInd b = poly[iB];
1643 const VertInd c = poly[iC];
1644
1645 // split pseudo-polygon in two parts and triangulate them
1646 // note: second part needs to be pushed on stack first to be processed first
1647
1648 // second part: points after the Delaunay point
1649 if(iB - iC > 1)
1650 {
1651 const TriInd iNext = addTriangle();
1652 iterations.push_back(make_tuple(iC, iB, iNext, iT, Index(1)));
1653 }
1654 else // pseudo-poly is reduced to a single outer edge
1655 {
1656 const Edge outerEdge(b, c);
1657 const TriInd outerTri = outerTris.at(outerEdge);
1658 if(outerTri != noNeighbor)
1659 {
1660 assert(outerTri != iT);
1661 t.neighbors[1] = outerTri;
1662 changeNeighbor(outerTri, c, b, iT);
1663 }
1664 else
1665 outerTris.at(outerEdge) = iT;
1666 }
1667 // first part: points before the Delaunay point
1668 if(iC - iA > 1)
1669 { // add next triangle and add another iteration
1670 const TriInd iNext = addTriangle();
1671 iterations.push_back(make_tuple(iA, iC, iNext, iT, Index(2)));
1672 }
1673 else
1674 { // pseudo-poly is reduced to a single outer edge
1675 const Edge outerEdge(c, a);
1676 const TriInd outerTri = outerTris.at(outerEdge);
1677 if(outerTri != noNeighbor)
1678 {
1679 assert(outerTri != iT);
1680 t.neighbors[2] = outerTri;
1681 changeNeighbor(outerTri, c, a, iT);
1682 }
1683 else
1684 outerTris.at(outerEdge) = iT;
1685 }
1686 // Finalize triangle
1687 // note: only when triangle is finalized to we add it as a neighbor to
1688 // parent to maintain triangulation topology consistency
1689 triangles[iParent].neighbors[iInParent] = iT;
1690 t.neighbors[0] = iParent;
1691 t.vertices = detail::arr3(a, b, c);
1692 setAdjacentTriangle(c, iT);
1693}
1694
1695template <typename T, typename TNearPointLocator>
1696IndexSizeType Triangulation<T, TNearPointLocator>::findDelaunayPoint(
1697 const std::vector<VertInd>& poly,
1698 const IndexSizeType iA,
1699 const IndexSizeType iB) const
1700{
1701 assert(iB - iA > 1);
1702 const V2d<T>& a = vertices[poly[iA]];
1703 const V2d<T>& b = vertices[poly[iB]];
1704 IndexSizeType out = iA + 1;
1705 const V2d<T>* c = &vertices[poly[out]]; // caching for better performance
1706 for(IndexSizeType i = iA + 1; i < iB; ++i)
1707 {
1708 const V2d<T>& v = vertices[poly[i]];
1709 if(isInCircumcircle(v, a, b, *c))
1710 {
1711 out = i;
1712 c = &v;
1713 }
1714 }
1715 assert(out > iA && out < iB); // point is between ends
1716 return out;
1717}
1718
1719template <typename T, typename TNearPointLocator>
1721 const std::vector<V2d<T> >& newVertices)
1722{
1723 return insertVertices(
1724 newVertices.begin(), newVertices.end(), getX_V2d<T>, getY_V2d<T>);
1725}
1726
1727template <typename T, typename TNearPointLocator>
1729{
1730 return m_vertTris.empty() && !vertices.empty();
1731}
1732
1733template <typename T, typename TNearPointLocator>
1734unordered_map<TriInd, LayerDepth>
1736 std::stack<TriInd> seeds,
1737 const LayerDepth layerDepth,
1738 std::vector<LayerDepth>& triDepths) const
1739{
1740 unordered_map<TriInd, LayerDepth> behindBoundary;
1741 while(!seeds.empty())
1742 {
1743 const TriInd iT = seeds.top();
1744 seeds.pop();
1745 triDepths[iT] = std::min(triDepths[iT], layerDepth);
1746 behindBoundary.erase(iT);
1747 const Triangle& t = triangles[iT];
1748 for(Index i(0); i < Index(3); ++i)
1749 {
1750 const Edge opEdge(t.vertices[ccw(i)], t.vertices[cw(i)]);
1751 const TriInd iN = t.neighbors[opoNbr(i)];
1752 if(iN == noNeighbor || triDepths[iN] <= layerDepth)
1753 continue;
1754 if(fixedEdges.count(opEdge))
1755 {
1756 const unordered_map<Edge, LayerDepth>::const_iterator cit =
1757 overlapCount.find(opEdge);
1758 const LayerDepth triDepth = cit == overlapCount.end()
1759 ? layerDepth + 1
1760 : layerDepth + cit->second + 1;
1761 behindBoundary[iN] = triDepth;
1762 continue;
1763 }
1764 seeds.push(iN);
1765 }
1766 }
1767 return behindBoundary;
1768}
1769
1770template <typename T, typename TNearPointLocator>
1771std::vector<LayerDepth>
1773{
1774 std::vector<LayerDepth> triDepths(
1775 triangles.size(), std::numeric_limits<LayerDepth>::max());
1776 std::stack<TriInd> seeds(TriDeque(1, m_vertTris[0]));
1777 LayerDepth layerDepth = 0;
1778 LayerDepth deepestSeedDepth = 0;
1779
1780 unordered_map<LayerDepth, TriIndUSet> seedsByDepth;
1781 do
1782 {
1783 const unordered_map<TriInd, LayerDepth>& newSeeds =
1784 peelLayer(seeds, layerDepth, triDepths);
1785
1786 seedsByDepth.erase(layerDepth);
1787 typedef unordered_map<TriInd, LayerDepth>::const_iterator Iter;
1788 for(Iter it = newSeeds.begin(); it != newSeeds.end(); ++it)
1789 {
1790 deepestSeedDepth = std::max(deepestSeedDepth, it->second);
1791 seedsByDepth[it->second].insert(it->first);
1792 }
1793 const TriIndUSet& nextLayerSeeds = seedsByDepth[layerDepth + 1];
1794 seeds = std::stack<TriInd>(
1795 TriDeque(nextLayerSeeds.begin(), nextLayerSeeds.end()));
1796 ++layerDepth;
1797 } while(!seeds.empty() || deepestSeedDepth > layerDepth);
1798
1799 return triDepths;
1800}
1801
1802template <typename T, typename TNearPointLocator>
1804 VertInd superGeomVertCount)
1805{
1806 for(VertInd iV = superGeomVertCount; iV < vertices.size(); ++iV)
1807 {
1808 insertVertex(iV);
1809 }
1810}
1811
1812template <typename T, typename TNearPointLocator>
1813void Triangulation<T, TNearPointLocator>::insertVertices_Randomized(
1814 VertInd superGeomVertCount)
1815{
1816 std::size_t vertexCount = vertices.size() - superGeomVertCount;
1817 std::vector<VertInd> ii(vertexCount);
1818 detail::iota(ii.begin(), ii.end(), superGeomVertCount);
1819 detail::random_shuffle(ii.begin(), ii.end());
1820 for(std::vector<VertInd>::iterator it = ii.begin(); it != ii.end(); ++it)
1821 {
1822 insertVertex(*it);
1823 }
1824}
1825
1826namespace detail
1827{
1828
1829// log2 implementation backwards compatible with pre c++11
1830template <typename T>
1831inline double log2_bc(T x)
1832{
1833#ifdef CDT_CXX11_IS_SUPPORTED
1834 return std::log2(x);
1835#else
1836 static double log2_constant = std::log(2.0);
1837 return std::log(static_cast<double>(x)) / log2_constant;
1838#endif
1839}
1840
1845inline std::size_t maxQueueLengthBFSKDTree(const std::size_t vertexCount)
1846{
1847 const int filledLayerPow2 =
1848 static_cast<int>(std::floor(log2_bc(vertexCount)) - 1);
1849 const std::size_t nodesInFilledTree =
1850 static_cast<std::size_t>(std::pow(2., filledLayerPow2 + 1) - 1);
1851 const std::size_t nodesInLastFilledLayer =
1852 static_cast<std::size_t>(std::pow(2., filledLayerPow2));
1853 const std::size_t nodesInLastLayer = vertexCount - nodesInFilledTree;
1854 return nodesInLastLayer >= nodesInLastFilledLayer
1855 ? nodesInLastFilledLayer + nodesInLastLayer -
1856 nodesInLastFilledLayer
1857 : nodesInLastFilledLayer;
1858}
1859
1860template <typename T>
1862{
1863public:
1864 FixedCapacityQueue(const std::size_t capacity)
1865 : m_vec(capacity)
1866 , m_front(m_vec.begin())
1867 , m_back(m_vec.begin())
1868 , m_size(0)
1869 {}
1870 bool empty() const
1871 {
1872 return m_size == 0;
1873 }
1874 const T& front() const
1875 {
1876 return *m_front;
1877 }
1878 void pop()
1879 {
1880 assert(m_size > 0);
1881 ++m_front;
1882 if(m_front == m_vec.end())
1883 m_front = m_vec.begin();
1884 --m_size;
1885 }
1886 void push(const T& t)
1887 {
1888 assert(m_size < m_vec.size());
1889 *m_back = t;
1890 ++m_back;
1891 if(m_back == m_vec.end())
1892 m_back = m_vec.begin();
1893 ++m_size;
1894 }
1895#ifdef CDT_CXX11_IS_SUPPORTED
1896 void push(const T&& t)
1897 {
1898 assert(m_size < m_vec.size());
1899 *m_back = t;
1900 ++m_back;
1901 if(m_back == m_vec.end())
1902 m_back = m_vec.begin();
1903 ++m_size;
1904 }
1905#endif
1906private:
1907 std::vector<T> m_vec;
1908 typename std::vector<T>::iterator m_front;
1909 typename std::vector<T>::iterator m_back;
1910 std::size_t m_size;
1911};
1912
1913template <typename T>
1915{
1916 const std::vector<V2d<T> >& m_vertices;
1917
1918public:
1919 less_than_x(const std::vector<V2d<T> >& vertices)
1920 : m_vertices(vertices)
1921 {}
1922 bool operator()(const VertInd a, const VertInd b) const
1923 {
1924 return m_vertices[a].x < m_vertices[b].x;
1925 }
1926};
1927
1928template <typename T>
1930{
1931 const std::vector<V2d<T> >& m_vertices;
1932
1933public:
1934 less_than_y(const std::vector<V2d<T> >& vertices)
1935 : m_vertices(vertices)
1936 {}
1937 bool operator()(const VertInd a, const VertInd b) const
1938 {
1939 return m_vertices[a].y < m_vertices[b].y;
1940 }
1941};
1942
1943} // namespace detail
1944
1945template <typename T, typename TNearPointLocator>
1947 VertInd superGeomVertCount,
1948 V2d<T> boxMin,
1949 V2d<T> boxMax)
1950{
1951 // calculate original indices
1952 const VertInd vertexCount =
1953 static_cast<IndexSizeType>(vertices.size()) - superGeomVertCount;
1954 if(vertexCount <= 0)
1955 return;
1956 std::vector<VertInd> ii(vertexCount);
1957 detail::iota(ii.begin(), ii.end(), superGeomVertCount);
1958
1959 typedef std::vector<VertInd>::iterator It;
1961 detail::maxQueueLengthBFSKDTree(vertexCount));
1962 queue.push(make_tuple(ii.begin(), ii.end(), boxMin, boxMax, VertInd(0)));
1963
1964 It first, last;
1965 V2d<T> newBoxMin, newBoxMax;
1966 VertInd parent, mid;
1967
1968 const detail::less_than_x<T> cmpX(vertices);
1969 const detail::less_than_y<T> cmpY(vertices);
1970
1971 while(!queue.empty())
1972 {
1973 tie(first, last, boxMin, boxMax, parent) = queue.front();
1974 queue.pop();
1975 assert(first != last);
1976
1977 const std::ptrdiff_t len = std::distance(first, last);
1978 if(len == 1)
1979 {
1980 insertVertex(*first, parent);
1981 continue;
1982 }
1983 const It midIt = first + len / 2;
1984 if(boxMax.x - boxMin.x >= boxMax.y - boxMin.y)
1985 {
1986 detail::portable_nth_element(first, midIt, last, cmpX);
1987 mid = *midIt;
1988 const T split = vertices[mid].x;
1989 newBoxMin.x = split;
1990 newBoxMin.y = boxMin.y;
1991 newBoxMax.x = split;
1992 newBoxMax.y = boxMax.y;
1993 }
1994 else
1995 {
1996 detail::portable_nth_element(first, midIt, last, cmpY);
1997 mid = *midIt;
1998 const T split = vertices[mid].y;
1999 newBoxMin.x = boxMin.x;
2000 newBoxMin.y = split;
2001 newBoxMax.x = boxMax.x;
2002 newBoxMax.y = split;
2003 }
2004 insertVertex(mid, parent);
2005 if(first != midIt)
2006 {
2007 queue.push(make_tuple(first, midIt, boxMin, newBoxMax, mid));
2008 }
2009 if(midIt + 1 != last)
2010 {
2011 queue.push(make_tuple(midIt + 1, last, newBoxMin, boxMax, mid));
2012 }
2013 }
2014}
2015
2016template <typename T, typename TNearPointLocator>
2017std::pair<TriInd, TriInd> Triangulation<T, TNearPointLocator>::edgeTriangles(
2018 const VertInd a,
2019 const VertInd b) const
2020{
2021 const TriInd triStart = m_vertTris[a];
2022 assert(triStart != noNeighbor);
2023 TriInd iT = triStart, iTNext = triStart;
2024 VertInd iV = noVertex;
2025 do
2026 {
2027 const Triangle& t = triangles[iT];
2028 tie(iTNext, iV) = t.next(a);
2029 assert(iTNext != noNeighbor);
2030 if(iV == b)
2031 {
2032 return std::make_pair(iT, iTNext);
2033 }
2034 iT = iTNext;
2035 } while(iT != triStart);
2036 return std::make_pair(invalidIndex, invalidIndex);
2037}
2038
2039template <typename T, typename TNearPointLocator>
2040bool Triangulation<T, TNearPointLocator>::hasEdge(
2041 const VertInd a,
2042 const VertInd b) const
2043{
2044 return edgeTriangles(a, b).first != invalidIndex;
2045}
2046
2047template <typename T, typename TNearPointLocator>
2048void Triangulation<T, TNearPointLocator>::setAdjacentTriangle(
2049 const VertInd v,
2050 const TriInd t)
2051{
2052 assert(t != noNeighbor);
2053 m_vertTris[v] = t;
2054 assert(
2055 triangles[t].vertices[0] == v || triangles[t].vertices[1] == v ||
2056 triangles[t].vertices[2] == v);
2057}
2058
2059template <typename T, typename TNearPointLocator>
2060void Triangulation<T, TNearPointLocator>::pivotVertexTriangleCW(const VertInd v)
2061{
2062 assert(m_vertTris[v] != noNeighbor);
2063 m_vertTris[v] = triangles[m_vertTris[v]].next(v).first;
2064 assert(m_vertTris[v] != noNeighbor);
2065 assert(
2066 triangles[m_vertTris[v]].vertices[0] == v ||
2067 triangles[m_vertTris[v]].vertices[1] == v ||
2068 triangles[m_vertTris[v]].vertices[2] == v);
2069}
2070
2071template <typename T, typename TNearPointLocator>
2072void Triangulation<T, TNearPointLocator>::tryAddVertexToLocator(const VertInd v)
2073{
2074 if(!m_nearPtLocator.empty()) // only if locator is initialized already
2075 m_nearPtLocator.addPoint(v, vertices);
2076}
2077
2078template <typename T, typename TNearPointLocator>
2079void Triangulation<T, TNearPointLocator>::tryInitNearestPointLocator()
2080{
2081 if(!vertices.empty() && m_nearPtLocator.empty())
2082 {
2083 m_nearPtLocator.initialize(vertices);
2084 }
2085}
2086
2087} // namespace CDT
Triangulation class.
array< T, 3 > arr3(const T &v0, const T &v1, const T &v2)
Needed for c++03 compatibility (no uniform initialization available)
std::size_t maxQueueLengthBFSKDTree(const std::size_t vertexCount)
Since KD-tree bulk load builds a balanced tree the maximum length of a queue can be pre-calculated: i...
Data structure representing a 2D constrained Delaunay triangulation.
void eraseOuterTriangles()
Erase triangles outside of constrained boundary using growing.
void conformToEdges(TEdgeIter first, TEdgeIter last, TGetEdgeVertexStart getStart, TGetEdgeVertexEnd getEnd)
Insert constraint edges into triangulation for Conforming Delaunay Triangulation (for example see fig...
std::vector< LayerDepth > calculateTriangleDepths() const
Calculate depth of each triangle in constraint triangulation.
bool isFinalized() const
Check if the triangulation was finalized with erase... method and super-triangle was removed.
void eraseOuterTrianglesAndHoles()
Erase triangles outside of constrained boundary and auto-detected holes.
void insertEdges(TEdgeIter first, TEdgeIter last, TGetEdgeVertexStart getStart, TGetEdgeVertexEnd getEnd)
Insert constraint edges into triangulation for Constrained Delaunay Triangulation (for example see fi...
void insertVertices(TVertexIter first, TVertexIter last, TGetVertexCoordX getX, TGetVertexCoordY getY)
Insert custom point-types specified by iterator range and X/Y-getters.
void eraseSuperTriangle()
Erase triangles adjacent to super triangle.
Triangulation()
Default constructor.
void initializedWithCustomSuperGeometry()
Call this method after directly setting custom super-geometry via vertices and triangles members.
unsigned short LayerDepth
Type used for storing layer depths for triangles.
Definition CDT.h:40
#define CDT_SOURCE_LOCATION
Macro for getting source location.
void flipEdge(TriInd iT, TriInd iTopo)
Flip an edge between two triangle.
TriIndVec & VertTrisInternal()
Access internal vertex adjacent triangles.
void removeTriangles(const TriIndUSet &removedTriangles)
Remove triangles with specified indices.
Namespace containing triangulation functionality.
std::vector< Edge > EdgeVec
Vector of edges.
Definition CDTUtils.h:245
CDT_EXPORT CDT_INLINE_IF_HEADER_ONLY VertInd opposedVertex(const Triangle &tri, TriInd iTopo)
Given two triangles, return vertex of first triangle opposed to the second.
Definition CDTUtils.hpp:255
unordered_set< Edge > EdgeUSet
Hash table of edges.
Definition CDTUtils.h:246
VertInd edge_get_v2(const Edge &e)
Get edge second vertex.
Definition CDTUtils.h:234
std::vector< TriInd > TriIndVec
Vector of triangle indices.
Definition CDTUtils.h:155
CDT_INLINE_IF_HEADER_ONLY Index edgeNeighborInd(const VerticesArr3 &vv, VertInd iVedge1, VertInd iVedge2)
Index of triangle's neighbor opposed to an edge.
Definition CDTUtils.hpp:193
VertInd edge_get_v1(const Edge &e)
Get edge first vertex.
Definition CDTUtils.h:228
CDT_EXPORT PtLineLocation::Enum classifyOrientation(T orientation, T orientationTolerance=T(0))
Classify value of orient2d predicate.
Definition CDTUtils.hpp:116
array< TriInd, 3 > NeighborsArr3
array of three neighbors
Definition CDTUtils.h:157
CDT_EXPORT T distance(const V2d< T > &a, const V2d< T > &b)
Distance between two 2D points.
Definition CDTUtils.hpp:302
IndexSizeType VertInd
Vertex index.
Definition CDTUtils.h:142
CDT_EXPORT Index cw(Index i)
Advance vertex or neighbor index clockwise.
Definition CDTUtils.hpp:80
array< VertInd, 3 > VerticesArr3
array of three vertex indices
Definition CDTUtils.h:156
CDT_INLINE_IF_HEADER_ONLY bool touchesSuperTriangle(const Triangle &t)
Check if any of triangle's vertices belongs to a super-triangle.
Definition CDTUtils.hpp:313
CDT_EXPORT CDT_INLINE_IF_HEADER_ONLY Index opoNbr(Index vertIndex)
Opposed neighbor index from vertex index.
Definition CDTUtils.hpp:158
CDT_EXPORT CDT_INLINE_IF_HEADER_ONLY Index vertexInd(const VerticesArr3 &vv, VertInd iV)
If triangle has a given vertex return vertex-index.
Definition CDTUtils.hpp:238
CDT_EXPORT T orient2D(const V2d< T > &p, const V2d< T > &v1, const V2d< T > &v2)
Orient p against line v1-v2 2D: robust geometric predicate.
Definition CDTUtils.hpp:99
CDT_EXPORT PtLineLocation::Enum locatePointLine(const V2d< T > &p, const V2d< T > &v1, const V2d< T > &v2, T orientationTolerance=T(0))
Check if point lies to the left of, to the right of, or on a line.
Definition CDTUtils.hpp:105
unordered_set< TriInd > TriIndUSet
Hash table of triangles.
Definition CDTUtils.h:247
CDT_EXPORT Index ccw(Index i)
Advance vertex or neighbor index counter-clockwise.
Definition CDTUtils.hpp:75
CDT_EXPORT Index edgeNeighbor(PtTriLocation::Enum location)
Neighbor index from a on-edge location.
Definition CDTUtils.hpp:92
unordered_map< TriInd, TriInd > TriIndUMap
Triangle hash map.
Definition CDTUtils.h:248
const T & getX_V2d(const V2d< T > &v)
X- coordinate getter for V2d.
Definition CDTUtils.h:106
unsigned char Index
Index in triangle.
Definition CDTUtils.h:140
CDT_EXPORT CDT_INLINE_IF_HEADER_ONLY TriInd opposedTriangle(const Triangle &tri, VertInd iVert)
Given triangle and a vertex find opposed triangle.
Definition CDTUtils.hpp:249
CDT_EXPORT bool isInCircumcircle(const V2d< T > &p, const V2d< T > &v1, const V2d< T > &v2, const V2d< T > &v3)
Test if point lies in a circumscribed circle of a triangle.
Definition CDTUtils.hpp:268
Edge RemapNoSuperTriangle(const Edge &e)
Remap removing super-triangle: subtract 3 from vertices.
CDT_EXPORT bool isOnEdge(PtTriLocation::Enum location)
Check if location is classified as on any of three edges.
Definition CDTUtils.hpp:85
IndexSizeType TriInd
Triangle index.
Definition CDTUtils.h:144
CDT_EXPORT PtTriLocation::Enum locatePointTriangle(const V2d< T > &p, const V2d< T > &v1, const V2d< T > &v2, const V2d< T > &v3)
Check if point a lies inside of, outside of, or on an edge of a triangle.
Definition CDTUtils.hpp:126
CDT_EXPORT CDT_INLINE_IF_HEADER_ONLY Index opposedVertexInd(const NeighborsArr3 &nn, TriInd iTopo)
Index of triangle's vertex opposed to a triangle.
Definition CDTUtils.hpp:227
const T & getY_V2d(const V2d< T > &v)
Y-coordinate getter for V2d.
Definition CDTUtils.h:113
Edge connecting two vertices: vertex with smaller index is always first.
Definition CDTUtils.h:209
VertInd v2() const
V2 getter.
Definition CDTUtils.hpp:62
VertInd v1() const
V1 getter.
Definition CDTUtils.hpp:57
@ TryResolve
attempt to resolve constraint edge intersections
@ NotAllowed
constraint edge intersections are not allowed
@ DontCheck
No checks: slightly faster but less safe.
Enum
The Enum itself.
@ SuperTriangle
conventional super-triangle
@ Custom
user-specified custom geometry (e.g., grid)
Triangulation triangle (counter-clockwise winding)
Definition CDTUtils.h:259
VerticesArr3 vertices
triangle's three vertices
Definition CDTUtils.h:260
static Triangle make(const array< VertInd, 3 > &vertices, const array< TriInd, 3 > &neighbors)
Factory method.
Definition CDTUtils.h:269
NeighborsArr3 neighbors
triangle's three neighbors
Definition CDTUtils.h:261
std::pair< TriInd, VertInd > next(const VertInd i) const
Next triangle adjacent to a vertex (clockwise)
Definition CDTUtils.h:277
2D vector
Definition CDTUtils.h:96
T y
Y-coordinate.
Definition CDTUtils.h:98
static V2d make(T x, T y)
Create vector from X and Y coordinates.
Definition CDTUtils.hpp:23
T x
X-coordinate.
Definition CDTUtils.h:97
@ Auto
Automatic insertion order optimized for better performance.