source: inundation/alpha_shape/alpha_shape.py @ 3330

Last change on this file since 3330 was 3027, checked in by duncan, 18 years ago

directory name change to avoid name conflicts.

File size: 25.3 KB
RevLine 
[1420]1"""Alpha shape
2Determine the shape of a set of points.
3
[1431]4From website by Kaspar Fischer:
[1433]5As mentionned in Edelsbrunner's and Muecke's paper, one can
[1420]6intuitively think of an alpha-shape as the following:
[1431]7
[1420]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
[1432]11block we can reach without bumping into chocolate pieces, even
[1420]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.
[1435]17
18Author: Vanessa Robins, ANU
[674]19"""
20
21import exceptions
[959]22from Numeric import array, Float, divide_safe, sqrt, product
[1381]23from load_mesh.loadASCII import import_points_file, export_boundary_file
[691]24import random
[674]25
[1431]26class AlphaError(exceptions.Exception):pass
27class PointError(AlphaError): pass
28class FlagError(AlphaError): pass
[674]29
30OUTPUT_FILE_TITLE = "# The alpha shape boundary defined by point index pairs of edges"
31INF = pow(10,9)
[1435]32EPSILON = 1.0e-12
[674]33
[959]34def alpha_shape_via_files(point_file, boundary_file, alpha= None):
[1420]35    """
36    Load a point file and return the alpha shape boundary as a boundary file.
[674]37   
[1420]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.
[1432]43    """
[1381]44    point_dict = import_points_file(point_file)
[674]45    points = point_dict['pointlist']
46   
[1432]47    AS = Alpha_Shape(points, alpha)
48    AS.write_boundary(boundary_file)
49   
[674]50
51class Alpha_Shape:
52
53    def __init__(self, points, alpha = None):
[1420]54        """
[1435]55        An Alpha_Shape requires input of a set of points. Other class routines
[1432]56        return the alpha shape boundary.
[674]57       
[1420]58        Inputs:       
[674]59          points: List of coordinate pairs [[x1, y1],[x2, y2]..]
[1420]60          alpha: The alpha value can be optionally specified.  If it is
61          not specified the optimum alpha value will be used.
[674]62        """
63        self._set_points(points)
[1432]64        self._alpha_shape_algorithm(alpha)
[674]65
66    def _set_points(self, points):
67        """
[1420]68        Create self.points array, do Error checking
69        Inputs:       
70          points: List of coordinate pairs [[x1, y1],[x2, y2]..]
[674]71        """
72        # print "setting points"
73        if len (points) <= 2:
74            raise PointError, "Too few points to find an alpha shape"
[691]75        if len(points)==3:
76            #check not in a straight line
[1433]77            # FIXME check points 1,2,3 if straingt, check if points 2,3,4, ect
[691]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"
[674]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        """
[1432]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.
[674]104        """
105        return self.boundary
106
[987]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):
[956]112        """
[1420]113        Use the flags to set constraints on the boundary:
[1435]114        raw_boundary    Return raw boundary i.e. the regular edges of the
115                        alpha shape.
[1432]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.
[956]121        """
122
[987]123        if raw_boundary:
[1432]124            # reset alpha shape boundary
[987]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:
[1432]129            #remove small holes
[987]130            self.boundary = self._remove_holes(boundary_points_fraction)
131        if smooth_indents:
[1432]132            #remove sharp triangular indents
[987]133            self.boundary = self._smooth_indents()
134        if expand_pinch:
[956]135            #deal with pinch-off
[987]136            self.boundary = self._expand_pinch()
[956]137       
138
[674]139    def get_delaunay(self):
140        """
141        """
142        return self.deltri
[691]143
144    def get_optimum_alpha(self):
145        """
146        """
147        return self.optimum_alpha
148
149    def get_alpha(self):
150        """
[1420]151        Return current alpha value.
[691]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]
[987]162        self._init_boundary_triangles()
[674]163   
164           
[1420]165    def _alpha_shape_algorithm(self, alpha=None):
166        """
167        Given a set of points (self.points) and an optional alpha value
[1432]168        determines the alpha shape boundary (stored in self.boundary,
[1420]169        accessed by get_boundary).
[674]170       
[1420]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        """
[674]175
176        #print "starting alpha shape algorithm"
177
[691]178        self.alpha = alpha
179
[1432]180        ## Build Delaunay triangulation
[3027]181        import mesh_engine.triang as triang
[674]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"
[691]193        #print "computing delaunay triangulation ... \n"
[1435]194        tridata = triang.genMesh(points,seglist,holelist,regionlist,
195                                 pointattlist,segattlist,trilist,mode)
[674]196        #print tridata
[1435]197        #print "point attlist: ", tridata['generatedpointattributelist'],"\n"
[674]198        #print "hull segments: ", tridata['generatedsegmentlist'], "\n"
199        self.deltri = tridata['generatedtrianglelist']
200        self.deltrinbr = tridata['generatedtriangleneighborlist']
[1435]201        self.hulledges = tridata['generatedsegmentlist'] 
[674]202
[691]203        #print "Number of delaunay triangles: ", len(self.deltri), "\n"
[674]204        #print "deltrinbrs: ", self.deltrinbr, "\n"
205
[1432]206        ## Build Alpha table
207        ## the following routines determine alpha thresholds for the
208        ## triangles, edges, and vertices of the delaunay triangulation
[674]209        self._tri_circumradius()
[691]210        # print "Largest circumradius ", max(self.triradius)
[674]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...
[1435]218            self.optimum_alpha = max([iv[0] for iv in self.vertexinterval \
219                                      if iv!=[] ])
[691]220            # print "optimum alpha ", self.optimum_alpha
[1435]221            alpha = self.optimum_alpha
222        self.alpha = alpha
[1433]223        reg_edge = self.get_regular_edges(self.alpha)
[674]224        self.boundary = [self.edge[k] for k in reg_edge]
225        #print "alpha boundary edges ", self.boundary
[1435]226        self._init_boundary_triangles()           
[674]227        return
228
229    def _tri_circumradius(self):
[1420]230        """
231        Compute circumradii of the delaunay triangles
232        """
233       
[674]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
[1432]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
[691]262        delta = 0.00000001
[1435]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
[691]267        while zeroind!=[]:
[1433]268            random.seed()
[691]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
[1435]284            zeroind = [k for k in range(len(denom)) if \
285                       (denom[k]< EPSILON and  denom[k] > -EPSILON)]
[674]286        try:
287            dx = divide_safe(y31*dist21 - y21*dist31,denom)
288            dy = divide_safe(x21*dist31 - x31*dist21,denom)
289        except ZeroDivisionError:
[1433]290            raise  AlphaError
[674]291           
292        self.triradius = 0.5*sqrt(dx*dx + dy*dy)
293        #print "triangle radii", self.triradius
294
295    def _edge_intervals(self):
[1420]296        """
297         for each edge, find triples
298         (length/2, min_adj_triradius, max_adj_triradius) if unattached
[1432]299         (min_adj_triradius, min_adj_triradius, max_adj_triradius) if attached.
300         An edge is attached if it is opposite an obtuse angle
[1420]301        """
302       
303        # It should be possible to rewrite this routine in an array-friendly
[1433]304        # form like _tri_circumradius()  if we need to speed things up.
305        # Hard to do though.
[674]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]
[1420]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]])
[674]317            elen = sqrt(dx*dx+dy*dy)
318            # really only need sign - not angle value:
[1420]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]])
[674]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],\
[1420]338                       min(self.triradius[t],self.triradius[trinbr[i]]),\
339                          max(self.triradius[t],self.triradius[trinbr[i]]) ])
[674]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):
[1420]353        """
354        for each vertex find pairs
355        (min_adj_triradius, max_adj_triradius)
356        """
[674]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 ]
[684]369            try:
370                vertexinterval[i] = [min(radii), max(radii)]
371                if vertexnbrs[i][-1]==-1:
372                    vertexinterval[i][1]=INF
373            except ValueError:
[1433]374                raise AlphaError
[1420]375
[674]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        """
[691]393        Given an alpha value,
[959]394        return the indices of edges on the boundary of the alpha-shape
[674]395        """
396        def reg_edge(k):
[1420]397            return self.edgeinterval[k][1]<=alpha and \
398                   self.edgeinterval[k][2]>alpha
[674]399
400        return filter(reg_edge, range(len(self.edgeinterval)))
401
[956]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):
[1420]408            return self.vertexinterval[k][0]<=alpha and \
409                   self.vertexinterval[k][1]>alpha
[956]410
411        return filter(exp_vert, range(len(self.vertexinterval)))       
412
[959]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
[987]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)
[674]436       
[987]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):
[956]451        """
[1364]452        Given the edges in self.boundary, finds the largest components.
453        The routine does this by implementing a union-find algorithm.
[956]454        """
[970]455
[987]456        #print "running _remove_holes \n"
457
458        bdry = self.boundary
[970]459       
[956]460        def findroot(i):
461            if vptr[i] < 0:
462                return i
463            k = findroot(vptr[i])
[1435]464            vptr[i] = k    # this produces "path compression" in the
465                           # union-find tree.
[956]466            return k
467
468
469
470        # get a list of unique vertex labels:
[959]471        verts = self._vertices_from_edges(bdry)
[956]472        #print "verts ", verts, "\n"
473
[1364]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
[1227]481        EMPTY = -max(verts)-len(verts)
[956]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
[1227]500                        vptr[rvr] = vptr[rvr]-1
[956]501                else:
502                    if (vptr[vr]==EMPTY):
503                        vptr[vr] = rvl
[1227]504                        vptr[rvl] = vptr[rvl]-1
[956]505                    else:
506                        if vptr[rvl] > vptr[rvr]:
[1227]507                            vptr[rvr] = vptr[rvr] + vptr[rvl]
[956]508                            vptr[rvl] = rvr
[1227]509                            vptr[vl] = rvr
[956]510                        else:
[1227]511                            vptr[rvl] = vptr[rvl] + vptr[rvr]
[956]512                            vptr[rvr] = rvl
[1227]513                            vptr[vr] = rvl
[956]514                #print "vptr: ", vptr, "\n"
515        # end edge loop
516
[1364]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"
[1227]521        # print "largest component has: ", -min(vptr), " points. \n"
[956]522        # discard the edges in the little components
[959]523        # (i.e. those components with less than 'small' fraction of bdry points)
524        cutoff = round(small*len(verts))
[1227]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)
[1364]529
530        # littleind has root indices for small components
[1435]531        littleind = [k for k in range(len(vptr)) if \
532                     (vptr[k]<0 and vptr[k]>-cutoff)] 
[956]533        if littleind:
[1364]534            # littlecomp has all vptr indices in the small components
[1435]535            littlecomp = [k for k in range(len(vptr)) if \
536                          findroot(k) in littleind]
[1364]537            # vdiscard has the vertex indices corresponding to vptr indices 
[956]538            vdiscard = [verts[k] for k in littlecomp] 
[1435]539            newbdry = [e for e in bdry if \
540                       not((e[0] in vdiscard) and (e[1] in vdiscard))]
[1364]541
542            newverts = self._vertices_from_edges(newbdry)
543            # recompute the boundary triangles
[987]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]
[1364]549                if (v0 in newverts or v1 in newverts or v2 in newverts):
[987]550                    newbt.append(bt)
[1364]551
[1435]552            self.boundarytriangle = newbt           
[956]553        else:
554            newbdry = bdry
[987]555   
[956]556        return newbdry
557
558
[987]559    def _smooth_indents(self):
[956]560        """
[1435]561        Given edges in bdry, test for acute-angle triangular indents
562        and remove them.
[956]563        """
[970]564       
[987]565        #print "running _smooth_indents \n"
[959]566
[987]567        bdry = self.boundary
568        bdrytri = self.boundarytriangle
569       
570        # find boundary triangles that have two edges in bdry
[1364]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
[987]576        b2etri = []
577        for ind in bdrytri:
578            bect = 0
[1364]579            v2ind = [0,1,2]
[959]580            for j in [0,1,2]:
[987]581                eda = (self.deltri[ind][(j+1)%3], self.deltri[ind][(j+2)%3])
[1435]582                edb = (self.deltri[ind][(j+2)%3], self.deltri[ind][(j+1)%3])
[987]583                if eda in bdry or edb in bdry:
584                    bect +=1
[1364]585                    v2ind.remove(j)
[987]586            if bect==2:
[1364]587                b2etri.append((ind,v2ind[0]))
[959]588
589        # test the bdrytri triangles for acute angles
[987]590        acutetri = []
591        for tind in b2etri:
[1364]592            tri = self.deltri[tind[0]]
593           
[1435]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..
[1364]602            if anglesign[tind[1]] > 0:
603                acutetri.append(tind[0])
[959]604
[987]605        #print "acute boundary triangles ", acutetri
[959]606
[1435]607        # adjust the bdry edges and triangles by adding
608        #in the acutetri triangles
[987]609        for pind in acutetri:
610            bdrytri.remove(pind) 
[959]611            tri = self.deltri[pind]
612            for i in [0,1,2]:
[987]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)
[959]620       
[987]621        #print "new boundary ", newbdry
622        return newbdry
[956]623
[987]624    def _expand_pinch(self):
[956]625        """
[970]626        Given edges in bdry, test for vertices with more than 2 incident edges.
[1432]627        Expand by adding back in associated triangles.
[956]628        """
[987]629        #print "running _expand_pinch \n"
[959]630
[987]631        bdry = self.boundary
632        bdrytri = self.boundarytriangle
633       
[959]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()
[1435]638        probv = [v[k] for k in range(len(v)) \
639                 if (v[k]!=v[k-1] and v.count(v[k])>2) ]
[987]640        #print "problem vertices: ", probv 
[959]641
[987]642        # find boundary triangles that have at least one vertex in probv
[970]643        probtri = []
[987]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)
[959]650
[987]651        #print "problem boundary triangle indices ", probtri
[959]652
[970]653        # "add in" the problem triangles
[987]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]))
[959]659
[987]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
[674]670#-------------------------------------------------------------
671if __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
[691]694        #print "about to call alpha shape routine \n"
[674]695        alpha_shape_via_files(point_file, boundary_file, alpha)
696       
Note: See TracBrowser for help on using the repository browser.