1 | """Alpha shape |
---|
2 | Determine the shape of a set of points. |
---|
3 | |
---|
4 | From website by Kaspar Fischer: |
---|
5 | As mentionned in Edelsbrunner's and Muecke's paper, one can |
---|
6 | intuitively think of an alpha-shape as the following: |
---|
7 | |
---|
8 | Imagine a huge mass of ice-cream making up the space and containing |
---|
9 | the points S as ``hard'' chocolate pieces. Using one of these |
---|
10 | sphere-formed ice-cream spoons we carve out all parts of the ice-cream |
---|
11 | block we can reach without bumping into chocolate pieces, even |
---|
12 | carving out holes in the inside (eg. parts not reachable by simply |
---|
13 | moving the spoon from the outside). We will eventually end up with a |
---|
14 | (not necessarily convex) object bounded by caps, arcs and points. If |
---|
15 | we now straighten all ``round'' faces to triangles and line segments, |
---|
16 | we have an intuitive description of what is called the alpha-shape. |
---|
17 | |
---|
18 | Author: Vanessa Robins, ANU |
---|
19 | """ |
---|
20 | |
---|
21 | import exceptions |
---|
22 | from Numeric import array, Float, divide_safe, sqrt, product |
---|
23 | from load_mesh.loadASCII import import_points_file, export_boundary_file |
---|
24 | import random |
---|
25 | |
---|
26 | class AlphaError(exceptions.Exception):pass |
---|
27 | class PointError(AlphaError): pass |
---|
28 | class FlagError(AlphaError): pass |
---|
29 | |
---|
30 | OUTPUT_FILE_TITLE = "# The alpha shape boundary defined by point index pairs of edges" |
---|
31 | INF = pow(10,9) |
---|
32 | EPSILON = 1.0e-12 |
---|
33 | |
---|
34 | def alpha_shape_via_files(point_file, boundary_file, alpha= None): |
---|
35 | """ |
---|
36 | Load a point file and return the alpha shape boundary as a boundary file. |
---|
37 | |
---|
38 | Inputs: |
---|
39 | point_file: File location of the input file, points format (.xya or .pts) |
---|
40 | boundary_file: File location of the generated output file |
---|
41 | alpha: The alpha value can be optionally specified. If it is not specified |
---|
42 | the optimum alpha value will be used. |
---|
43 | """ |
---|
44 | point_dict = import_points_file(point_file) |
---|
45 | points = point_dict['pointlist'] |
---|
46 | |
---|
47 | AS = Alpha_Shape(points, alpha) |
---|
48 | AS.write_boundary(boundary_file) |
---|
49 | |
---|
50 | |
---|
51 | class Alpha_Shape: |
---|
52 | |
---|
53 | def __init__(self, points, alpha = None): |
---|
54 | """ |
---|
55 | An Alpha_Shape requires input of a set of points. Other class routines |
---|
56 | return the alpha shape boundary. |
---|
57 | |
---|
58 | Inputs: |
---|
59 | points: List of coordinate pairs [[x1, y1],[x2, y2]..] |
---|
60 | alpha: The alpha value can be optionally specified. If it is |
---|
61 | not specified the optimum alpha value will be used. |
---|
62 | """ |
---|
63 | self._set_points(points) |
---|
64 | self._alpha_shape_algorithm(alpha) |
---|
65 | |
---|
66 | def _set_points(self, points): |
---|
67 | """ |
---|
68 | Create self.points array, do Error checking |
---|
69 | Inputs: |
---|
70 | points: List of coordinate pairs [[x1, y1],[x2, y2]..] |
---|
71 | """ |
---|
72 | # print "setting points" |
---|
73 | if len (points) <= 2: |
---|
74 | raise PointError, "Too few points to find an alpha shape" |
---|
75 | if len(points)==3: |
---|
76 | #check not in a straight line |
---|
77 | # FIXME check points 1,2,3 if straingt, check if points 2,3,4, ect |
---|
78 | x01 = points[0][0] - points[1][0] |
---|
79 | y01 = points[0][1] - points[1][1] |
---|
80 | x12 = points[1][0] - points[2][0] |
---|
81 | y12 = points[1][1] - points[2][1] |
---|
82 | crossprod = x01*y12 - x12*y01 |
---|
83 | if crossprod==0: |
---|
84 | raise PointError, "Three points on a straight line" |
---|
85 | |
---|
86 | #Convert input to Numeric arrays |
---|
87 | self.points = array(points).astype(Float) |
---|
88 | |
---|
89 | |
---|
90 | def write_boundary(self,file_name): |
---|
91 | """ |
---|
92 | Write the boundary to a file |
---|
93 | """ |
---|
94 | #print " this info will be in the file" |
---|
95 | export_boundary_file(file_name, self.get_boundary(), |
---|
96 | OUTPUT_FILE_TITLE, delimiter = ',') |
---|
97 | |
---|
98 | def get_boundary(self): |
---|
99 | """ |
---|
100 | Return a list of tuples. |
---|
101 | Each tuple represents a segment in the boundary |
---|
102 | by the index of its two end points. |
---|
103 | The list of tuples represents the alpha shape boundary. |
---|
104 | """ |
---|
105 | return self.boundary |
---|
106 | |
---|
107 | def set_boundary_type(self,raw_boundary=True, |
---|
108 | remove_holes=False, |
---|
109 | smooth_indents=False, |
---|
110 | expand_pinch=False, |
---|
111 | boundary_points_fraction=0.2): |
---|
112 | """ |
---|
113 | Use the flags to set constraints on the boundary: |
---|
114 | raw_boundary Return raw boundary i.e. the regular edges of the |
---|
115 | alpha shape. |
---|
116 | remove_holes filter to remove small holes |
---|
117 | (small is defined by boundary_points_fraction ) |
---|
118 | smooth_indents remove sharp triangular indents in boundary |
---|
119 | expand_pinch test for pinch-off and correct |
---|
120 | i.e. a boundary vertex with more than two edges. |
---|
121 | """ |
---|
122 | |
---|
123 | if raw_boundary: |
---|
124 | # reset alpha shape boundary |
---|
125 | reg_edge = self.get_regular_edges(self.alpha) |
---|
126 | self.boundary = [self.edge[k] for k in reg_edge] |
---|
127 | self._init_boundary_triangles() |
---|
128 | if remove_holes: |
---|
129 | #remove small holes |
---|
130 | self.boundary = self._remove_holes(boundary_points_fraction) |
---|
131 | if smooth_indents: |
---|
132 | #remove sharp triangular indents |
---|
133 | self.boundary = self._smooth_indents() |
---|
134 | if expand_pinch: |
---|
135 | #deal with pinch-off |
---|
136 | self.boundary = self._expand_pinch() |
---|
137 | |
---|
138 | |
---|
139 | def get_delaunay(self): |
---|
140 | """ |
---|
141 | """ |
---|
142 | return self.deltri |
---|
143 | |
---|
144 | def get_optimum_alpha(self): |
---|
145 | """ |
---|
146 | """ |
---|
147 | return self.optimum_alpha |
---|
148 | |
---|
149 | def get_alpha(self): |
---|
150 | """ |
---|
151 | Return current alpha value. |
---|
152 | """ |
---|
153 | return self.alpha |
---|
154 | |
---|
155 | def set_alpha(self,alpha): |
---|
156 | """ |
---|
157 | Set alpha and update alpha-boundary. |
---|
158 | """ |
---|
159 | self.alpha = alpha |
---|
160 | reg_edge = self.get_regular_edges(alpha) |
---|
161 | self.boundary = [self.edge[k] for k in reg_edge] |
---|
162 | self._init_boundary_triangles() |
---|
163 | |
---|
164 | |
---|
165 | def _alpha_shape_algorithm(self, alpha=None): |
---|
166 | """ |
---|
167 | Given a set of points (self.points) and an optional alpha value |
---|
168 | determines the alpha shape boundary (stored in self.boundary, |
---|
169 | accessed by get_boundary). |
---|
170 | |
---|
171 | Inputs: |
---|
172 | alpha: The alpha value can be optionally specified. If it is |
---|
173 | not specified the optimum alpha value will be used. |
---|
174 | """ |
---|
175 | |
---|
176 | #print "starting alpha shape algorithm" |
---|
177 | |
---|
178 | self.alpha = alpha |
---|
179 | |
---|
180 | ## Build Delaunay triangulation |
---|
181 | import triang |
---|
182 | points = [] |
---|
183 | seglist = [] |
---|
184 | holelist = [] |
---|
185 | regionlist = [] |
---|
186 | pointattlist = [] |
---|
187 | segattlist = [] |
---|
188 | trilist = [] |
---|
189 | |
---|
190 | points = [(pt[0], pt[1]) for pt in self.points] |
---|
191 | pointattlist = [ [] for i in range(len(points)) ] |
---|
192 | mode = "Qzcn" |
---|
193 | #print "computing delaunay triangulation ... \n" |
---|
194 | tridata = triang.genMesh(points,seglist,holelist,regionlist, |
---|
195 | pointattlist,segattlist,trilist,mode) |
---|
196 | #print tridata |
---|
197 | #print "point attlist: ", tridata['generatedpointattributelist'],"\n" |
---|
198 | #print "hull segments: ", tridata['generatedsegmentlist'], "\n" |
---|
199 | self.deltri = tridata['generatedtrianglelist'] |
---|
200 | self.deltrinbr = tridata['generatedtriangleneighborlist'] |
---|
201 | self.hulledges = tridata['generatedsegmentlist'] |
---|
202 | |
---|
203 | #print "Number of delaunay triangles: ", len(self.deltri), "\n" |
---|
204 | #print "deltrinbrs: ", self.deltrinbr, "\n" |
---|
205 | |
---|
206 | ## Build Alpha table |
---|
207 | ## the following routines determine alpha thresholds for the |
---|
208 | ## triangles, edges, and vertices of the delaunay triangulation |
---|
209 | self._tri_circumradius() |
---|
210 | # print "Largest circumradius ", max(self.triradius) |
---|
211 | self._edge_intervals() |
---|
212 | self._vertex_intervals() |
---|
213 | |
---|
214 | if alpha==None: |
---|
215 | # Find optimum alpha |
---|
216 | # Ken Clarkson's hull program uses smallest alpha so that |
---|
217 | # every vertex is non-singular so... |
---|
218 | self.optimum_alpha = max([iv[0] for iv in self.vertexinterval \ |
---|
219 | if iv!=[] ]) |
---|
220 | # print "optimum alpha ", self.optimum_alpha |
---|
221 | alpha = self.optimum_alpha |
---|
222 | self.alpha = alpha |
---|
223 | reg_edge = self.get_regular_edges(self.alpha) |
---|
224 | self.boundary = [self.edge[k] for k in reg_edge] |
---|
225 | #print "alpha boundary edges ", self.boundary |
---|
226 | self._init_boundary_triangles() |
---|
227 | return |
---|
228 | |
---|
229 | def _tri_circumradius(self): |
---|
230 | """ |
---|
231 | Compute circumradii of the delaunay triangles |
---|
232 | """ |
---|
233 | |
---|
234 | x = self.points[:,0] |
---|
235 | y = self.points[:,1] |
---|
236 | ind1 = [self.deltri[j][0] for j in range(len(self.deltri))] |
---|
237 | ind2 = [self.deltri[j][1] for j in range(len(self.deltri))] |
---|
238 | ind3 = [self.deltri[j][2] for j in range(len(self.deltri))] |
---|
239 | |
---|
240 | x1 = array([ x[j] for j in ind1 ]) |
---|
241 | y1 = array([ y[j] for j in ind1 ]) |
---|
242 | x2 = array([ x[j] for j in ind2 ]) |
---|
243 | y2 = array([ y[j] for j in ind2 ]) |
---|
244 | x3 = array([ x[j] for j in ind3 ]) |
---|
245 | y3 = array([ y[j] for j in ind3 ]) |
---|
246 | |
---|
247 | x21 = x2-x1 |
---|
248 | x31 = x3-x1 |
---|
249 | y21 = y2-y1 |
---|
250 | y31 = y3-y1 |
---|
251 | |
---|
252 | dist21 = x21*x21 + y21*y21 |
---|
253 | dist31 = x31*x31 + y31*y31 |
---|
254 | |
---|
255 | denom = x21*y31 - x31*y21 |
---|
256 | #print "denom = ", denom |
---|
257 | |
---|
258 | # dx/2, dy/2 give circumcenter relative to x1,y1. |
---|
259 | # dx = (y31*dist21 - y21*dist31)/denom |
---|
260 | # dy = (x21*dist31 - x31*dist21)/denom |
---|
261 | # first need to check for near-zero values of denom |
---|
262 | delta = 0.00000001 |
---|
263 | zeroind = [k for k in range(len(denom)) if \ |
---|
264 | (denom[k]< EPSILON and denom[k] > -EPSILON)] |
---|
265 | # if some denom values are close to zero, |
---|
266 | # we perturb the associated vertices and recalculate |
---|
267 | while zeroind!=[]: |
---|
268 | random.seed() |
---|
269 | print "Warning: degenerate triangles found in alpha_shape.py, results may be inaccurate." |
---|
270 | for d in zeroind: |
---|
271 | x1[d] = x1[d]+delta*(random.random()-0.5) |
---|
272 | x2[d] = x2[d]+delta*(random.random()-0.5) |
---|
273 | x3[d] = x3[d]+delta*(random.random()-0.5) |
---|
274 | y1[d] = y1[d]+delta*(random.random()-0.5) |
---|
275 | y2[d] = y2[d]+delta*(random.random()-0.5) |
---|
276 | y3[d] = y3[d]+delta*(random.random()-0.5) |
---|
277 | x21 = x2-x1 |
---|
278 | x31 = x3-x1 |
---|
279 | y21 = y2-y1 |
---|
280 | y31 = y3-y1 |
---|
281 | dist21 = x21*x21 + y21*y21 |
---|
282 | dist31 = x31*x31 + y31*y31 |
---|
283 | denom = x21*y31 - x31*y21 |
---|
284 | zeroind = [k for k in range(len(denom)) if \ |
---|
285 | (denom[k]< EPSILON and denom[k] > -EPSILON)] |
---|
286 | try: |
---|
287 | dx = divide_safe(y31*dist21 - y21*dist31,denom) |
---|
288 | dy = divide_safe(x21*dist31 - x31*dist21,denom) |
---|
289 | except ZeroDivisionError: |
---|
290 | raise AlphaError |
---|
291 | |
---|
292 | self.triradius = 0.5*sqrt(dx*dx + dy*dy) |
---|
293 | #print "triangle radii", self.triradius |
---|
294 | |
---|
295 | def _edge_intervals(self): |
---|
296 | """ |
---|
297 | for each edge, find triples |
---|
298 | (length/2, min_adj_triradius, max_adj_triradius) if unattached |
---|
299 | (min_adj_triradius, min_adj_triradius, max_adj_triradius) if attached. |
---|
300 | An edge is attached if it is opposite an obtuse angle |
---|
301 | """ |
---|
302 | |
---|
303 | # It should be possible to rewrite this routine in an array-friendly |
---|
304 | # form like _tri_circumradius() if we need to speed things up. |
---|
305 | # Hard to do though. |
---|
306 | |
---|
307 | edges = [] |
---|
308 | edgenbrs = [] |
---|
309 | edgeinterval = [] |
---|
310 | for t in range(len(self.deltri)): |
---|
311 | tri = self.deltri[t] |
---|
312 | trinbr = self.deltrinbr[t] |
---|
313 | dx = array([self.points[tri[(i+1)%3],0] - |
---|
314 | self.points[tri[(i+2)%3],0] for i in [0,1,2]]) |
---|
315 | dy = array([self.points[tri[(i+1)%3],1] - |
---|
316 | self.points[tri[(i+2)%3],1] for i in [0,1,2]]) |
---|
317 | elen = sqrt(dx*dx+dy*dy) |
---|
318 | # really only need sign - not angle value: |
---|
319 | anglesign = array([(-dx[(i+1)%3]*dx[(i+2)%3]- |
---|
320 | dy[(i+1)%3]*dy[(i+2)%3]) for i in [0,1,2]]) |
---|
321 | |
---|
322 | #print "dx ",dx,"\n" |
---|
323 | #print "dy ",dy,"\n" |
---|
324 | #print "edge lengths of triangle ",t,"\t",elen,"\n" |
---|
325 | #print "angles ",angle,"\n" |
---|
326 | |
---|
327 | for i in [0,1,2]: |
---|
328 | j = (i+1)%3 |
---|
329 | k = (i+2)%3 |
---|
330 | if trinbr[i]==-1: |
---|
331 | edges.append((tri[j], tri[k])) |
---|
332 | edgenbrs.append((t, -1)) |
---|
333 | edgeinterval.append([0.5*elen[i], self.triradius[t], INF]) |
---|
334 | elif (tri[j]<tri[k]): |
---|
335 | edges.append((tri[j], tri[k])) |
---|
336 | edgenbrs.append((t, trinbr[i])) |
---|
337 | edgeinterval.append([0.5*elen[i],\ |
---|
338 | min(self.triradius[t],self.triradius[trinbr[i]]),\ |
---|
339 | max(self.triradius[t],self.triradius[trinbr[i]]) ]) |
---|
340 | else: |
---|
341 | continue |
---|
342 | if anglesign[i] < 0: |
---|
343 | edgeinterval[-1][0] = edgeinterval[-1][1] |
---|
344 | |
---|
345 | self.edge = edges |
---|
346 | self.edgenbr = edgenbrs |
---|
347 | self.edgeinterval = edgeinterval |
---|
348 | #print "edges: ",edges, "\n" |
---|
349 | #print "edge nbrs:", edgenbrs ,"\n" |
---|
350 | #print "edge intervals: ",edgeinterval , "\n" |
---|
351 | |
---|
352 | def _vertex_intervals(self): |
---|
353 | """ |
---|
354 | for each vertex find pairs |
---|
355 | (min_adj_triradius, max_adj_triradius) |
---|
356 | """ |
---|
357 | nv = len(self.points) |
---|
358 | vertexnbrs = [ [] for i in range(nv)] |
---|
359 | vertexinterval = [ [] for i in range(nv)] |
---|
360 | for t in range(len(self.deltri)): |
---|
361 | for j in self.deltri[t]: |
---|
362 | vertexnbrs[j].append(t) |
---|
363 | for h in range(len(self.hulledges)): |
---|
364 | for j in self.hulledges[h]: |
---|
365 | vertexnbrs[j].append(-1) |
---|
366 | |
---|
367 | for i in range(nv): |
---|
368 | radii = [ self.triradius[t] for t in vertexnbrs[i] if t>=0 ] |
---|
369 | try: |
---|
370 | vertexinterval[i] = [min(radii), max(radii)] |
---|
371 | if vertexnbrs[i][-1]==-1: |
---|
372 | vertexinterval[i][1]=INF |
---|
373 | except ValueError: |
---|
374 | raise AlphaError |
---|
375 | |
---|
376 | self.vertexnbr = vertexnbrs |
---|
377 | self.vertexinterval = vertexinterval |
---|
378 | #print "vertex nbrs ", vertexnbrs, "\n" |
---|
379 | #print "vertex intervals ",vertexinterval, "\n" |
---|
380 | |
---|
381 | def get_alpha_triangles(self,alpha): |
---|
382 | """ |
---|
383 | Given an alpha value, |
---|
384 | return indices of triangles in the alpha-shape |
---|
385 | """ |
---|
386 | def tri_rad_lta(k): |
---|
387 | return self.triradius[k]<=alpha |
---|
388 | |
---|
389 | return filter(tri_rad_lta, range(len(self.triradius))) |
---|
390 | |
---|
391 | def get_regular_edges(self,alpha): |
---|
392 | """ |
---|
393 | Given an alpha value, |
---|
394 | return the indices of edges on the boundary of the alpha-shape |
---|
395 | """ |
---|
396 | def reg_edge(k): |
---|
397 | return self.edgeinterval[k][1]<=alpha and \ |
---|
398 | self.edgeinterval[k][2]>alpha |
---|
399 | |
---|
400 | return filter(reg_edge, range(len(self.edgeinterval))) |
---|
401 | |
---|
402 | def get_exposed_vertices(self,alpha): |
---|
403 | """ |
---|
404 | Given an alpha value, |
---|
405 | return the vertices on the boundary of the alpha-shape |
---|
406 | """ |
---|
407 | def exp_vert(k): |
---|
408 | return self.vertexinterval[k][0]<=alpha and \ |
---|
409 | self.vertexinterval[k][1]>alpha |
---|
410 | |
---|
411 | return filter(exp_vert, range(len(self.vertexinterval))) |
---|
412 | |
---|
413 | def _vertices_from_edges(self,elist): |
---|
414 | """ |
---|
415 | Returns the list of unique vertex labels from edges in elist |
---|
416 | """ |
---|
417 | |
---|
418 | v1 = [elist[k][0] for k in range(len(elist))] |
---|
419 | v2 = [elist[k][1] for k in range(len(elist))] |
---|
420 | v = v1+v2 |
---|
421 | v.sort() |
---|
422 | vertices = [v[k] for k in range(len(v)) if v[k]!=v[k-1] ] |
---|
423 | return vertices |
---|
424 | |
---|
425 | def _init_boundary_triangles(self): |
---|
426 | """ |
---|
427 | Creates the initial list of triangle indices |
---|
428 | exterior to and touching the boundary of the alpha shape |
---|
429 | """ |
---|
430 | def tri_rad_gta(k): |
---|
431 | return self.triradius[k]>self.alpha |
---|
432 | |
---|
433 | extrind = filter(tri_rad_gta, range(len(self.triradius))) |
---|
434 | |
---|
435 | bv = self._vertices_from_edges(self.boundary) |
---|
436 | |
---|
437 | btri = [] |
---|
438 | for et in extrind: |
---|
439 | v0 = self.deltri[et][0] |
---|
440 | v1 = self.deltri[et][1] |
---|
441 | v2 = self.deltri[et][2] |
---|
442 | if v0 in bv or v1 in bv or v2 in bv: |
---|
443 | btri.append(et) |
---|
444 | |
---|
445 | self.boundarytriangle = btri |
---|
446 | |
---|
447 | #print "exterior triangles: ", extrind |
---|
448 | |
---|
449 | |
---|
450 | def _remove_holes(self,small): |
---|
451 | """ |
---|
452 | Given the edges in self.boundary, finds the largest components. |
---|
453 | The routine does this by implementing a union-find algorithm. |
---|
454 | """ |
---|
455 | |
---|
456 | #print "running _remove_holes \n" |
---|
457 | |
---|
458 | bdry = self.boundary |
---|
459 | |
---|
460 | def findroot(i): |
---|
461 | if vptr[i] < 0: |
---|
462 | return i |
---|
463 | k = findroot(vptr[i]) |
---|
464 | vptr[i] = k # this produces "path compression" in the |
---|
465 | # union-find tree. |
---|
466 | return k |
---|
467 | |
---|
468 | |
---|
469 | |
---|
470 | # get a list of unique vertex labels: |
---|
471 | verts = self._vertices_from_edges(bdry) |
---|
472 | #print "verts ", verts, "\n" |
---|
473 | |
---|
474 | # vptr represents the union-find tree. |
---|
475 | # if vptr[i] = EMPTY < 0, the vertex verts[i] has not been visited yet |
---|
476 | # if vptr[i] = j > 0, then j verts[j] is the parent of verts[i]. |
---|
477 | # if vptr[i] = n < 0, then verts[i] is a root vertex and |
---|
478 | # represents a connected component of n vertices. |
---|
479 | |
---|
480 | #initialise vptr to negative number outside range |
---|
481 | EMPTY = -max(verts)-len(verts) |
---|
482 | vptr = [EMPTY for k in range(len(verts))] |
---|
483 | #print "vptr init: ", vptr, "\n" |
---|
484 | |
---|
485 | #add edges and maintain union tree |
---|
486 | for i in range(len(bdry)): |
---|
487 | #print "edge ",i,"\t",bdry[i] |
---|
488 | vl = verts.index(bdry[i][0]) |
---|
489 | vr = verts.index(bdry[i][1]) |
---|
490 | rvl = findroot(vl) |
---|
491 | rvr = findroot(vr) |
---|
492 | #print "roots: ",rvl, rvr |
---|
493 | if not(rvl==rvr): |
---|
494 | if (vptr[vl]==EMPTY): |
---|
495 | if (vptr[vr]==EMPTY): |
---|
496 | vptr[vl] = -2 |
---|
497 | vptr[vr] = vl |
---|
498 | else: |
---|
499 | vptr[vl] = rvr |
---|
500 | vptr[rvr] = vptr[rvr]-1 |
---|
501 | else: |
---|
502 | if (vptr[vr]==EMPTY): |
---|
503 | vptr[vr] = rvl |
---|
504 | vptr[rvl] = vptr[rvl]-1 |
---|
505 | else: |
---|
506 | if vptr[rvl] > vptr[rvr]: |
---|
507 | vptr[rvr] = vptr[rvr] + vptr[rvl] |
---|
508 | vptr[rvl] = rvr |
---|
509 | vptr[vl] = rvr |
---|
510 | else: |
---|
511 | vptr[rvl] = vptr[rvl] + vptr[rvr] |
---|
512 | vptr[rvr] = rvl |
---|
513 | vptr[vr] = rvl |
---|
514 | #print "vptr: ", vptr, "\n" |
---|
515 | # end edge loop |
---|
516 | |
---|
517 | if vptr.count(EMPTY): |
---|
518 | raise FlagError, "We didn't hit all the vertices in the boundary" |
---|
519 | |
---|
520 | # print "number of vertices in the connected components: ", [-v for v in vptr if v<0], "\n" |
---|
521 | # print "largest component has: ", -min(vptr), " points. \n" |
---|
522 | # discard the edges in the little components |
---|
523 | # (i.e. those components with less than 'small' fraction of bdry points) |
---|
524 | cutoff = round(small*len(verts)) |
---|
525 | # print "cutoff component size is ", cutoff, "\n" |
---|
526 | largest_component = -min(vptr) |
---|
527 | if cutoff > largest_component: |
---|
528 | cutoff = round((1-small)*largest_component) |
---|
529 | |
---|
530 | # littleind has root indices for small components |
---|
531 | littleind = [k for k in range(len(vptr)) if \ |
---|
532 | (vptr[k]<0 and vptr[k]>-cutoff)] |
---|
533 | if littleind: |
---|
534 | # littlecomp has all vptr indices in the small components |
---|
535 | littlecomp = [k for k in range(len(vptr)) if \ |
---|
536 | findroot(k) in littleind] |
---|
537 | # vdiscard has the vertex indices corresponding to vptr indices |
---|
538 | vdiscard = [verts[k] for k in littlecomp] |
---|
539 | newbdry = [e for e in bdry if \ |
---|
540 | not((e[0] in vdiscard) and (e[1] in vdiscard))] |
---|
541 | |
---|
542 | newverts = self._vertices_from_edges(newbdry) |
---|
543 | # recompute the boundary triangles |
---|
544 | newbt = [] |
---|
545 | for bt in self.boundarytriangle: |
---|
546 | v0 = self.deltri[bt][0] |
---|
547 | v1 = self.deltri[bt][1] |
---|
548 | v2 = self.deltri[bt][2] |
---|
549 | if (v0 in newverts or v1 in newverts or v2 in newverts): |
---|
550 | newbt.append(bt) |
---|
551 | |
---|
552 | self.boundarytriangle = newbt |
---|
553 | else: |
---|
554 | newbdry = bdry |
---|
555 | |
---|
556 | return newbdry |
---|
557 | |
---|
558 | |
---|
559 | def _smooth_indents(self): |
---|
560 | """ |
---|
561 | Given edges in bdry, test for acute-angle triangular indents |
---|
562 | and remove them. |
---|
563 | """ |
---|
564 | |
---|
565 | #print "running _smooth_indents \n" |
---|
566 | |
---|
567 | bdry = self.boundary |
---|
568 | bdrytri = self.boundarytriangle |
---|
569 | |
---|
570 | # find boundary triangles that have two edges in bdry |
---|
571 | # v2ind has the place index relative to the triangle deltri[ind] |
---|
572 | # for the bdry vertex where the two edges meet. |
---|
573 | |
---|
574 | verts = self._vertices_from_edges(bdry) |
---|
575 | |
---|
576 | b2etri = [] |
---|
577 | for ind in bdrytri: |
---|
578 | bect = 0 |
---|
579 | v2ind = [0,1,2] |
---|
580 | for j in [0,1,2]: |
---|
581 | eda = (self.deltri[ind][(j+1)%3], self.deltri[ind][(j+2)%3]) |
---|
582 | edb = (self.deltri[ind][(j+2)%3], self.deltri[ind][(j+1)%3]) |
---|
583 | if eda in bdry or edb in bdry: |
---|
584 | bect +=1 |
---|
585 | v2ind.remove(j) |
---|
586 | if bect==2: |
---|
587 | b2etri.append((ind,v2ind[0])) |
---|
588 | |
---|
589 | # test the bdrytri triangles for acute angles |
---|
590 | acutetri = [] |
---|
591 | for tind in b2etri: |
---|
592 | tri = self.deltri[tind[0]] |
---|
593 | |
---|
594 | dx = array([self.points[tri[(i+1)%3],0] - \ |
---|
595 | self.points[tri[(i+2)%3],0] for i in [0,1,2]]) |
---|
596 | dy = array([self.points[tri[(i+1)%3],1] - \ |
---|
597 | self.points[tri[(i+2)%3],1] for i in [0,1,2]]) |
---|
598 | anglesign = array([(-dx[(i+1)%3]*dx[(i+2)%3]-\ |
---|
599 | dy[(i+1)%3]*dy[(i+2)%3]) for i in [0,1,2]]) |
---|
600 | # record any triangle that has an acute angle spanned by |
---|
601 | #two edges along the boundary.. |
---|
602 | if anglesign[tind[1]] > 0: |
---|
603 | acutetri.append(tind[0]) |
---|
604 | |
---|
605 | #print "acute boundary triangles ", acutetri |
---|
606 | |
---|
607 | # adjust the bdry edges and triangles by adding |
---|
608 | #in the acutetri triangles |
---|
609 | for pind in acutetri: |
---|
610 | bdrytri.remove(pind) |
---|
611 | tri = self.deltri[pind] |
---|
612 | for i in [0,1,2]: |
---|
613 | bdry.append((tri[(i+1)%3], tri[(i+2)%3])) |
---|
614 | |
---|
615 | newbdry = [] |
---|
616 | for ed in bdry: |
---|
617 | numed = bdry.count(ed)+bdry.count((ed[1],ed[0])) |
---|
618 | if numed%2 == 1: |
---|
619 | newbdry.append(ed) |
---|
620 | |
---|
621 | #print "new boundary ", newbdry |
---|
622 | return newbdry |
---|
623 | |
---|
624 | def _expand_pinch(self): |
---|
625 | """ |
---|
626 | Given edges in bdry, test for vertices with more than 2 incident edges. |
---|
627 | Expand by adding back in associated triangles. |
---|
628 | """ |
---|
629 | #print "running _expand_pinch \n" |
---|
630 | |
---|
631 | bdry = self.boundary |
---|
632 | bdrytri = self.boundarytriangle |
---|
633 | |
---|
634 | v1 = [bdry[k][0] for k in range(len(bdry))] |
---|
635 | v2 = [bdry[k][1] for k in range(len(bdry))] |
---|
636 | v = v1+v2 |
---|
637 | v.sort() |
---|
638 | probv = [v[k] for k in range(len(v)) \ |
---|
639 | if (v[k]!=v[k-1] and v.count(v[k])>2) ] |
---|
640 | #print "problem vertices: ", probv |
---|
641 | |
---|
642 | # find boundary triangles that have at least one vertex in probv |
---|
643 | probtri = [] |
---|
644 | for ind in bdrytri: |
---|
645 | v0 = self.deltri[ind][0] |
---|
646 | v1 = self.deltri[ind][1] |
---|
647 | v2 = self.deltri[ind][2] |
---|
648 | if v0 in probv or v1 in probv or v2 in probv: |
---|
649 | probtri.append(ind) |
---|
650 | |
---|
651 | #print "problem boundary triangle indices ", probtri |
---|
652 | |
---|
653 | # "add in" the problem triangles |
---|
654 | for pind in probtri: |
---|
655 | bdrytri.remove(pind) |
---|
656 | tri = self.deltri[pind] |
---|
657 | for i in [0,1,2]: |
---|
658 | bdry.append((tri[(i+1)%3], tri[(i+2)%3])) |
---|
659 | |
---|
660 | newbdry = [] |
---|
661 | for ed in bdry: |
---|
662 | numed = bdry.count(ed)+bdry.count((ed[1],ed[0])) |
---|
663 | if numed%2 == 1: |
---|
664 | newbdry.append(ed) |
---|
665 | |
---|
666 | #print "new boundary ", newbdry |
---|
667 | return newbdry |
---|
668 | |
---|
669 | |
---|
670 | #------------------------------------------------------------- |
---|
671 | if __name__ == "__main__": |
---|
672 | """ |
---|
673 | Load in a data point file. |
---|
674 | Determine the alpha shape boundary |
---|
675 | Save the boundary to a file. |
---|
676 | |
---|
677 | usage: alpha_shape.py point_file.xya boundary_file.bnd [alpha] |
---|
678 | |
---|
679 | The alpha value is optional. |
---|
680 | """ |
---|
681 | |
---|
682 | import os, sys |
---|
683 | usage = "usage: %s point_file.xya boundary_file.bnd [alpha]"%os.path.basename(sys.argv[0]) |
---|
684 | if len(sys.argv) < 3: |
---|
685 | print usage |
---|
686 | else: |
---|
687 | point_file = sys.argv[1] |
---|
688 | boundary_file = sys.argv[2] |
---|
689 | if len(sys.argv) > 4: |
---|
690 | alpha = sys.argv[3] |
---|
691 | else: |
---|
692 | alpha = None |
---|
693 | |
---|
694 | #print "about to call alpha shape routine \n" |
---|
695 | alpha_shape_via_files(point_file, boundary_file, alpha) |
---|
696 | |
---|