source: inundation/ga/storm_surge/alpha_shape/alpha_shape.py @ 684

Last change on this file since 684 was 684, checked in by duncan, 20 years ago

update from Vanessa

File size: 11.2 KB
Line 
1"""Alpha shape   
2"""
3
4import exceptions
5from Numeric import array, Float, divide_safe, sqrt
6from load_mesh.loadASCII import load_xya_file, export_boundary_file
7
8class PointError(exceptions.Exception): pass
9
10OUTPUT_FILE_TITLE = "# The alpha shape boundary defined by point index pairs of edges"
11INF = pow(10,9)
12
13def alpha_shape_via_files(point_file, boundary_file, alpha= None):
14   
15    from load_mesh.loadASCII import load_xya_file
16   
17    point_dict = load_xya_file(point_file)
18    points = point_dict['pointlist']
19    #title_string = point_dict['title']
20   
21    alpha = Alpha_Shape(points)
22    alpha.write_boundary(boundary_file)
23
24class Alpha_Shape:
25
26    def __init__(self, points, alpha = None):
27
28        """
29        Inputs:
30       
31          points: List of coordinate pairs [[x1, y1],[x2, y2]..]
32
33          alpha: alpha shape parameter
34         
35        """
36        self._set_points(points)
37        self._alpha_shape_algorithm(alpha)
38           
39
40    def _set_points(self, points):
41        """
42        """
43        # print "setting points"
44        if len (points) <= 2:
45            raise PointError, "Too few points to find an alpha shape"
46       
47        #Convert input to Numeric arrays
48        self.points = array(points).astype(Float)
49
50   
51    def write_boundary(self,file_name):
52        """
53        Write the boundary to a file
54        """
55        #print " this info will be in the file"
56        export_boundary_file(file_name, self.get_boundary(),
57                             OUTPUT_FILE_TITLE, delimiter = ',')
58   
59    def get_boundary(self):
60        """
61        """
62        #print "getting boundary"
63        return self.boundary
64
65    def get_delaunay(self):
66        """
67        """
68        return self.deltri
69   
70           
71    def _alpha_shape_algorithm(self,alpha):
72       
73        # A python style guide suggests using _[function name] to
74        # specify internal functions
75        # - this is the first time I've used it though - DSG
76
77        #print "starting alpha shape algorithm"
78
79        # Build Delaunay triangulation
80        import triang
81        points = []
82        seglist = []
83        holelist = []
84        regionlist = []
85        pointattlist = []
86        segattlist = []
87        trilist = []
88
89        points = [(pt[0], pt[1]) for pt in self.points]
90        pointattlist = [ [] for i in range(len(points)) ] 
91        mode = "Qzcn"
92        print "computing delaunay triangulation ... \n" 
93        tridata = triang.genMesh(points,seglist,holelist,regionlist,pointattlist,segattlist,trilist,mode)
94        #print tridata
95        #print "point attribute list: ", tridata['generatedpointattributelist'],"\n"
96        #print "hull segments: ", tridata['generatedsegmentlist'], "\n"
97        self.deltri = tridata['generatedtrianglelist']
98        self.deltrinbr = tridata['generatedtriangleneighborlist']
99        self.hulledges = tridata['generatedsegmentlist']                             
100
101        print "Number of delaunay triangles: ", len(self.deltri), "\n"
102        #print "deltrinbrs: ", self.deltrinbr, "\n"
103
104        # Build Alpha table
105        print "Building alpha table ... \n" 
106        self._tri_circumradius()
107        print "Largest circumradius ", max(self.triradius)
108        self._edge_intervals()
109        self._vertex_intervals()
110
111        if alpha==None:
112            # Find optimum alpha
113            # Ken Clarkson's hull program uses smallest alpha so that
114            # every vertex is non-singular so...
115            self.optimum_alpha = max([ interval[0] for interval in self.vertexinterval])
116            print "optimum alpha ", self.optimum_alpha
117
118            #alpha_tri = self.get_alpha_triangles(self.optimum_alpha)
119            #print "alpha shape triangles ", alpha_tri
120
121            reg_edge = self.get_regular_edges(self.optimum_alpha)
122            #print "alpha boundary edges", reg_edge
123        else:
124            reg_edge = self.get_regular_edges(alpha)
125   
126        self.boundary = [self.edge[k] for k in reg_edge]
127        #print "alpha boundary edges ", self.boundary
128           
129        # this produces a baaad boundary
130        #boundary = []
131        #for point_index in range(len(self.points)-1):
132        #    boundary.append([point_index,point_index +1])
133        #boundary.append([len(self.points)-1,0])
134        #self.boundary = boundary
135        return
136
137    def _tri_circumradius(self):
138
139        # compute circumradii of the delaunay triangles
140        x = self.points[:,0]
141        y = self.points[:,1]
142        ind1 = [self.deltri[j][0] for j in range(len(self.deltri))]
143        ind2 = [self.deltri[j][1] for j in range(len(self.deltri))]
144        ind3 = [self.deltri[j][2] for j in range(len(self.deltri))]
145
146        x1 = array([ x[j] for j in ind1 ])
147        y1 = array([ y[j] for j in ind1 ])
148        x2 = array([ x[j] for j in ind2 ])
149        y2 = array([ y[j] for j in ind2 ])
150        x3 = array([ x[j] for j in ind3 ])
151        y3 = array([ y[j] for j in ind3 ])
152
153        x21 = x2-x1
154        x31 = x3-x1
155        y21 = y2-y1
156        y31 = y3-y1
157
158        dist21 = x21*x21 + y21*y21
159        dist31 = x31*x31 + y31*y31
160
161        denom = x21*y31 - x31*y21
162        #print "denom = ", denom
163               
164        # dx/2, dy/2 give circumcenter relative to x1,y1.
165        #dx = (y31*dist21 - y21*dist31)/denom
166        #dy = (x21*dist31 - x31*dist21)/denom
167        # from Numeric import divide_safe
168        try:
169            dx = divide_safe(y31*dist21 - y21*dist31,denom)
170            dy = divide_safe(x21*dist31 - x31*dist21,denom)
171        except ZeroDivisionError:
172            print "Found some degenerate triangles."
173            raise
174            # really need to take care of cases when denom = 0.
175            # eg: something like:
176            # print "Handling degenerate triangles."
177            # ... add some random displacments to trouble points?
178           
179        self.triradius = 0.5*sqrt(dx*dx + dy*dy)
180        #print "triangle radii", self.triradius
181        return
182
183    def _edge_intervals(self):
184        # for each edge, find triples
185        # (length/2, min_adj_triradius, max_adj_triradius) if unattached
186        # (min_adj_triradius, min_adj_triradius, max_adj_triradius) if attached.
187
188        # It should be possible to rewrite this routine in an array-friendly form
189        # like _tri_circumradius()  if we need to speed things up.
190
191        edges = []
192        edgenbrs = []
193        edgeinterval = []
194        for t in range(len(self.deltri)):
195            tri = self.deltri[t]
196            trinbr = self.deltrinbr[t]
197            dx = array([self.points[tri[(i+1)%3],0] - self.points[tri[(i+2)%3],0] for i in [0,1,2]])
198            dy = array([self.points[tri[(i+1)%3],1] - self.points[tri[(i+2)%3],1] for i in [0,1,2]])
199            elen = sqrt(dx*dx+dy*dy)
200            #angle = array([(-dx[(i+1)%3]*dx[(i+2)%3]-dy[(i+1)%3]*dy[(i+2)%3])/(elen[(i+1)%3]*elen[(i+2)%3]) for i in [0,1,2]])
201            #angle = arccos(angle)
202            # really only need sign - not angle value:
203            anglesign = array([(-dx[(i+1)%3]*dx[(i+2)%3]-dy[(i+1)%3]*dy[(i+2)%3]) for i in [0,1,2]])
204           
205            #print "dx ",dx,"\n"
206            #print "dy ",dy,"\n"
207            #print "edge lengths of triangle ",t,"\t",elen,"\n"
208            #print "angles ",angle,"\n"
209                         
210            for i in [0,1,2]:
211                j = (i+1)%3
212                k = (i+2)%3
213                if trinbr[i]==-1:
214                    edges.append((tri[j], tri[k]))
215                    edgenbrs.append((t, -1))
216                    edgeinterval.append([0.5*elen[i], self.triradius[t], INF])
217                elif (tri[j]<tri[k]):
218                    edges.append((tri[j], tri[k]))
219                    edgenbrs.append((t, trinbr[i]))
220                    edgeinterval.append([0.5*elen[i],\
221                                         min(self.triradius[t],self.triradius[trinbr[i]]),\
222                                             max(self.triradius[t],self.triradius[trinbr[i]]) ])
223                else:
224                    continue
225                #if angle[i]>pi/2:
226                if anglesign[i] < 0: 
227                    edgeinterval[-1][0] = edgeinterval[-1][1]
228                   
229        self.edge = edges
230        self.edgenbr = edgenbrs
231        self.edgeinterval = edgeinterval
232        #print "edges: ",edges, "\n"
233        #print "edge nbrs:", edgenbrs ,"\n"
234        #print "edge intervals: ",edgeinterval , "\n"
235
236    def _vertex_intervals(self):
237        # for each vertex find pairs
238        # (min_adj_triradius, max_adj_triradius)
239        nv = len(self.points)
240        vertexnbrs = [ [] for i in range(nv)]
241        vertexinterval = [ [] for i in range(nv)] 
242        for t in range(len(self.deltri)):
243            for j in self.deltri[t]:
244                vertexnbrs[j].append(t)
245        for h in range(len(self.hulledges)):
246            for j in self.hulledges[h]:
247                vertexnbrs[j].append(-1)
248       
249        for i in range(nv):
250            radii = [ self.triradius[t] for t in vertexnbrs[i] if t>=0 ]
251            try:
252                vertexinterval[i] = [min(radii), max(radii)]
253                if vertexnbrs[i][-1]==-1:
254                    vertexinterval[i][1]=INF
255            except ValueError:
256                print "Vertex %i is isolated?"%i
257                print "coords: ",self.points[i]
258                print "Vertex nbrs: ", vertexnbrs[i]
259                print "nbr radii: ",radii
260                vertexinterval[i] = [0,INF]
261                pass
262
263        self.vertexnbr = vertexnbrs
264        self.vertexinterval = vertexinterval
265        #print "vertex nbrs ", vertexnbrs, "\n"
266        #print "vertex intervals ",vertexinterval, "\n"
267       
268    def get_alpha_triangles(self,alpha):
269        """
270        Given an alpha value,
271        return indices of triangles in the alpha-shape
272        """
273        def tri_rad_lta(k):
274            return self.triradius[k]<=alpha
275
276        return filter(tri_rad_lta, range(len(self.triradius)))
277
278    def get_regular_edges(self,alpha):
279        """
280        Given and alpha value,
281        return the edges on the boundary of the alpha-shape
282        """
283        def reg_edge(k):
284            return self.edgeinterval[k][1]<=alpha and self.edgeinterval[k][2]>alpha
285
286        return filter(reg_edge, range(len(self.edgeinterval)))
287
288       
289   
290#-------------------------------------------------------------
291if __name__ == "__main__":
292    """
293    Load in a data point file.
294    Determine the alpha shape boundary
295    Save the boundary to a file.
296
297    usage: alpha_shape.py point_file.xya boundary_file.bnd [alpha]
298   
299    The alpha value is optional.
300    """
301   
302    import os, sys
303    usage = "usage: %s point_file.xya boundary_file.bnd [alpha]"%os.path.basename(sys.argv[0])
304    # I made up the .bnd affix. Other ideas welcome. -DSG
305    if len(sys.argv) < 3:
306        print usage
307    else:
308        point_file = sys.argv[1]
309        boundary_file = sys.argv[2]
310        if len(sys.argv) > 4:
311            alpha = sys.argv[3]
312        else:
313            alpha = None
314
315        print "about to call alpha shape routine \n"
316        alpha_shape_via_files(point_file, boundary_file, alpha)
317       
Note: See TracBrowser for help on using the repository browser.