source: anuga_core/source/anuga/alpha_shape/alpha_shape.py @ 7340

Last change on this file since 7340 was 7317, checked in by rwilson, 16 years ago

Replaced 'print' statements with log.critical() calls.

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