forked from jmespadero/pyDelaunay2D
-
Notifications
You must be signed in to change notification settings - Fork 0
/
delaunay2D.py
233 lines (198 loc) · 9.13 KB
/
delaunay2D.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
# -*- coding: ascii -*-
"""
Simple structured Delaunay triangulation in 2D with Bowyer-Watson algorithm.
Written by Jose M. Espadero ( http://github.com/jmespadero/pyDelaunay2D )
Based on code from Ayron Catteau. Published at http://github.com/ayron/delaunay
Just pretend to be simple and didactic. The only requisite is numpy.
Robust checks disabled by default. May not work in degenerate set of points.
"""
import numpy as np
from math import sqrt
class Delaunay2D:
"""
Class to compute a Delaunay triangulation in 2D
ref: http://en.wikipedia.org/wiki/Bowyer-Watson_algorithm
ref: http://www.geom.uiuc.edu/~samuelp/del_project.html
"""
def __init__(self, center=(0, 0), radius=9999):
""" Init and create a new frame to contain the triangulation
center -- Optional position for the center of the frame. Default (0,0)
radius -- Optional distance from corners to the center.
"""
center = np.asarray(center)
# Create coordinates for the corners of the frame
self.coords = [center+radius*np.array((-1, -1)),
center+radius*np.array((+1, -1)),
center+radius*np.array((+1, +1)),
center+radius*np.array((-1, +1))]
# Create two dicts to store triangle neighbours and circumcircles.
self.triangles = {}
self.circles = {}
# Create two CCW triangles for the frame
T1 = (0, 1, 3)
T2 = (2, 3, 1)
self.triangles[T1] = [T2, None, None]
self.triangles[T2] = [T1, None, None]
# Compute circumcenters and circumradius for each triangle
for t in self.triangles:
self.circles[t] = self.circumcenter(t)
def circumcenter(self, tri):
"""Compute circumcenter and circumradius of a triangle in 2D.
Uses an extension of the method described here:
http://www.ics.uci.edu/~eppstein/junkyard/circumcenter.html
"""
pts = np.asarray([self.coords[v] for v in tri])
pts2 = np.dot(pts, pts.T)
A = np.bmat([[2 * pts2, [[1],
[1],
[1]]],
[[[1, 1, 1, 0]]]])
b = np.hstack((np.sum(pts * pts, axis=1), [1]))
x = np.linalg.solve(A, b)
bary_coords = x[:-1]
center = np.dot(bary_coords, pts)
# radius = np.linalg.norm(pts[0] - center) # euclidean distance
radius = np.sum(np.square(pts[0] - center)) # squared distance
return (center, radius)
def inCircleFast(self, tri, p):
"""Check if point p is inside of precomputed circumcircle of tri.
"""
center, radius = self.circles[tri]
return np.sum(np.square(center - p)) <= radius
def inCircleRobust(self, tri, p):
"""Check if point p is inside of circumcircle around the triangle tri.
This is a robust predicate, slower than compare distance to centers
ref: http://www.cs.cmu.edu/~quake/robust.html
"""
m1 = np.asarray([self.coords[v] - p for v in tri])
m2 = np.sum(np.square(m1), axis=1).reshape((3, 1))
m = np.hstack((m1, m2)) # The 3x3 matrix to check
return np.linalg.det(m) <= 0
def addPoint(self, p):
"""Add a point to the current DT, and refine it using Bowyer-Watson.
"""
p = np.asarray(p)
idx = len(self.coords)
# print("coords[", idx,"] ->",p)
self.coords.append(p)
# Search the triangle(s) whose circumcircle contains p
bad_triangles = []
for T in self.triangles:
# Choose one method: inCircleRobust(T, p) or inCircleFast(T, p)
if self.inCircleFast(T, p):
bad_triangles.append(T)
# Find the CCW boundary (star shape) of the bad triangles,
# expressed as a list of edges (point pairs) and the opposite
# triangle to each edge.
boundary = []
# Choose a "random" triangle and edge
T = bad_triangles[0]
edge = 0
# get the opposite triangle of this edge
while True:
# Check if edge of triangle T is on the boundary...
# if opposite triangle of this edge is external to the list
tri_op = self.triangles[T][edge]
if tri_op not in bad_triangles:
# Insert edge and external triangle into boundary list
boundary.append((T[(edge+1) % 3], T[(edge-1) % 3], tri_op))
# Move to next CCW edge in this triangle
edge = (edge + 1) % 3
# Check if boundary is a closed loop
if boundary[0][0] == boundary[-1][1]:
break
else:
# Move to next CCW edge in opposite triangle
edge = (self.triangles[tri_op].index(T) + 1) % 3
T = tri_op
# Remove triangles too near of point p of our solution
for T in bad_triangles:
del self.triangles[T]
del self.circles[T]
# Retriangle the hole left by bad_triangles
new_triangles = []
for (e0, e1, tri_op) in boundary:
# Create a new triangle using point p and edge extremes
T = (idx, e0, e1)
# Store circumcenter and circumradius of the triangle
self.circles[T] = self.circumcenter(T)
# Set opposite triangle of the edge as neighbour of T
self.triangles[T] = [tri_op, None, None]
# Try to set T as neighbour of the opposite triangle
if tri_op:
# search the neighbour of tri_op that use edge (e1, e0)
for i, neigh in enumerate(self.triangles[tri_op]):
if neigh:
if e1 in neigh and e0 in neigh:
# change link to use our new triangle
self.triangles[tri_op][i] = T
# Add triangle to a temporal list
new_triangles.append(T)
# Link the new triangles each another
N = len(new_triangles)
for i, T in enumerate(new_triangles):
self.triangles[T][1] = new_triangles[(i+1) % N] # next
self.triangles[T][2] = new_triangles[(i-1) % N] # previous
def exportTriangles(self):
"""Export the current list of Delaunay triangles
"""
# Filter out triangles with any vertex in the extended BBox
return [(a-4, b-4, c-4)
for (a, b, c) in self.triangles if a > 3 and b > 3 and c > 3]
def exportCircles(self):
"""Export the circumcircles as a list of (center, radius)
"""
# Remember to compute circumcircles if not done before
# for t in self.triangles:
# self.circles[t] = self.circumcenter(t)
# Filter out triangles with any vertex in the extended BBox
# Do sqrt of radius before of return
return [(self.circles[(a, b, c)][0], sqrt(self.circles[(a, b, c)][1]))
for (a, b, c) in self.triangles if a > 3 and b > 3 and c > 3]
def exportDT(self):
"""Export the current set of Delaunay coordinates and triangles.
"""
# Filter out coordinates in the extended BBox
coord = self.coords[4:]
# Filter out triangles with any vertex in the extended BBox
tris = [(a-4, b-4, c-4)
for (a, b, c) in self.triangles if a > 3 and b > 3 and c > 3]
return coord, tris
def exportExtendedDT(self):
"""Export the Extended Delaunay Triangulation (with the frame vertex).
"""
return self.coords, list(self.triangles)
def exportVoronoiRegions(self):
"""Export coordinates and regions of Voronoi diagram as indexed data.
"""
# Remember to compute circumcircles if not done before
# for t in self.triangles:
# self.circles[t] = self.circumcenter(t)
useVertex = {i: [] for i in range(len(self.coords))}
vor_coors = []
index = {}
# Build a list of coordinates and one index per triangle/region
for tidx, (a, b, c) in enumerate(sorted(self.triangles)):
vor_coors.append(self.circles[(a, b, c)][0])
# Insert triangle, rotating it so the key is the "last" vertex
useVertex[a] += [(b, c, a)]
useVertex[b] += [(c, a, b)]
useVertex[c] += [(a, b, c)]
# Set tidx as the index to use with this triangle
index[(a, b, c)] = tidx
index[(c, a, b)] = tidx
index[(b, c, a)] = tidx
# init regions per coordinate dictionary
regions = {}
# Sort each region in a coherent order, and substitude each triangle
# by its index
for i in range(4, len(self.coords)):
v = useVertex[i][0][0] # Get a vertex of a triangle
r = []
for _ in range(len(useVertex[i])):
# Search the triangle beginning with vertex v
t = [t for t in useVertex[i] if t[0] == v][0]
r.append(index[t]) # Add the index of this triangle to region
v = t[1] # Choose the next vertex to search
regions[i-4] = r # Store region.
return vor_coors, regions