CDT  1.3.0
C++ library for constrained Delaunay triangulation
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 Triangle& t = triangles[iT];
140 if(t.vertices[0] < 3 || t.vertices[1] < 3 || t.vertices[2] < 3)
141 toErase.insert(iT);
142 }
143 finalizeTriangulation(toErase);
144}
145
146template <typename T, typename TNearPointLocator>
148{
149 // make dummy triangles adjacent to super-triangle's vertices
150 assert(m_vertTris[0] != noNeighbor);
151 const std::stack<TriInd> seed(std::deque<TriInd>(1, m_vertTris[0]));
152 const TriIndUSet toErase = growToBoundary(seed);
153 finalizeTriangulation(toErase);
154}
155
156template <typename T, typename TNearPointLocator>
158{
159 const std::vector<LayerDepth> triDepths = calculateTriangleDepths();
160 TriIndUSet toErase;
161 toErase.reserve(triangles.size());
162 for(std::size_t iT = 0; iT != triangles.size(); ++iT)
163 {
164 if(triDepths[iT] % 2 == 0)
165 toErase.insert(static_cast<TriInd>(iT));
166 }
167 finalizeTriangulation(toErase);
168}
169
172{
173 return Edge(VertInd(e.v1() - 3), VertInd(e.v2() - 3));
174}
175
176template <typename T, typename TNearPointLocator>
178 const TriIndUSet& removedTriangles)
179{
180 if(removedTriangles.empty())
181 return;
182 // remove triangles and calculate triangle index mapping
183 TriIndUMap triIndMap;
184 for(TriInd iT(0), iTnew(0); iT < TriInd(triangles.size()); ++iT)
185 {
186 if(removedTriangles.count(iT))
187 continue;
188 triIndMap[iT] = iTnew;
189 triangles[iTnew] = triangles[iT];
190 iTnew++;
191 }
192 triangles.erase(triangles.end() - removedTriangles.size(), triangles.end());
193 // adjust triangles' neighbors
194 for(TriInd iT(0); iT < triangles.size(); ++iT)
195 {
196 Triangle& t = triangles[iT];
197 // update neighbors to account for removed triangles
198 NeighborsArr3& nn = t.neighbors;
199 for(NeighborsArr3::iterator n = nn.begin(); n != nn.end(); ++n)
200 {
201 if(removedTriangles.count(*n))
202 {
203 *n = noNeighbor;
204 }
205 else if(*n != noNeighbor)
206 {
207 *n = triIndMap[*n];
208 }
209 }
210 }
211}
212template <typename T, typename TNearPointLocator>
214{
215 return m_vertTris;
216}
217
218template <typename T, typename TNearPointLocator>
220 const TriIndUSet& removedTriangles)
221{
222 eraseDummies();
223 m_vertTris = TriIndVec();
224 // remove super-triangle
225 if(m_superGeomType == SuperGeometryType::SuperTriangle)
226 {
227 vertices.erase(vertices.begin(), vertices.begin() + 3);
228 // Edge re-mapping
229 { // fixed edges
230 EdgeUSet updatedFixedEdges;
231 typedef CDT::EdgeUSet::const_iterator It;
232 for(It e = fixedEdges.begin(); e != fixedEdges.end(); ++e)
233 {
234 updatedFixedEdges.insert(RemapNoSuperTriangle(*e));
235 }
236 fixedEdges = updatedFixedEdges;
237 }
238 { // overlap count
239 unordered_map<Edge, BoundaryOverlapCount> updatedOverlapCount;
240 typedef unordered_map<Edge, BoundaryOverlapCount>::const_iterator
241 It;
242 for(It it = overlapCount.begin(); it != overlapCount.end(); ++it)
243 {
244 updatedOverlapCount.insert(std::make_pair(
245 RemapNoSuperTriangle(it->first), it->second));
246 }
247 overlapCount = updatedOverlapCount;
248 }
249 { // split edges mapping
250 unordered_map<Edge, EdgeVec> updatedPieceToOriginals;
251 typedef unordered_map<Edge, EdgeVec>::const_iterator It;
252 for(It it = pieceToOriginals.begin(); it != pieceToOriginals.end();
253 ++it)
254 {
255 EdgeVec ee = it->second;
256 for(EdgeVec::iterator eeIt = ee.begin(); eeIt != ee.end();
257 ++eeIt)
258 {
259 *eeIt = RemapNoSuperTriangle(*eeIt);
260 }
261 updatedPieceToOriginals.insert(
262 std::make_pair(RemapNoSuperTriangle(it->first), ee));
263 }
264 pieceToOriginals = updatedPieceToOriginals;
265 }
266 }
267 // remove other triangles
268 removeTriangles(removedTriangles);
269 // adjust triangle vertices: account for removed super-triangle
270 if(m_superGeomType == SuperGeometryType::SuperTriangle)
271 {
272 for(TriangleVec::iterator t = triangles.begin(); t != triangles.end();
273 ++t)
274 {
275 VerticesArr3& vv = t->vertices;
276 for(VerticesArr3::iterator v = vv.begin(); v != vv.end(); ++v)
277 {
278 *v -= 3;
279 }
280 }
281 }
282}
283
284template <typename T, typename TNearPointLocator>
286{
287 m_nearPtLocator.initialize(vertices);
288 m_nTargetVerts = vertices.size();
289 m_superGeomType = SuperGeometryType::Custom;
290}
291
292template <typename T, typename TNearPointLocator>
294 std::stack<TriInd> seeds) const
295{
296 TriIndUSet traversed;
297 while(!seeds.empty())
298 {
299 const TriInd iT = seeds.top();
300 seeds.pop();
301 traversed.insert(iT);
302 const Triangle& t = triangles[iT];
303 for(Index i(0); i < Index(3); ++i)
304 {
305 const Edge opEdge(t.vertices[ccw(i)], t.vertices[cw(i)]);
306 if(fixedEdges.count(opEdge))
307 continue;
308 const TriInd iN = t.neighbors[opoNbr(i)];
309 if(iN != noNeighbor && traversed.count(iN) == 0)
310 seeds.push(iN);
311 }
312 }
313 return traversed;
314}
315
316template <typename T, typename TNearPointLocator>
317void Triangulation<T, TNearPointLocator>::makeDummy(const TriInd iT)
318{
319 m_dummyTris.push_back(iT);
320}
321
322template <typename T, typename TNearPointLocator>
323TriInd Triangulation<T, TNearPointLocator>::addTriangle(const Triangle& t)
324{
325 if(m_dummyTris.empty())
326 {
327 triangles.push_back(t);
328 return TriInd(triangles.size() - 1);
329 }
330 const TriInd nxtDummy = m_dummyTris.back();
331 m_dummyTris.pop_back();
332 triangles[nxtDummy] = t;
333 return nxtDummy;
334}
335
336template <typename T, typename TNearPointLocator>
337TriInd Triangulation<T, TNearPointLocator>::addTriangle()
338{
339 if(m_dummyTris.empty())
340 {
341 const Triangle dummy = {
342 {noVertex, noVertex, noVertex},
343 {noNeighbor, noNeighbor, noNeighbor}};
344 triangles.push_back(dummy);
345 return TriInd(triangles.size() - 1);
346 }
347 const TriInd nxtDummy = m_dummyTris.back();
348 m_dummyTris.pop_back();
349 return nxtDummy;
350}
351
352template <typename T, typename TNearPointLocator>
354 const std::vector<Edge>& edges)
355{
356 insertEdges(edges.begin(), edges.end(), edge_get_v1, edge_get_v2);
357}
358
359template <typename T, typename TNearPointLocator>
361 const std::vector<Edge>& edges)
362{
363 conformToEdges(edges.begin(), edges.end(), edge_get_v1, edge_get_v2);
364}
365
366template <typename T, typename TNearPointLocator>
368{
369 if(!fixedEdges.insert(edge).second)
370 {
371 ++overlapCount[edge]; // if edge is already fixed increment the counter
372 }
373}
374
375namespace detail
376{
377
378// add element to 'to' if not already in 'to'
379template <typename T, typename Allocator1>
380void insert_unique(std::vector<T, Allocator1>& to, const T& elem)
381{
382 if(std::find(to.begin(), to.end(), elem) == to.end())
383 {
384 to.push_back(elem);
385 }
386}
387
388// add elements of 'from' that are not present in 'to' to 'to'
389template <typename T, typename Allocator1, typename Allocator2>
390void insert_unique(
391 std::vector<T, Allocator1>& to,
392 const std::vector<T, Allocator2>& from)
393{
394 typedef typename std::vector<T, Allocator2>::const_iterator Cit;
395 to.reserve(to.size() + from.size());
396 for(Cit cit = from.begin(); cit != from.end(); ++cit)
397 {
398 insert_unique(to, *cit);
399 }
400}
401
402} // namespace detail
403
404template <typename T, typename TNearPointLocator>
405void Triangulation<T, TNearPointLocator>::fixEdge(
406 const Edge& edge,
407 const Edge& originalEdge)
408{
409 fixEdge(edge);
410 if(edge != originalEdge)
411 detail::insert_unique(pieceToOriginals[edge], originalEdge);
412}
413
414template <typename T, typename TNearPointLocator>
415void Triangulation<T, TNearPointLocator>::fixEdge(
416 const Edge& edge,
417 const BoundaryOverlapCount overlaps)
418{
419 fixedEdges.insert(edge);
420 overlapCount[edge] = overlaps; // override overlap counter
421}
422
423namespace detail
424{
425
426template <typename T>
427T lerp(const T& a, const T& b, const T t)
428{
429 return (T(1) - t) * a + t * b;
430}
431
432// Precondition: ab and cd intersect normally
433template <typename T>
434V2d<T> intersectionPosition(
435 const V2d<T>& a,
436 const V2d<T>& b,
437 const V2d<T>& c,
438 const V2d<T>& d)
439{
440 using namespace predicates::adaptive;
441 // interpolate point on the shorter segment
442 if(distanceSquared(a, b) < distanceSquared(c, d))
443 {
444 const T a_cd = orient2d(c.x, c.y, d.x, d.y, a.x, a.y);
445 const T b_cd = orient2d(c.x, c.y, d.x, d.y, b.x, b.y);
446 const T t = a_cd / (a_cd - b_cd);
447 return V2d<T>::make(lerp(a.x, b.x, t), lerp(a.y, b.y, t));
448 }
449 else
450 {
451 const T c_ab = orient2d(a.x, a.y, b.x, b.y, c.x, c.y);
452 const T d_ab = orient2d(a.x, a.y, b.x, b.y, d.x, d.y);
453 const T t = c_ab / (c_ab - d_ab);
454 return V2d<T>::make(lerp(c.x, d.x, t), lerp(c.y, d.y, t));
455 }
456}
457
458} // namespace detail
459
460template <typename T, typename TNearPointLocator>
461IndexSizeType Triangulation<T, TNearPointLocator>::previousEdgeEncounter(
462 const VertInd v,
463 const std::vector<VertInd>& poly) const
464{
465 // check if vertex was encountered earlier in the pseudo-polygons:
466 // it's triangle was reset
467 if(m_vertTris[v] != noNeighbor)
468 return invalidIndex;
469 // rewind back to the previous entry of the point
470 IndexSizeType i = poly.size() - 1;
471 while(poly[i] != v)
472 {
473 assert(i > 0);
474 --i;
475 }
476 // was the edge or only the vertex encountered before?
477 // check if previous 'next' is same as current 'previous'
478 if(poly[i + 1] == poly.back())
479 return i; // edge was encountered earlier
480 return invalidIndex; // edge was not encountered earlier
481}
482
483template <typename T, typename TNearPointLocator>
484void Triangulation<T, TNearPointLocator>::insertEdgeIteration(
485 const Edge edge,
486 const Edge originalEdge,
487 EdgeVec& remaining,
488 std::vector<TriangulatePseudopolygonTask>& tppIterations)
489{
490 const VertInd iA = edge.v1();
491 VertInd iB = edge.v2();
492 if(iA == iB) // edge connects a vertex to itself
493 return;
494
495 if(hasEdge(iA, iB))
496 {
497 fixEdge(edge, originalEdge);
498 return;
499 }
500
501 const V2d<T>& a = vertices[iA];
502 const V2d<T>& b = vertices[iB];
503 const T distanceTolerance =
504 m_minDistToConstraintEdge == T(0)
505 ? T(0)
506 : m_minDistToConstraintEdge * distance(a, b);
507
508 TriInd iT;
509 // Note: 'L' is left and 'R' is right of the inserted constraint edge
510 VertInd iVL, iVR;
511 tie(iT, iVL, iVR) = intersectedTriangle(iA, a, b, distanceTolerance);
512 // if one of the triangle vertices is on the edge, move edge start
513 if(iT == noNeighbor)
514 {
515 const Edge edgePart(iA, iVL);
516 fixEdge(edgePart, originalEdge);
517 remaining.push_back(Edge(iVL, iB));
518 return;
519 }
520 Triangle t = triangles[iT];
521 std::vector<TriInd> intersected(1, iT);
522 std::vector<VertInd> polyL, polyR;
523 std::vector<TriInd> outerTrisL, outerTrisR;
524 polyL.reserve(2);
525 polyL.push_back(iA);
526 polyL.push_back(iVL);
527 outerTrisL.push_back(edgeNeighbor(t, iA, iVL));
528 polyR.reserve(2);
529 polyR.push_back(iA);
530 polyR.push_back(iVR);
531 outerTrisR.push_back(edgeNeighbor(t, iA, iVR));
532 removeAdjacentTriangle(iVL);
533 removeAdjacentTriangle(iVR);
534 VertInd iV = iA;
535
536 while(!t.containsVertex(iB))
537 {
538 const TriInd iTopo = opposedTriangle(t, iV);
539 const Triangle& tOpo = triangles[iTopo];
540 const VertInd iVopo = opposedVertex(tOpo, iT);
541
542 // Resolve intersection between two constraint edges if needed
543 if(m_intersectingEdgesStrategy ==
545 fixedEdges.count(Edge(iVL, iVR)))
546 {
547 const VertInd iNewVert = static_cast<VertInd>(vertices.size());
548
549 // split constraint edge that already exists in triangulation
550 const Edge splitEdge(iVL, iVR);
551 const Edge half1(iVL, iNewVert);
552 const Edge half2(iNewVert, iVR);
553 const BoundaryOverlapCount overlaps = overlapCount[splitEdge];
554 // remove the edge that will be split
555 fixedEdges.erase(splitEdge);
556 overlapCount.erase(splitEdge);
557 // add split edge's halves
558 fixEdge(half1, overlaps);
559 fixEdge(half2, overlaps);
560 // maintain piece-to-original mapping
561 EdgeVec newOriginals(1, splitEdge);
562 const unordered_map<Edge, EdgeVec>::const_iterator originalsIt =
563 pieceToOriginals.find(splitEdge);
564 if(originalsIt != pieceToOriginals.end())
565 { // edge being split was split before: pass-through originals
566 newOriginals = originalsIt->second;
567 pieceToOriginals.erase(originalsIt);
568 }
569 detail::insert_unique(pieceToOriginals[half1], newOriginals);
570 detail::insert_unique(pieceToOriginals[half2], newOriginals);
571 // set adjacent triangles again
572 for(IndexSizeType i = 0; i < outerTrisL.size(); ++i)
573 setAdjacentTriangle(polyL[i + 1], outerTrisL[i]);
574 for(IndexSizeType i = 0; i < outerTrisR.size(); ++i)
575 setAdjacentTriangle(polyR[i + 1], outerTrisR[i]);
576 // add a new point at the intersection of two constraint edges
577 const V2d<T> newV = detail::intersectionPosition(
578 vertices[iA], vertices[iB], vertices[iVL], vertices[iVR]);
579 addNewVertex(newV, noNeighbor);
580 std::stack<TriInd> triStack =
581 insertVertexOnEdge(iNewVert, iT, iTopo);
582 tryAddVertexToLocator(iNewVert);
583 ensureDelaunayByEdgeFlips(newV, iNewVert, triStack);
584 // TODO: is it's possible to re-use pseudo-polygons
585 // for inserting [iA, iNewVert] edge half?
586 remaining.push_back(Edge(iA, iNewVert));
587 remaining.push_back(Edge(iNewVert, iB));
588 return;
589 }
590
591 const PtLineLocation::Enum loc =
592 locatePointLine(vertices[iVopo], a, b, distanceTolerance);
593 if(loc == PtLineLocation::Left)
594 {
595 // hanging edge check
596 const IndexSizeType prev = previousEdgeEncounter(iVopo, polyL);
597 if(prev != invalidIndex)
598 {
599 outerTrisL[prev] = noNeighbor; // previous entry of hanging edge
600 outerTrisL.push_back(noNeighbor);
601 }
602 else
603 {
604 outerTrisL.push_back(edgeNeighbor(tOpo, polyL.back(), iVopo));
605 }
606 polyL.push_back(iVopo);
607 removeAdjacentTriangle(iVopo);
608 iV = iVL;
609 iVL = iVopo;
610 }
611 else if(loc == PtLineLocation::Right)
612 {
613 // hanging edge check
614 const IndexSizeType prev = previousEdgeEncounter(iVopo, polyR);
615 if(prev != invalidIndex)
616 {
617 outerTrisR[prev] = noNeighbor; // previous entry of hanging edge
618 outerTrisR.push_back(noNeighbor);
619 }
620 else
621 {
622 outerTrisR.push_back(edgeNeighbor(tOpo, polyR.back(), iVopo));
623 }
624 polyR.push_back(iVopo);
625 removeAdjacentTriangle(iVopo);
626 iV = iVR;
627 iVR = iVopo;
628 }
629 else // encountered point on the edge
630 iB = iVopo;
631
632 intersected.push_back(iTopo);
633 iT = iTopo;
634 t = triangles[iT];
635 }
636 outerTrisL.push_back(edgeNeighbor(t, polyL.back(), iB));
637 outerTrisR.push_back(edgeNeighbor(t, polyR.back(), iB));
638 polyL.push_back(iB);
639 polyR.push_back(iB);
640
641 assert(!intersected.empty());
642 // make sure start/end vertices have a valid adjacent triangle
643 // that is not intersected by an edge
644 if(m_vertTris[iA] == intersected.front())
645 pivotVertexTriangleCW(iA);
646 if(m_vertTris[iB] == intersected.back())
647 pivotVertexTriangleCW(iB);
648 // reset adjacent triangle for pseudo-polygons' vertices
649 // typedef std::vector<VertInd>::const_iterator PolyCit;
650 // for(PolyCit it = polyL.begin() + 1; it != polyL.end() - 1; ++it)
651 // removeAdjacentTriangle(*it);
652 // for(PolyCit it = polyR.begin() + 1; it != polyR.end() - 1; ++it)
653 // removeAdjacentTriangle(*it);
654 // Remove intersected triangles
655 typedef std::vector<TriInd>::const_iterator TriIndCit;
656 for(TriIndCit it = intersected.begin(); it != intersected.end(); ++it)
657 makeDummy(*it);
658 // Triangulate pseudo-polygons on both sides
659 std::reverse(polyR.begin(), polyR.end());
660 std::reverse(outerTrisR.begin(), outerTrisR.end());
661 const TriInd iTL = addTriangle();
662 const TriInd iTR = addTriangle();
663 triangulatePseudopolygon(polyL, outerTrisL, iTL, iTR, tppIterations);
664 triangulatePseudopolygon(polyR, outerTrisR, iTR, iTL, tppIterations);
665
666 if(iB != edge.v2()) // encountered point on the edge
667 {
668 // fix edge part
669 const Edge edgePart(iA, iB);
670 fixEdge(edgePart, originalEdge);
671 remaining.push_back(Edge(iB, edge.v2()));
672 return;
673 }
674 else
675 {
676 fixEdge(edge, originalEdge);
677 }
678}
679
680template <typename T, typename TNearPointLocator>
681void Triangulation<T, TNearPointLocator>::insertEdge(
682 Edge edge,
683 const Edge originalEdge,
684 EdgeVec& remaining,
685 std::vector<TriangulatePseudopolygonTask>& tppIterations)
686{
687 // use iteration over recursion to avoid stack overflows
688 remaining.clear();
689 remaining.push_back(edge);
690 while(!remaining.empty())
691 {
692 edge = remaining.back();
693 remaining.pop_back();
694 insertEdgeIteration(edge, originalEdge, remaining, tppIterations);
695 }
696}
697
698template <typename T, typename TNearPointLocator>
699void Triangulation<T, TNearPointLocator>::conformToEdgeIteration(
700 Edge edge,
701 const EdgeVec& originals,
702 BoundaryOverlapCount overlaps,
703 std::vector<ConformToEdgeTask>& remaining)
704{
705 const VertInd iA = edge.v1();
706 VertInd iB = edge.v2();
707 if(iA == iB) // edge connects a vertex to itself
708 return;
709
710 if(hasEdge(iA, iB))
711 {
712 overlaps > 0 ? fixEdge(edge, overlaps) : fixEdge(edge);
713 // avoid marking edge as a part of itself
714 if(!originals.empty() && edge != originals.front())
715 {
716 detail::insert_unique(pieceToOriginals[edge], originals);
717 }
718 return;
719 }
720
721 const V2d<T>& a = vertices[iA];
722 const V2d<T>& b = vertices[iB];
723 const T distanceTolerance =
724 m_minDistToConstraintEdge == T(0)
725 ? T(0)
726 : m_minDistToConstraintEdge * distance(a, b);
727 TriInd iT;
728 VertInd iVleft, iVright;
729 tie(iT, iVleft, iVright) = intersectedTriangle(iA, a, b, distanceTolerance);
730 // if one of the triangle vertices is on the edge, move edge start
731 if(iT == noNeighbor)
732 {
733 const Edge edgePart(iA, iVleft);
734 overlaps > 0 ? fixEdge(edgePart, overlaps) : fixEdge(edgePart);
735 detail::insert_unique(pieceToOriginals[edgePart], originals);
736#ifdef CDT_CXX11_IS_SUPPORTED
737 remaining.emplace_back(Edge(iVleft, iB), originals, overlaps);
738#else
739 remaining.push_back(make_tuple(Edge(iVleft, iB), originals, overlaps));
740#endif
741 return;
742 }
743
744 VertInd iV = iA;
745 Triangle t = triangles[iT];
746 while(std::find(t.vertices.begin(), t.vertices.end(), iB) ==
747 t.vertices.end())
748 {
749 const TriInd iTopo = opposedTriangle(t, iV);
750 const Triangle& tOpo = triangles[iTopo];
751 const VertInd iVopo = opposedVertex(tOpo, iT);
752 const V2d<T> vOpo = vertices[iVopo];
753
754 // Resolve intersection between two constraint edges if needed
755 if(m_intersectingEdgesStrategy ==
757 fixedEdges.count(Edge(iVleft, iVright)))
758 {
759 const VertInd iNewVert = static_cast<VertInd>(vertices.size());
760
761 // split constraint edge that already exists in triangulation
762 const Edge splitEdge(iVleft, iVright);
763 const Edge half1(iVleft, iNewVert);
764 const Edge half2(iNewVert, iVright);
765
766 const unordered_map<Edge, BoundaryOverlapCount>::const_iterator
767 splitEdgeOverlapsIt = overlapCount.find(splitEdge);
768 const BoundaryOverlapCount splitEdgeOverlaps =
769 splitEdgeOverlapsIt != overlapCount.end()
770 ? splitEdgeOverlapsIt->second
771 : 0;
772 // remove the edge that will be split and add split edge's
773 // halves
774 fixedEdges.erase(splitEdge);
775 if(splitEdgeOverlaps > 0)
776 {
777 overlapCount.erase(splitEdgeOverlapsIt);
778 fixEdge(half1, splitEdgeOverlaps);
779 fixEdge(half2, splitEdgeOverlaps);
780 }
781 else
782 {
783 fixEdge(half1);
784 fixEdge(half2);
785 }
786 // maintain piece-to-original mapping
787 EdgeVec newOriginals(1, splitEdge);
788 const unordered_map<Edge, EdgeVec>::const_iterator originalsIt =
789 pieceToOriginals.find(splitEdge);
790 if(originalsIt != pieceToOriginals.end())
791 { // edge being split was split before: pass-through originals
792 newOriginals = originalsIt->second;
793 pieceToOriginals.erase(originalsIt);
794 }
795 detail::insert_unique(pieceToOriginals[half1], newOriginals);
796 detail::insert_unique(pieceToOriginals[half2], newOriginals);
797
798 // add a new point at the intersection of two constraint edges
799 const V2d<T> newV = detail::intersectionPosition(
800 vertices[iA],
801 vertices[iB],
802 vertices[iVleft],
803 vertices[iVright]);
804 addNewVertex(newV, noNeighbor);
805 std::stack<TriInd> triStack =
806 insertVertexOnEdge(iNewVert, iT, iTopo);
807 tryAddVertexToLocator(iNewVert);
808 ensureDelaunayByEdgeFlips(newV, iNewVert, triStack);
809#ifdef CDT_CXX11_IS_SUPPORTED
810 remaining.emplace_back(Edge(iNewVert, iB), originals, overlaps);
811 remaining.emplace_back(Edge(iA, iNewVert), originals, overlaps);
812#else
813 remaining.push_back(
814 make_tuple(Edge(iNewVert, iB), originals, overlaps));
815 remaining.push_back(
816 make_tuple(Edge(iA, iNewVert), originals, overlaps));
817#endif
818 return;
819 }
820
821 iT = iTopo;
822 t = triangles[iT];
823
824 const PtLineLocation::Enum loc =
825 locatePointLine(vOpo, a, b, distanceTolerance);
826 if(loc == PtLineLocation::Left)
827 {
828 iV = iVleft;
829 iVleft = iVopo;
830 }
831 else if(loc == PtLineLocation::Right)
832 {
833 iV = iVright;
834 iVright = iVopo;
835 }
836 else // encountered point on the edge
837 iB = iVopo;
838 }
839
840 // encountered one or more points on the edge: add remaining edge part
841 if(iB != edge.v2())
842 {
843#ifdef CDT_CXX11_IS_SUPPORTED
844 remaining.emplace_back(Edge(iB, edge.v2()), originals, overlaps);
845#else
846 remaining.push_back(
847 make_tuple(Edge(iB, edge.v2()), originals, overlaps));
848#endif
849 }
850
851 // add mid-point to triangulation
852 const VertInd iMid = static_cast<VertInd>(vertices.size());
853 const V2d<T>& start = vertices[iA];
854 const V2d<T>& end = vertices[iB];
855 addNewVertex(
856 V2d<T>::make((start.x + end.x) / T(2), (start.y + end.y) / T(2)),
857 noNeighbor);
858 const std::vector<Edge> flippedFixedEdges =
859 insertVertex_FlipFixedEdges(iMid);
860
861#ifdef CDT_CXX11_IS_SUPPORTED
862 remaining.emplace_back(Edge(iMid, iB), originals, overlaps);
863 remaining.emplace_back(Edge(iA, iMid), originals, overlaps);
864#else
865 remaining.push_back(make_tuple(Edge(iMid, iB), originals, overlaps));
866 remaining.push_back(make_tuple(Edge(iA, iMid), originals, overlaps));
867#endif
868
869 // re-introduce fixed edges that were flipped
870 // and make sure overlap count is preserved
871 for(std::vector<Edge>::const_iterator it = flippedFixedEdges.begin();
872 it != flippedFixedEdges.end();
873 ++it)
874 {
875 const Edge& flippedFixedEdge = *it;
876 fixedEdges.erase(flippedFixedEdge);
877
878 BoundaryOverlapCount prevOverlaps = 0;
879 const unordered_map<Edge, BoundaryOverlapCount>::const_iterator
880 overlapsIt = overlapCount.find(flippedFixedEdge);
881 if(overlapsIt != overlapCount.end())
882 {
883 prevOverlaps = overlapsIt->second;
884 overlapCount.erase(overlapsIt);
885 }
886 // override overlapping boundaries count when re-inserting an edge
887 EdgeVec prevOriginals(1, flippedFixedEdge);
888 const unordered_map<Edge, EdgeVec>::const_iterator originalsIt =
889 pieceToOriginals.find(flippedFixedEdge);
890 if(originalsIt != pieceToOriginals.end())
891 {
892 prevOriginals = originalsIt->second;
893 }
894#ifdef CDT_CXX11_IS_SUPPORTED
895 remaining.emplace_back(flippedFixedEdge, prevOriginals, prevOverlaps);
896#else
897 remaining.push_back(
898 make_tuple(flippedFixedEdge, prevOriginals, prevOverlaps));
899#endif
900 }
901}
902
903template <typename T, typename TNearPointLocator>
904void Triangulation<T, TNearPointLocator>::conformToEdge(
905 Edge edge,
906 EdgeVec originals,
907 BoundaryOverlapCount overlaps,
908 std::vector<ConformToEdgeTask>& remaining)
909{
910 // use iteration over recursion to avoid stack overflows
911 remaining.clear();
912#ifdef CDT_CXX11_IS_SUPPORTED
913 remaining.emplace_back(edge, originals, overlaps);
914#else
915 remaining.push_back(make_tuple(edge, originals, overlaps));
916#endif
917 while(!remaining.empty())
918 {
919 tie(edge, originals, overlaps) = remaining.back();
920 remaining.pop_back();
921 conformToEdgeIteration(edge, originals, overlaps, remaining);
922 }
923}
924
935template <typename T, typename TNearPointLocator>
936tuple<TriInd, VertInd, VertInd>
937Triangulation<T, TNearPointLocator>::intersectedTriangle(
938 const VertInd iA,
939 const V2d<T>& a,
940 const V2d<T>& b,
941 const T orientationTolerance) const
942{
943 const TriInd startTri = m_vertTris[iA];
944 TriInd iT = startTri;
945 do
946 {
947 const Triangle t = triangles[iT];
948 const Index i = vertexInd(t.vertices, iA);
949 const VertInd iP2 = t.vertices[ccw(i)];
950 const T orientP2 = orient2D(vertices[iP2], a, b);
951 const PtLineLocation::Enum locP2 = classifyOrientation(orientP2);
952 if(locP2 == PtLineLocation::Right)
953 {
954 const VertInd iP1 = t.vertices[cw(i)];
955 const T orientP1 = orient2D(vertices[iP1], a, b);
956 const PtLineLocation::Enum locP1 = classifyOrientation(orientP1);
957 if(locP1 == PtLineLocation::OnLine)
958 {
959 return make_tuple(noNeighbor, iP1, iP1);
960 }
961 if(locP1 == PtLineLocation::Left)
962 {
963 if(orientationTolerance)
964 {
965 T closestOrient;
966 VertInd iClosestP;
967 if(std::abs(orientP1) <= std::abs(orientP2))
968 {
969 closestOrient = orientP1;
970 iClosestP = iP1;
971 }
972 else
973 {
974 closestOrient = orientP2;
975 iClosestP = iP2;
976 }
978 closestOrient, orientationTolerance) ==
979 PtLineLocation::OnLine)
980 {
981 return make_tuple(noNeighbor, iClosestP, iClosestP);
982 }
983 }
984 return make_tuple(iT, iP1, iP2);
985 }
986 }
987 iT = t.next(iA).first;
988 } while(iT != startTri);
989 throw std::runtime_error("Could not find vertex triangle intersected by "
990 "edge. Note: can be caused by duplicate points.");
991}
992
993template <typename T, typename TNearPointLocator>
994void Triangulation<T, TNearPointLocator>::addSuperTriangle(const Box2d<T>& box)
995{
996 m_nTargetVerts = 3;
997 m_superGeomType = SuperGeometryType::SuperTriangle;
998
999 const V2d<T> center = {
1000 (box.min.x + box.max.x) / T(2), (box.min.y + box.max.y) / T(2)};
1001 const T w = box.max.x - box.min.x;
1002 const T h = box.max.y - box.min.y;
1003 T r = std::sqrt(w * w + h * h) / T(2); // incircle radius
1004 r *= T(1.1);
1005 const T R = T(2) * r; // excircle radius
1006 const T shiftX = R * std::sqrt(T(3)) / T(2); // R * cos(30 deg)
1007 const V2d<T> posV1 = {center.x - shiftX, center.y - r};
1008 const V2d<T> posV2 = {center.x + shiftX, center.y - r};
1009 const V2d<T> posV3 = {center.x, center.y + R};
1010 addNewVertex(posV1, TriInd(0));
1011 addNewVertex(posV2, TriInd(0));
1012 addNewVertex(posV3, TriInd(0));
1013 const Triangle superTri = {
1014 {VertInd(0), VertInd(1), VertInd(2)},
1015 {noNeighbor, noNeighbor, noNeighbor}};
1016 addTriangle(superTri);
1017 if(m_vertexInsertionOrder != VertexInsertionOrder::Auto)
1018 {
1019 m_nearPtLocator.initialize(vertices);
1020 }
1021}
1022
1023template <typename T, typename TNearPointLocator>
1024void Triangulation<T, TNearPointLocator>::addNewVertex(
1025 const V2d<T>& pos,
1026 const TriInd iT)
1027{
1028 vertices.push_back(pos);
1029 m_vertTris.push_back(iT);
1030}
1031
1032template <typename T, typename TNearPointLocator>
1033std::vector<Edge>
1034Triangulation<T, TNearPointLocator>::insertVertex_FlipFixedEdges(
1035 const VertInd iV1)
1036{
1037 std::vector<Edge> flippedFixedEdges;
1038
1039 const V2d<T>& v1 = vertices[iV1];
1040 const VertInd startVertex = m_nearPtLocator.nearPoint(v1, vertices);
1041 array<TriInd, 2> trisAt = walkingSearchTrianglesAt(v1, startVertex);
1042 std::stack<TriInd> triStack =
1043 trisAt[1] == noNeighbor ? insertVertexInsideTriangle(iV1, trisAt[0])
1044 : insertVertexOnEdge(iV1, trisAt[0], trisAt[1]);
1045
1046 TriInd iTopo, n1, n2, n3, n4;
1047 VertInd iV2, iV3, iV4;
1048 while(!triStack.empty())
1049 {
1050 const TriInd iT = triStack.top();
1051 triStack.pop();
1052
1053 edgeFlipInfo(iT, iV1, iTopo, iV2, iV3, iV4, n1, n2, n3, n4);
1054 if(iTopo != noNeighbor && isFlipNeeded(v1, iV1, iV2, iV3, iV4))
1055 {
1056 // if flipped edge is fixed, remember it
1057 const Edge flippedEdge(iV2, iV4);
1058 if(!fixedEdges.empty() &&
1059 fixedEdges.find(flippedEdge) != fixedEdges.end())
1060 {
1061 flippedFixedEdges.push_back(flippedEdge);
1062 }
1063
1064 flipEdge(iT, iTopo, iV1, iV2, iV3, iV4, n1, n2, n3, n4);
1065 triStack.push(iT);
1066 triStack.push(iTopo);
1067 }
1068 }
1069
1070 tryAddVertexToLocator(iV1);
1071 return flippedFixedEdges;
1072}
1073
1074template <typename T, typename TNearPointLocator>
1075void Triangulation<T, TNearPointLocator>::insertVertex(
1076 const VertInd iVert,
1077 const VertInd walkStart)
1078{
1079 const V2d<T>& v = vertices[iVert];
1080 const array<TriInd, 2> trisAt = walkingSearchTrianglesAt(v, walkStart);
1081 std::stack<TriInd> triStack =
1082 trisAt[1] == noNeighbor
1083 ? insertVertexInsideTriangle(iVert, trisAt[0])
1084 : insertVertexOnEdge(iVert, trisAt[0], trisAt[1]);
1085 ensureDelaunayByEdgeFlips(v, iVert, triStack);
1086}
1087
1088template <typename T, typename TNearPointLocator>
1089void Triangulation<T, TNearPointLocator>::insertVertex(const VertInd iVert)
1090{
1091 const V2d<T>& v = vertices[iVert];
1092 const VertInd walkStart = m_nearPtLocator.nearPoint(v, vertices);
1093 insertVertex(iVert, walkStart);
1094 tryAddVertexToLocator(iVert);
1095}
1096
1097template <typename T, typename TNearPointLocator>
1098void Triangulation<T, TNearPointLocator>::ensureDelaunayByEdgeFlips(
1099 const V2d<T>& v1,
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(v1, 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 V2d<T>& v,
1225 const VertInd iV1,
1226 const VertInd iV2,
1227 const VertInd iV3,
1228 const VertInd iV4) const
1229{
1230 if(fixedEdges.count(Edge(iV2, iV4)))
1231 return false; // flip not needed if the original edge is fixed
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(v, v3, v4);
1247 if(iV4 < 3)
1248 return locatePointLine(v4, v2, v3) ==
1249 locatePointLine(v, 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 return locatePointLine(v2, v, v4) == locatePointLine(v3, v, v4);
1257 if(iV4 < 3)
1258 return locatePointLine(v4, v2, v) == locatePointLine(v3, v2, v);
1259 return false; // original edge does not touch super-triangle
1260 }
1261 // flip-candidate edge does not touch super-triangle
1262 if(iV2 < 3)
1263 return locatePointLine(v2, v3, v4) == locatePointLine(v, v3, v4);
1264 if(iV4 < 3)
1265 return locatePointLine(v4, v2, v3) == locatePointLine(v, v2, v3);
1266 }
1267 return isInCircumcircle(v, v2, v3, v4);
1268}
1269
1270/* Flip edge between T and Topo:
1271 *
1272 * v4 | - old edge
1273 * /|\ ~ - new edge
1274 * / | \
1275 * n3 / T' \ n4
1276 * / | \
1277 * / | \
1278 * T -> v1 ~~~~~~~~ v3 <- Topo
1279 * \ | /
1280 * \ | /
1281 * n1 \Topo'/ n2
1282 * \ | /
1283 * \|/
1284 * v2
1285 */
1286template <typename T, typename TNearPointLocator>
1288 const TriInd iT,
1289 const TriInd iTopo,
1290 const VertInd v1,
1291 const VertInd v2,
1292 const VertInd v3,
1293 const VertInd v4,
1294 const TriInd n1,
1295 const TriInd n2,
1296 const TriInd n3,
1297 const TriInd n4)
1298{
1299 // change vertices and neighbors
1300 using detail::arr3;
1301 triangles[iT] = Triangle::make(arr3(v4, v1, v3), arr3(n3, iTopo, n4));
1302 triangles[iTopo] = Triangle::make(arr3(v2, v3, v1), arr3(n2, iT, n1));
1303 // adjust neighboring triangles and vertices
1304 changeNeighbor(n1, iT, iTopo);
1305 changeNeighbor(n4, iTopo, iT);
1306 // only adjust adjacent triangles if triangulation is not finalized:
1307 // can happen when called from outside on an already finalized
1308 // triangulation
1309 if(!isFinalized())
1310 {
1311 setAdjacentTriangle(v4, iT);
1312 setAdjacentTriangle(v2, iTopo);
1313 }
1314}
1315
1316/* Insert point into triangle: split into 3 triangles:
1317 * - create 2 new triangles
1318 * - re-use old triangle for the 3rd
1319 * v3
1320 * / | \
1321 * / | \ <-- original triangle (t)
1322 * / | \
1323 * n3 / | \ n2
1324 * /newT2|newT1\
1325 * / v \
1326 * / __/ \__ \
1327 * / __/ \__ \
1328 * / _/ t' \_ \
1329 * v1 ___________________ v2
1330 * n1
1331 */
1332template <typename T, typename TNearPointLocator>
1333std::stack<TriInd>
1334Triangulation<T, TNearPointLocator>::insertVertexInsideTriangle(
1335 VertInd v,
1336 TriInd iT)
1337{
1338 const TriInd iNewT1 = addTriangle();
1339 const TriInd iNewT2 = addTriangle();
1340
1341 Triangle& t = triangles[iT];
1342 const array<VertInd, 3> vv = t.vertices;
1343 const array<TriInd, 3> nn = t.neighbors;
1344 const VertInd v1 = vv[0], v2 = vv[1], v3 = vv[2];
1345 const TriInd n1 = nn[0], n2 = nn[1], n3 = nn[2];
1346 // make two new triangles and convert current triangle to 3rd new
1347 // triangle
1348 using detail::arr3;
1349 triangles[iNewT1] = Triangle::make(arr3(v2, v3, v), arr3(n2, iNewT2, iT));
1350 triangles[iNewT2] = Triangle::make(arr3(v3, v1, v), arr3(n3, iT, iNewT1));
1351 t = Triangle::make(arr3(v1, v2, v), arr3(n1, iNewT1, iNewT2));
1352 // adjust adjacent triangles
1353 setAdjacentTriangle(v, iT);
1354 setAdjacentTriangle(v3, iNewT1);
1355 // change triangle neighbor's neighbors to new triangles
1356 changeNeighbor(n2, iT, iNewT1);
1357 changeNeighbor(n3, iT, iNewT2);
1358 // return newly added triangles
1359 std::stack<TriInd> newTriangles;
1360 newTriangles.push(iT);
1361 newTriangles.push(iNewT1);
1362 newTriangles.push(iNewT2);
1363 return newTriangles;
1364}
1365
1366/* Inserting a point on the edge between two triangles
1367 * T1 (top) v1
1368 * /|\
1369 * n1 / | \ n4
1370 * / | \
1371 * / T1' | Tnew1\
1372 * v2-------v-------v4
1373 * \ T2' | Tnew2/
1374 * \ | /
1375 * n2 \ | / n3
1376 * \|/
1377 * T2 (bottom) v3
1378 */
1379template <typename T, typename TNearPointLocator>
1380std::stack<TriInd> Triangulation<T, TNearPointLocator>::insertVertexOnEdge(
1381 VertInd v,
1382 TriInd iT1,
1383 TriInd iT2)
1384{
1385 const TriInd iTnew1 = addTriangle();
1386 const TriInd iTnew2 = addTriangle();
1387
1388 Triangle& t1 = triangles[iT1];
1389 Triangle& t2 = triangles[iT2];
1390 Index i = opposedVertexInd(t1.neighbors, iT2);
1391 const VertInd v1 = t1.vertices[i];
1392 const VertInd v2 = t1.vertices[ccw(i)];
1393 const TriInd n1 = t1.neighbors[i];
1394 const TriInd n4 = t1.neighbors[cw(i)];
1395 i = opposedVertexInd(t2.neighbors, iT1);
1396 const VertInd v3 = t2.vertices[i];
1397 const VertInd v4 = t2.vertices[ccw(i)];
1398 const TriInd n3 = t2.neighbors[i];
1399 const TriInd n2 = t2.neighbors[cw(i)];
1400 // add new triangles and change existing ones
1401 using detail::arr3;
1402 t1 = Triangle::make(arr3(v, v1, v2), arr3(iTnew1, n1, iT2));
1403 t2 = Triangle::make(arr3(v, v2, v3), arr3(iT1, n2, iTnew2));
1404 triangles[iTnew1] = Triangle::make(arr3(v, v4, v1), arr3(iTnew2, n4, iT1));
1405 triangles[iTnew2] = Triangle::make(arr3(v, v3, v4), arr3(iT2, n3, iTnew1));
1406 // adjust adjacent triangles
1407 setAdjacentTriangle(v, iT1);
1408 setAdjacentTriangle(v4, iTnew1);
1409 // adjust neighboring triangles and vertices
1410 changeNeighbor(n4, iT1, iTnew1);
1411 changeNeighbor(n3, iT2, iTnew2);
1412 // return newly added triangles
1413 std::stack<TriInd> newTriangles;
1414 newTriangles.push(iT1);
1415 newTriangles.push(iTnew2);
1416 newTriangles.push(iT2);
1417 newTriangles.push(iTnew1);
1418 return newTriangles;
1419}
1420
1421template <typename T, typename TNearPointLocator>
1422array<TriInd, 2>
1423Triangulation<T, TNearPointLocator>::trianglesAt(const V2d<T>& pos) const
1424{
1425 array<TriInd, 2> out = {noNeighbor, noNeighbor};
1426 for(TriInd i = TriInd(0); i < TriInd(triangles.size()); ++i)
1427 {
1428 const Triangle& t = triangles[i];
1429 const V2d<T>& v1 = vertices[t.vertices[0]];
1430 const V2d<T>& v2 = vertices[t.vertices[1]];
1431 const V2d<T>& v3 = vertices[t.vertices[2]];
1432 const PtTriLocation::Enum loc = locatePointTriangle(pos, v1, v2, v3);
1433 if(loc == PtTriLocation::Outside)
1434 continue;
1435 out[0] = i;
1436 if(isOnEdge(loc))
1437 out[1] = t.neighbors[edgeNeighbor(loc)];
1438 return out;
1439 }
1440 throw std::runtime_error("No triangle was found at position");
1441}
1442
1443template <typename T, typename TNearPointLocator>
1444TriInd Triangulation<T, TNearPointLocator>::walkTriangles(
1445 const VertInd startVertex,
1446 const V2d<T>& pos) const
1447{
1448 // begin walk in search of triangle at pos
1449 TriInd currTri = m_vertTris[startVertex];
1450 bool found = false;
1451 detail::SplitMix64RandGen prng;
1452 while(!found)
1453 {
1454 const Triangle& t = triangles[currTri];
1455 found = true;
1456 // stochastic offset to randomize which edge we check first
1457 const Index offset(prng() % 3);
1458 for(Index i_(0); i_ < Index(3); ++i_)
1459 {
1460 const Index i((i_ + offset) % 3);
1461 const V2d<T>& vStart = vertices[t.vertices[i]];
1462 const V2d<T>& vEnd = vertices[t.vertices[ccw(i)]];
1463 const PtLineLocation::Enum edgeCheck =
1464 locatePointLine(pos, vStart, vEnd);
1465 const TriInd iN = t.neighbors[i];
1466 if(edgeCheck == PtLineLocation::Right && iN != noNeighbor)
1467 {
1468 found = false;
1469 currTri = t.neighbors[i];
1470 break;
1471 }
1472 }
1473 }
1474 return currTri;
1475}
1476
1477template <typename T, typename TNearPointLocator>
1478array<TriInd, 2> Triangulation<T, TNearPointLocator>::walkingSearchTrianglesAt(
1479 const V2d<T>& pos,
1480 const VertInd startVertex) const
1481{
1482 array<TriInd, 2> out = {noNeighbor, noNeighbor};
1483 const TriInd iT = walkTriangles(startVertex, pos);
1484 // Finished walk, locate point in current triangle
1485 const Triangle& t = triangles[iT];
1486 const V2d<T>& v1 = vertices[t.vertices[0]];
1487 const V2d<T>& v2 = vertices[t.vertices[1]];
1488 const V2d<T>& v3 = vertices[t.vertices[2]];
1489 const PtTriLocation::Enum loc = locatePointTriangle(pos, v1, v2, v3);
1490 if(loc == PtTriLocation::Outside)
1491 throw std::runtime_error("No triangle was found at position");
1492 out[0] = iT;
1493 if(isOnEdge(loc))
1494 out[1] = t.neighbors[edgeNeighbor(loc)];
1495 return out;
1496}
1497
1498/* Flip edge between T and Topo:
1499 *
1500 * v4 | - old edge
1501 * /|\ ~ - new edge
1502 * / | \
1503 * n3 / T' \ n4
1504 * / | \
1505 * / | \
1506 * T -> v1~~~~~~~~~v3 <- Topo
1507 * \ | /
1508 * \ | /
1509 * n1 \Topo'/ n2
1510 * \ | /
1511 * \|/
1512 * v2
1513 */
1514template <typename T, typename TNearPointLocator>
1516 const TriInd iT,
1517 const TriInd iTopo)
1518{
1519 Triangle& t = triangles[iT];
1520 Triangle& tOpo = triangles[iTopo];
1521 const array<TriInd, 3>& triNs = t.neighbors;
1522 const array<TriInd, 3>& triOpoNs = tOpo.neighbors;
1523 const array<VertInd, 3>& triVs = t.vertices;
1524 const array<VertInd, 3>& triOpoVs = tOpo.vertices;
1525 // find vertices and neighbors
1526 Index i = opposedVertexInd(t.neighbors, iTopo);
1527 const VertInd v1 = triVs[i];
1528 const VertInd v2 = triVs[ccw(i)];
1529 const TriInd n1 = triNs[i];
1530 const TriInd n3 = triNs[cw(i)];
1531 i = opposedVertexInd(tOpo.neighbors, iT);
1532 const VertInd v3 = triOpoVs[i];
1533 const VertInd v4 = triOpoVs[ccw(i)];
1534 const TriInd n4 = triOpoNs[i];
1535 const TriInd n2 = triOpoNs[cw(i)];
1536 // change vertices and neighbors
1537 using detail::arr3;
1538 t = Triangle::make(arr3(v4, v1, v3), arr3(n3, iTopo, n4));
1539 tOpo = Triangle::make(arr3(v2, v3, v1), arr3(n2, iT, n1));
1540 // adjust neighboring triangles and vertices
1541 changeNeighbor(n1, iT, iTopo);
1542 changeNeighbor(n4, iTopo, iT);
1543 // only adjust adjacent triangles if triangulation is not finalized:
1544 // can happen when called from outside on an already finalized
1545 // triangulation
1546 if(!isFinalized())
1547 {
1548 setAdjacentTriangle(v4, iT);
1549 setAdjacentTriangle(v2, iTopo);
1550 }
1551}
1552
1553template <typename T, typename TNearPointLocator>
1555 const TriInd iT,
1556 const TriInd oldNeighbor,
1557 const TriInd newNeighbor)
1558{
1559 if(iT == noNeighbor)
1560 return;
1561 NeighborsArr3& nn = triangles[iT].neighbors;
1562 assert(
1563 nn[0] == oldNeighbor || nn[1] == oldNeighbor || nn[2] == oldNeighbor);
1564 if(nn[0] == oldNeighbor)
1565 nn[0] = newNeighbor;
1566 else if(nn[1] == oldNeighbor)
1567 nn[1] = newNeighbor;
1568 else
1569 nn[2] = newNeighbor;
1570}
1571
1572template <typename T, typename TNearPointLocator>
1573void Triangulation<T, TNearPointLocator>::changeNeighbor(
1574 const TriInd iT,
1575 const VertInd iVedge1,
1576 const VertInd iVedge2,
1577 const TriInd newNeighbor)
1578{
1579 assert(iT != noNeighbor);
1580 Triangle& t = triangles[iT];
1581 t.neighbors[edgeNeighborInd(t.vertices, iVedge1, iVedge2)] = newNeighbor;
1582}
1583
1584template <typename T, typename TNearPointLocator>
1585void Triangulation<T, TNearPointLocator>::triangulatePseudopolygon(
1586 const std::vector<VertInd>& poly,
1587 const std::vector<TriInd>& outerTris,
1588 const TriInd iT,
1589 const TriInd iN,
1590 std::vector<TriangulatePseudopolygonTask>& iterations)
1591{
1592 assert(poly.size() > 2);
1593 // note: uses interation instead of recursion to avoid stack overflows
1594 iterations.clear();
1595 iterations.push_back(make_tuple(0, poly.size() - 1, iT, iN, Index(0)));
1596 while(!iterations.empty())
1597 {
1598 triangulatePseudopolygonIteration(poly, outerTris, iterations);
1599 }
1600}
1601
1602template <typename T, typename TNearPointLocator>
1603void Triangulation<T, TNearPointLocator>::triangulatePseudopolygonIteration(
1604 const std::vector<VertInd>& poly,
1605 const std::vector<TriInd>& outerTris,
1606 std::vector<TriangulatePseudopolygonTask>& iterations)
1607{
1608 IndexSizeType iA, iB;
1609 TriInd iT, iParent;
1610 Index iInParent;
1611 assert(!iterations.empty());
1612 tie(iA, iB, iT, iParent, iInParent) = iterations.back();
1613 iterations.pop_back();
1614 assert(iB - iA > 1 && iT != noNeighbor && iParent != noNeighbor);
1615 Triangle& t = triangles[iT];
1616 // find Delaunay point
1617 const IndexSizeType iC = findDelaunayPoint(poly, iA, iB);
1618
1619 const VertInd a = poly[iA];
1620 const VertInd b = poly[iB];
1621 const VertInd c = poly[iC];
1622 // split pseudo-polygon in two parts and triangulate them
1623 //
1624 // note: first part needs to be pushed on stack last
1625 // in order to be processed first
1626 //
1627 // second part: points after the Delaunay point
1628 if(iB - iC > 1)
1629 {
1630 const TriInd iNext = addTriangle();
1631 iterations.push_back(make_tuple(iC, iB, iNext, iT, Index(1)));
1632 }
1633 else // pseudo-poly is reduced to a single outer edge
1634 {
1635 const TriInd outerTri = outerTris[iC];
1636 if(outerTri != noNeighbor)
1637 {
1638 assert(outerTri != iT);
1639 t.neighbors[1] = outerTri;
1640 changeNeighbor(outerTri, c, b, iT);
1641 }
1642 }
1643 // first part: points before the Delaunay point
1644 if(iC - iA > 1)
1645 { // add next triangle and add another iteration
1646 const TriInd iNext = addTriangle();
1647 iterations.push_back(make_tuple(iA, iC, iNext, iT, Index(2)));
1648 }
1649 else
1650 { // pseudo-poly is reduced to a single outer edge
1651 const TriInd outerTri =
1652 outerTris[iA] != noNeighbor ? outerTris[iA] : m_vertTris[c];
1653 if(outerTri != noNeighbor)
1654 {
1655 assert(outerTri != iT);
1656 t.neighbors[2] = outerTri;
1657 changeNeighbor(outerTri, c, a, iT);
1658 }
1659 }
1660 // Finalize triangle
1661 // note: only when triangle is finalized to we add it as a neighbor to
1662 // parent to maintain triangulation topology consistency
1663 triangles[iParent].neighbors[iInParent] = iT;
1664 t.neighbors[0] = iParent;
1665 t.vertices = detail::arr3(a, b, c);
1666 // needs to be done at the end not to affect finding edge triangles
1667 setAdjacentTriangle(c, iT);
1668}
1669
1670template <typename T, typename TNearPointLocator>
1671IndexSizeType Triangulation<T, TNearPointLocator>::findDelaunayPoint(
1672 const std::vector<VertInd>& poly,
1673 const IndexSizeType iA,
1674 const IndexSizeType iB) const
1675{
1676 assert(iB - iA > 1);
1677 const V2d<T>& a = vertices[poly[iA]];
1678 const V2d<T>& b = vertices[poly[iB]];
1679 IndexSizeType out = iA + 1;
1680 const V2d<T>* c = &vertices[poly[out]]; // caching for better performance
1681 for(IndexSizeType i = iA + 1; i < iB; ++i)
1682 {
1683 const V2d<T>& v = vertices[poly[i]];
1684 if(isInCircumcircle(v, a, b, *c))
1685 {
1686 out = i;
1687 c = &v;
1688 }
1689 }
1690 assert(out > iA && out < iB); // point is between ends
1691 return out;
1692}
1693
1694template <typename T, typename TNearPointLocator>
1696 const std::vector<V2d<T> >& newVertices)
1697{
1698 return insertVertices(
1699 newVertices.begin(), newVertices.end(), getX_V2d<T>, getY_V2d<T>);
1700}
1701
1702template <typename T, typename TNearPointLocator>
1704{
1705 return m_vertTris.empty() && !vertices.empty();
1706}
1707
1708template <typename T, typename TNearPointLocator>
1709unordered_map<TriInd, LayerDepth>
1711 std::stack<TriInd> seeds,
1712 const LayerDepth layerDepth,
1713 std::vector<LayerDepth>& triDepths) const
1714{
1715 unordered_map<TriInd, LayerDepth> behindBoundary;
1716 while(!seeds.empty())
1717 {
1718 const TriInd iT = seeds.top();
1719 seeds.pop();
1720 triDepths[iT] = layerDepth;
1721 behindBoundary.erase(iT);
1722 const Triangle& t = triangles[iT];
1723 for(Index i(0); i < Index(3); ++i)
1724 {
1725 const Edge opEdge(t.vertices[ccw(i)], t.vertices[cw(i)]);
1726 const TriInd iN = t.neighbors[opoNbr(i)];
1727 if(iN == noNeighbor || triDepths[iN] <= layerDepth)
1728 continue;
1729 if(fixedEdges.count(opEdge))
1730 {
1731 const unordered_map<Edge, LayerDepth>::const_iterator cit =
1732 overlapCount.find(opEdge);
1733 const LayerDepth triDepth = cit == overlapCount.end()
1734 ? layerDepth + 1
1735 : layerDepth + cit->second + 1;
1736 behindBoundary[iN] = triDepth;
1737 continue;
1738 }
1739 seeds.push(iN);
1740 }
1741 }
1742 return behindBoundary;
1743}
1744
1745template <typename T, typename TNearPointLocator>
1746std::vector<LayerDepth>
1748{
1749 std::vector<LayerDepth> triDepths(
1750 triangles.size(), std::numeric_limits<LayerDepth>::max());
1751 std::stack<TriInd> seeds(TriDeque(1, m_vertTris[0]));
1752 LayerDepth layerDepth = 0;
1753 LayerDepth deepestSeedDepth = 0;
1754
1755 unordered_map<LayerDepth, TriIndUSet> seedsByDepth;
1756 do
1757 {
1758 const unordered_map<TriInd, LayerDepth>& newSeeds =
1759 peelLayer(seeds, layerDepth, triDepths);
1760
1761 seedsByDepth.erase(layerDepth);
1762 typedef unordered_map<TriInd, LayerDepth>::const_iterator Iter;
1763 for(Iter it = newSeeds.begin(); it != newSeeds.end(); ++it)
1764 {
1765 deepestSeedDepth = std::max(deepestSeedDepth, it->second);
1766 seedsByDepth[it->second].insert(it->first);
1767 }
1768 const TriIndUSet& nextLayerSeeds = seedsByDepth[layerDepth + 1];
1769 seeds = std::stack<TriInd>(
1770 TriDeque(nextLayerSeeds.begin(), nextLayerSeeds.end()));
1771 ++layerDepth;
1772 } while(!seeds.empty() || deepestSeedDepth > layerDepth);
1773
1774 return triDepths;
1775}
1776
1777template <typename T, typename TNearPointLocator>
1779 VertInd superGeomVertCount)
1780{
1781 for(VertInd iV = superGeomVertCount; iV < vertices.size(); ++iV)
1782 {
1783 insertVertex(iV);
1784 }
1785}
1786
1787template <typename T, typename TNearPointLocator>
1788void Triangulation<T, TNearPointLocator>::insertVertices_Randomized(
1789 VertInd superGeomVertCount)
1790{
1791 std::size_t vertexCount = vertices.size() - superGeomVertCount;
1792 std::vector<VertInd> ii(vertexCount);
1793 detail::iota(ii.begin(), ii.end(), superGeomVertCount);
1794 detail::random_shuffle(ii.begin(), ii.end());
1795 for(std::vector<VertInd>::iterator it = ii.begin(); it != ii.end(); ++it)
1796 {
1797 insertVertex(*it);
1798 }
1799}
1800
1801namespace detail
1802{
1803
1804// log2 implementation backwards compatible with pre c++11
1805template <typename T>
1806inline double log2_bc(T x)
1807{
1808#ifdef CDT_CXX11_IS_SUPPORTED
1809 return std::log2(x);
1810#else
1811 static double log2_constant = std::log(2.0);
1812 return std::log(static_cast<double>(x)) / log2_constant;
1813#endif
1814}
1815
1820inline std::size_t maxQueueLengthBFSKDTree(const std::size_t vertexCount)
1821{
1822 int filledLayerPow2 = std::floor(log2_bc(vertexCount)) - 1;
1823 std::size_t nodesInFilledTree = std::pow(2., filledLayerPow2 + 1) - 1;
1824 std::size_t nodesInLastFilledLayer = std::pow(2., filledLayerPow2);
1825 std::size_t nodesInLastLayer = vertexCount - nodesInFilledTree;
1826 return nodesInLastLayer >= nodesInLastFilledLayer
1827 ? nodesInLastFilledLayer + nodesInLastLayer -
1828 nodesInLastFilledLayer
1829 : nodesInLastFilledLayer;
1830}
1831
1832template <typename T>
1834{
1835public:
1836 FixedCapacityQueue(const std::size_t capacity)
1837 : m_vec(capacity)
1838 , m_front(m_vec.begin())
1839 , m_back(m_vec.begin())
1840 , m_size(0)
1841 {}
1842 bool empty() const
1843 {
1844 return m_size == 0;
1845 }
1846 const T& front() const
1847 {
1848 return *m_front;
1849 }
1850 void pop()
1851 {
1852 assert(m_size > 0);
1853 ++m_front;
1854 if(m_front == m_vec.end())
1855 m_front = m_vec.begin();
1856 --m_size;
1857 }
1858 void push(const T& t)
1859 {
1860 assert(m_size < m_vec.size());
1861 *m_back = t;
1862 ++m_back;
1863 if(m_back == m_vec.end())
1864 m_back = m_vec.begin();
1865 ++m_size;
1866 }
1867#ifdef CDT_CXX11_IS_SUPPORTED
1868 void push(const T&& t)
1869 {
1870 assert(m_size < m_vec.size());
1871 *m_back = t;
1872 ++m_back;
1873 if(m_back == m_vec.end())
1874 m_back = m_vec.begin();
1875 ++m_size;
1876 }
1877#endif
1878private:
1879 std::vector<T> m_vec;
1880 typename std::vector<T>::iterator m_front;
1881 typename std::vector<T>::iterator m_back;
1882 std::size_t m_size;
1883};
1884
1885template <typename T>
1887{
1888 less_than_x(const std::vector<V2d<T> >& vertices)
1889 : vertices(vertices)
1890 {}
1891 bool operator()(const VertInd a, const VertInd b) const
1892 {
1893 return vertices[a].x < vertices[b].x;
1894 }
1895 const std::vector<V2d<T> >& vertices;
1896};
1897
1898template <typename T>
1900{
1901 less_than_y(const std::vector<V2d<T> >& vertices)
1902 : vertices(vertices)
1903 {}
1904 bool operator()(const VertInd a, const VertInd b) const
1905 {
1906 return vertices[a].y < vertices[b].y;
1907 }
1908 const std::vector<V2d<T> >& vertices;
1909};
1910
1911} // namespace detail
1912
1913template <typename T, typename TNearPointLocator>
1915 VertInd superGeomVertCount,
1916 V2d<T> boxMin,
1917 V2d<T> boxMax)
1918{
1919 // calculate original indices
1920 std::size_t vertexCount = vertices.size() - superGeomVertCount;
1921 std::vector<VertInd> ii(vertexCount);
1922 detail::iota(ii.begin(), ii.end(), superGeomVertCount);
1923
1924 typedef std::vector<VertInd>::iterator It;
1926 detail::maxQueueLengthBFSKDTree(vertexCount));
1927 queue.push(make_tuple(ii.begin(), ii.end(), boxMin, boxMax, VertInd(0)));
1928
1929 It first, last;
1930 V2d<T> newBoxMin, newBoxMax;
1931 VertInd parent, mid;
1932
1933 const detail::less_than_x<T> cmpX(vertices);
1934 const detail::less_than_y<T> cmpY(vertices);
1935
1936 while(!queue.empty())
1937 {
1938 tie(first, last, boxMin, boxMax, parent) = queue.front();
1939 queue.pop();
1940 assert(first != last);
1941
1942 const std::ptrdiff_t len = std::distance(first, last);
1943 if(len == 1)
1944 {
1945 insertVertex(*first, parent);
1946 continue;
1947 }
1948 const It midIt = first + len / 2;
1949 if(boxMax.x - boxMin.x >= boxMax.y - boxMin.y)
1950 {
1951 detail::portable_nth_element(first, midIt, last, cmpX);
1952 mid = *midIt;
1953 const T split = vertices[mid].x;
1954 newBoxMin.x = split;
1955 newBoxMin.y = boxMin.y;
1956 newBoxMax.x = split;
1957 newBoxMax.y = boxMax.y;
1958 }
1959 else
1960 {
1961 detail::portable_nth_element(first, midIt, last, cmpY);
1962 mid = *midIt;
1963 const T split = vertices[mid].y;
1964 newBoxMin.x = boxMin.x;
1965 newBoxMin.y = split;
1966 newBoxMax.x = boxMax.x;
1967 newBoxMax.y = split;
1968 }
1969 insertVertex(mid, parent);
1970 if(first != midIt)
1971 {
1972 queue.push(make_tuple(first, midIt, boxMin, newBoxMax, mid));
1973 }
1974 if(midIt + 1 != last)
1975 {
1976 queue.push(make_tuple(midIt + 1, last, newBoxMin, boxMax, mid));
1977 }
1978 }
1979}
1980
1981template <typename T, typename TNearPointLocator>
1982bool Triangulation<T, TNearPointLocator>::hasEdge(
1983 const VertInd a,
1984 const VertInd b) const
1985{
1986 const TriInd triStart = m_vertTris[a];
1987 assert(triStart != noNeighbor);
1988 TriInd iT = triStart;
1989 VertInd iV = noVertex;
1990 do
1991 {
1992 const Triangle& t = triangles[iT];
1993 tie(iT, iV) = t.next(a);
1994 assert(iT != noNeighbor);
1995 if(iV == b)
1996 return true;
1997 } while(iT != triStart);
1998 return false;
1999}
2000
2001template <typename T, typename TNearPointLocator>
2002void Triangulation<T, TNearPointLocator>::setAdjacentTriangle(
2003 const VertInd v,
2004 const TriInd t)
2005{
2006 assert(t != noNeighbor);
2007 m_vertTris[v] = t;
2008 assert(
2009 triangles[t].vertices[0] == v || triangles[t].vertices[1] == v ||
2010 triangles[t].vertices[2] == v);
2011}
2012
2013template <typename T, typename TNearPointLocator>
2014void Triangulation<T, TNearPointLocator>::pivotVertexTriangleCW(const VertInd v)
2015{
2016 assert(m_vertTris[v] != noNeighbor);
2017 m_vertTris[v] = triangles[m_vertTris[v]].next(v).first;
2018 assert(m_vertTris[v] != noNeighbor);
2019 assert(
2020 triangles[m_vertTris[v]].vertices[0] == v ||
2021 triangles[m_vertTris[v]].vertices[1] == v ||
2022 triangles[m_vertTris[v]].vertices[2] == v);
2023}
2024
2025template <typename T, typename TNearPointLocator>
2026void Triangulation<T, TNearPointLocator>::removeAdjacentTriangle(
2027 const VertInd v)
2028{
2029 m_vertTris[v] = noNeighbor;
2030}
2031
2032template <typename T, typename TNearPointLocator>
2033void Triangulation<T, TNearPointLocator>::tryAddVertexToLocator(const VertInd v)
2034{
2035 if(!m_nearPtLocator.empty()) // only if locator is initialized already
2036 m_nearPtLocator.addPoint(v, vertices);
2037}
2038
2039template <typename T, typename TNearPointLocator>
2040void Triangulation<T, TNearPointLocator>::tryInitNearestPointLocator()
2041{
2042 if(!vertices.empty() && m_nearPtLocator.empty())
2043 {
2044 m_nearPtLocator.initialize(vertices);
2045 }
2046}
2047
2048} // 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)
Ensure that triangulation conforms to constraints (fixed edges)
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 constraints (custom-type fixed edges) into triangulation.
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
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:243
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:183
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:114
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:290
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_EXPORT CDT_INLINE_IF_HEADER_ONLY Index opoNbr(Index vertIndex)
Opposed neighbor index from vertex index.
Definition: CDTUtils.hpp:150
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:226
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:97
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:103
CDT_EXPORT T distanceSquared(const V2d< T > &a, const V2d< T > &b)
Squared distance between two 2D points.
Definition: CDTUtils.hpp:296
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:90
unordered_map< TriInd, TriInd > TriIndUMap
Triangle hash map.
Definition: CDTUtils.h:248
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:237
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:256
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:124
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:215
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
@ Ignore
constraint edge intersections are not checked
Definition: Triangulation.h:81
@ Resolve
constraint edge intersections are resolved
Definition: Triangulation.h:82
Enum
The Enum itself.
Definition: Triangulation.h:64
@ SuperTriangle
conventional super-triangle
Definition: Triangulation.h:65
@ Custom
user-specified custom geometry (e.g., grid)
Definition: Triangulation.h:66
Triangulation triangle (counter-clockwise winding)
Definition: CDTUtils.h:257
VerticesArr3 vertices
triangle's three vertices
Definition: CDTUtils.h:258
static Triangle make(const array< VertInd, 3 > &vertices, const array< TriInd, 3 > &neighbors)
Factory method.
Definition: CDTUtils.h:267
NeighborsArr3 neighbors
triangle's three neighbors
Definition: CDTUtils.h:259
std::pair< TriInd, VertInd > next(const VertInd i) const
Next triangle adjacent to a vertex (clockwise)
Definition: CDTUtils.h:275
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.
Definition: Triangulation.h:50