source: anuga_core/source/anuga/utilities/polygon.py @ 7516

Last change on this file since 7516 was 7516, checked in by rwilson, 15 years ago

Added polygon fill - UNTESTED!

File size: 46.5 KB
Line 
1#!/usr/bin/env python
2
3"""Polygon manipulations"""
4
5import numpy as num
6
7from math import sqrt
8from anuga.utilities.numerical_tools import ensure_numeric
9from anuga.geospatial_data.geospatial_data import ensure_absolute, Geospatial_data
10from anuga.config import netcdf_float
11import anuga.utilities.log as log
12
13
14##
15# @brief Determine whether a point is on a line segment.
16# @param point (x, y) of point in question (tuple, list or array).
17# @param line ((x1,y1), (x2,y2)) for line (tuple, list or array).
18# @param rtol Relative error for 'close'.
19# @param atol Absolute error for 'close'.
20# @return True or False.
21def point_on_line(point, line, rtol=1.0e-5, atol=1.0e-8):
22    """Determine whether a point is on a line segment
23
24    Input:
25        point is given by [x, y]
26        line is given by [x0, y0], [x1, y1]] or
27        the equivalent 2x2 numeric array with each row corresponding to a point.
28
29    Output:
30
31    Note: Line can be degenerate and function still works to discern coinciding
32          points from non-coinciding.
33    """
34
35    point = ensure_numeric(point)
36    line = ensure_numeric(line)
37
38    res = _point_on_line(point[0], point[1],
39                         line[0,0], line[0,1],
40                         line[1,0], line[1,1],
41                         rtol, atol)
42
43    return bool(res)
44
45
46######
47# Result functions used in intersection() below for collinear lines.
48# (p0,p1) defines line 0, (p2,p3) defines line 1.
49######
50
51# result functions for possible states
52def lines_dont_coincide(p0,p1,p2,p3):               return (3, None)
53def lines_0_fully_included_in_1(p0,p1,p2,p3):       return (2,
54                                                            num.array([p0,p1]))
55def lines_1_fully_included_in_0(p0,p1,p2,p3):       return (2,
56                                                            num.array([p2,p3]))
57def lines_overlap_same_direction(p0,p1,p2,p3):      return (2,
58                                                            num.array([p0,p3]))
59def lines_overlap_same_direction2(p0,p1,p2,p3):     return (2,
60                                                            num.array([p2,p1]))
61def lines_overlap_opposite_direction(p0,p1,p2,p3):  return (2,
62                                                            num.array([p0,p2]))
63def lines_overlap_opposite_direction2(p0,p1,p2,p3): return (2,
64                                                            num.array([p3,p1]))
65
66# this function called when an impossible state is found
67def lines_error(p1, p2, p3, p4):
68    raise RuntimeError, ('INTERNAL ERROR: p1=%s, p2=%s, p3=%s, p4=%s'
69                         % (str(p1), str(p2), str(p3), str(p4)))
70
71#                     0s1    0e1    1s0    1e0   # line 0 starts on 1, 0 ends 1, 1 starts 0, 1 ends 0
72collinear_result = { (False, False, False, False): lines_dont_coincide,
73                     (False, False, False, True ): lines_error,
74                     (False, False, True,  False): lines_error,
75                     (False, False, True,  True ): lines_1_fully_included_in_0,
76                     (False, True,  False, False): lines_error,
77                     (False, True,  False, True ): lines_overlap_opposite_direction2,
78                     (False, True,  True,  False): lines_overlap_same_direction2,
79                     (False, True,  True,  True ): lines_1_fully_included_in_0,
80                     (True,  False, False, False): lines_error,
81                     (True,  False, False, True ): lines_overlap_same_direction,
82                     (True,  False, True,  False): lines_overlap_opposite_direction,
83                     (True,  False, True,  True ): lines_1_fully_included_in_0,
84                     (True,  True,  False, False): lines_0_fully_included_in_1,
85                     (True,  True,  False, True ): lines_0_fully_included_in_1,
86                     (True,  True,  True,  False): lines_0_fully_included_in_1,
87                     (True,  True,  True,  True ): lines_0_fully_included_in_1
88                   }
89
90##
91# @brief Finds intersection point of two line segments.
92# @param line0 First line ((x1,y1), (x2,y2)).
93# @param line1 Second line ((x1,y1), (x2,y2)).
94# @param rtol Relative error for 'close'.
95# @param atol Absolute error for 'close'.
96# @return (status, value) where:
97#             status = 0 - no intersection, value set to None
98#                      1 - intersection found, value=(x,y)
99#                      2 - lines collienar, overlap, value=overlap segment
100#                      3 - lines collinear, no overlap, value is None
101#                      4 - lines parallel, value is None
102def intersection(line0, line1, rtol=1.0e-5, atol=1.0e-8):
103    """Returns intersecting point between two line segments.
104
105    However, if parallel lines coincide partly (i.e. share a common segment),
106    the line segment where lines coincide is returned
107
108    Inputs:
109        line0, line1: Each defined by two end points as in: [[x0, y0], [x1, y1]]
110                      A line can also be a 2x2 numpy array with each row
111                      corresponding to a point.
112
113    Output:
114        status, value - where status and value is interpreted as follows:
115        status == 0: no intersection, value set to None.
116        status == 1: intersection point found and returned in value as [x,y].
117        status == 2: Collinear overlapping lines found.
118                     Value takes the form [[x0,y0], [x1,y1]].
119        status == 3: Collinear non-overlapping lines. Value set to None.
120        status == 4: Lines are parallel. Value set to None.
121    """
122
123    # FIXME (Ole): Write this in C
124
125    line0 = ensure_numeric(line0, num.float)
126    line1 = ensure_numeric(line1, num.float)
127
128    x0 = line0[0,0]; y0 = line0[0,1]
129    x1 = line0[1,0]; y1 = line0[1,1]
130
131    x2 = line1[0,0]; y2 = line1[0,1]
132    x3 = line1[1,0]; y3 = line1[1,1]
133
134    denom = (y3-y2)*(x1-x0) - (x3-x2)*(y1-y0)
135    u0 = (x3-x2)*(y0-y2) - (y3-y2)*(x0-x2)
136    u1 = (x2-x0)*(y1-y0) - (y2-y0)*(x1-x0)
137
138    if num.allclose(denom, 0.0, rtol=rtol, atol=atol):
139        # Lines are parallel - check if they are collinear
140        if num.allclose([u0, u1], 0.0, rtol=rtol, atol=atol):
141            # We now know that the lines are collinear
142            state_tuple = (point_on_line([x0, y0], line1, rtol=rtol, atol=atol),
143                           point_on_line([x1, y1], line1, rtol=rtol, atol=atol),
144                           point_on_line([x2, y2], line0, rtol=rtol, atol=atol),
145                           point_on_line([x3, y3], line0, rtol=rtol, atol=atol))
146
147            return collinear_result[state_tuple]([x0,y0], [x1,y1],
148                                                 [x2,y2], [x3,y3])
149        else:
150            # Lines are parallel but aren't collinear
151            return 4, None #FIXME (Ole): Add distance here instead of None
152    else:
153        # Lines are not parallel, check if they intersect
154        u0 = u0/denom
155        u1 = u1/denom
156
157        x = x0 + u0*(x1-x0)
158        y = y0 + u0*(y1-y0)
159
160        # Sanity check - can be removed to speed up if needed
161        assert num.allclose(x, x2 + u1*(x3-x2), rtol=rtol, atol=atol)
162        assert num.allclose(y, y2 + u1*(y3-y2), rtol=rtol, atol=atol)
163
164        # Check if point found lies within given line segments
165        if 0.0 <= u0 <= 1.0 and 0.0 <= u1 <= 1.0:
166            # We have intersection
167            return 1, num.array([x, y])
168        else:
169            # No intersection
170            return 0, None
171
172##
173# @brief Finds intersection point of two line segments.
174# @param line0 First line ((x1,y1), (x2,y2)).
175# @param line1 Second line ((x1,y1), (x2,y2)).
176# @return (status, value) where:
177#             status = 0 - no intersection, value set to None
178#                      1 - intersection found, value=(x,y)
179#                      2 - lines collienar, overlap, value=overlap segment
180#                      3 - lines collinear, no overlap, value is None
181#                      4 - lines parallel, value is None
182# @note Wrapper for C function.
183def NEW_C_intersection(line0, line1):
184    """Returns intersecting point between two line segments.
185
186    However, if parallel lines coincide partly (i.e. share a common segment),
187    the line segment where lines coincide is returned
188
189    Inputs:
190        line0, line1: Each defined by two end points as in: [[x0, y0], [x1, y1]]
191                      A line can also be a 2x2 numpy array with each row
192                      corresponding to a point.
193
194    Output:
195        status, value - where status and value is interpreted as follows:
196        status == 0: no intersection, value set to None.
197        status == 1: intersection point found and returned in value as [x,y].
198        status == 2: Collinear overlapping lines found.
199                     Value takes the form [[x0,y0], [x1,y1]].
200        status == 3: Collinear non-overlapping lines. Value set to None.
201        status == 4: Lines are parallel. Value set to None.
202    """
203
204    line0 = ensure_numeric(line0, num.float)
205    line1 = ensure_numeric(line1, num.float)
206
207    status, value = _intersection(line0[0,0], line0[0,1],
208                                  line0[1,0], line0[1,1],
209                                  line1[0,0], line1[0,1],
210                                  line1[1,0], line1[1,1])
211
212    return status, value
213
214def is_inside_triangle(point, triangle, 
215                       closed=True, 
216                       rtol=1.0e-12,
217                       atol=1.0e-12,                     
218                       check_inputs=True, 
219                       verbose=False):
220    """Determine if one point is inside a triangle
221   
222    This uses the barycentric method:
223   
224    Triangle is A, B, C
225    Point P can then be written as
226   
227    P = A + alpha * (C-A) + beta * (B-A)
228    or if we let
229    v=P-A, v0=C-A, v1=B-A   
230   
231    v = alpha*v0 + beta*v1
232
233    Dot this equation by v0 and v1 to get two:
234   
235    dot(v0, v) = alpha*dot(v0, v0) + beta*dot(v0, v1)
236    dot(v1, v) = alpha*dot(v1, v0) + beta*dot(v1, v1)   
237   
238    or if a_ij = dot(v_i, v_j) and b_i = dot(v_i, v)
239    the matrix equation:
240   
241    a_00 a_01   alpha     b_0
242                       =
243    a_10 a_11   beta      b_1
244   
245    Solving for alpha and beta yields:
246   
247    alpha = (b_0*a_11 - b_1*a_01)/denom
248    beta =  (b_1*a_00 - b_0*a_10)/denom
249   
250    with denom = a_11*a_00 - a_10*a_01
251   
252    The point is in the triangle whenever
253    alpha and beta and their sums are in the unit interval.
254   
255    rtol and atol will determine how close the point has to be to the edge
256    before it is deemed to be on the edge.
257   
258    """
259
260    triangle = ensure_numeric(triangle)       
261    point = ensure_numeric(point, num.float)   
262   
263    if check_inputs is True:
264        msg = 'is_inside_triangle must be invoked with one point only'
265        assert num.allclose(point.shape, [2]), msg
266   
267   
268    # Use C-implementation
269    return bool(_is_inside_triangle(point, triangle, int(closed), rtol, atol))
270   
271   
272
273    # FIXME (Ole): The rest of this function has been made
274    # obsolete by the C extension.
275   
276    # Quickly reject points that are clearly outside
277    if point[0] < min(triangle[:,0]): return False 
278    if point[0] > max(triangle[:,0]): return False   
279    if point[1] < min(triangle[:,1]): return False
280    if point[1] > max(triangle[:,1]): return False       
281
282
283    # Start search   
284    A = triangle[0, :]
285    B = triangle[1, :]
286    C = triangle[2, :]
287   
288    # Now check if point lies wholly inside triangle
289    v0 = C-A
290    v1 = B-A       
291       
292    a00 = num.inner(v0, v0)
293    a10 = a01 = num.inner(v0, v1)
294    a11 = num.inner(v1, v1)
295   
296    denom = a11*a00 - a01*a10
297   
298    if abs(denom) > 0.0:
299        v = point-A
300        b0 = num.inner(v0, v)       
301        b1 = num.inner(v1, v)           
302   
303        alpha = (b0*a11 - b1*a01)/denom
304        beta = (b1*a00 - b0*a10)/denom       
305
306        if (alpha > 0.0) and (beta > 0.0) and (alpha+beta < 1.0):
307            return True 
308
309
310    if closed is True:
311        # Check if point lies on one of the edges
312       
313        for X, Y in [[A,B], [B,C], [C,A]]:
314            res = _point_on_line(point[0], point[1],
315                                 X[0], X[1],
316                                 Y[0], Y[1],
317                                 rtol, atol)
318           
319            if res:
320                return True
321               
322    return False
323
324def is_inside_polygon_quick(point, polygon, closed=True, verbose=False):
325    """Determine if one point is inside a polygon
326    Both point and polygon are assumed to be numeric arrays or lists and
327    no georeferencing etc or other checks will take place.
328   
329    As such it is faster than is_inside_polygon
330    """
331
332    # FIXME(Ole): This function isn't being used
333    polygon = ensure_numeric(polygon, num.float)
334    points = ensure_numeric(point, num.float) # Convert point to array of points
335    points = num.ascontiguousarray(points[num.newaxis, :])
336    msg = ('is_inside_polygon() must be invoked with one point only.\n'
337           'I got %s and converted to %s' % (str(point), str(points.shape)))
338    assert points.shape[0] == 1 and points.shape[1] == 2, msg
339   
340    indices = num.zeros(1, num.int)
341
342    count = _separate_points_by_polygon(points, polygon, indices,
343                                        int(closed), int(verbose))
344
345    return count > 0
346
347
348def is_inside_polygon(point, polygon, closed=True, verbose=False):
349    """Determine if one point is inside a polygon
350
351    See inside_polygon for more details
352    """
353
354    indices = inside_polygon(point, polygon, closed, verbose)
355
356    if indices.shape[0] == 1:
357        return True
358    elif indices.shape[0] == 0:
359        return False
360    else:
361        msg = 'is_inside_polygon must be invoked with one point only'
362        raise msg
363
364##
365# @brief Determine which of a set of points are inside a polygon.
366# @param points A set of points (tuple, list or array).
367# @param polygon A set of points defining a polygon (tuple, list or array).
368# @param closed True if points on boundary are considered 'inside' polygon.
369# @param verbose True if this function is to be verbose.
370# @return A list of indices of points inside the polygon.
371def inside_polygon(points, polygon, closed=True, verbose=False):
372    """Determine points inside a polygon
373
374       Functions inside_polygon and outside_polygon have been defined in
375       terms of separate_by_polygon which will put all inside indices in
376       the first part of the indices array and outside indices in the last
377
378       See separate_points_by_polygon for documentation
379
380       points and polygon can be a geospatial instance,
381       a list or a numeric array
382    """
383
384    try:
385        points = ensure_absolute(points)
386    except NameError, e:
387        raise NameError, e
388    except:
389        # If this fails it is going to be because the points can't be
390        # converted to a numeric array.
391        msg = 'Points could not be converted to numeric array' 
392        raise Exception, msg
393
394    polygon = ensure_absolute(polygon)       
395    try:
396        polygon = ensure_absolute(polygon)
397    except NameError, e:
398        raise NameError, e
399    except:
400        # If this fails it is going to be because the points can't be
401        # converted to a numeric array.
402        msg = ('Polygon %s could not be converted to numeric array'
403               % (str(polygon)))
404        raise Exception, msg
405
406    if len(points.shape) == 1:
407        # Only one point was passed in. Convert to array of points
408        points = num.reshape(points, (1,2))
409
410    indices, count = separate_points_by_polygon(points, polygon,
411                                                closed=closed,
412                                                verbose=verbose)
413
414    # Return indices of points inside polygon
415    return indices[:count]
416
417##
418# @brief Determine if one point is outside a polygon.
419# @param point The point of interest.
420# @param polygon The polygon to test inclusion in.
421# @param closed True if points on boundary are considered 'inside' polygon.
422# @param verbose True if this function is to be verbose.
423# @return True if point is outside the polygon.
424# @note Uses inside_polygon() to do the work.
425def is_outside_polygon(point, polygon, closed=True, verbose=False,
426                       points_geo_ref=None, polygon_geo_ref=None):
427    """Determine if one point is outside a polygon
428
429    See outside_polygon for more details
430    """
431
432    indices = outside_polygon(point, polygon, closed, verbose)
433
434    if indices.shape[0] == 1:
435        return True
436    elif indices.shape[0] == 0:
437        return False
438    else:
439        msg = 'is_outside_polygon must be invoked with one point only'
440        raise Exception, msg
441
442##
443# @brief Determine which of a set of points are outside a polygon.
444# @param points A set of points (tuple, list or array).
445# @param polygon A set of points defining a polygon (tuple, list or array).
446# @param closed True if points on boundary are considered 'inside' polygon.
447# @param verbose True if this function is to be verbose.
448# @return A list of indices of points outside the polygon.
449def outside_polygon(points, polygon, closed = True, verbose = False):
450    """Determine points outside a polygon
451
452       Functions inside_polygon and outside_polygon have been defined in
453       terms of separate_by_polygon which will put all inside indices in
454       the first part of the indices array and outside indices in the last
455
456       See separate_points_by_polygon for documentation
457    """
458
459    try:
460        points = ensure_numeric(points, num.float)
461    except NameError, e:
462        raise NameError, e
463    except:
464        msg = 'Points could not be converted to numeric array'
465        raise Exception, msg
466
467    try:
468        polygon = ensure_numeric(polygon, num.float)
469    except NameError, e:
470        raise NameError, e
471    except:
472        msg = 'Polygon could not be converted to numeric array'
473        raise Exception, msg
474
475    if len(points.shape) == 1:
476        # Only one point was passed in. Convert to array of points
477        points = num.reshape(points, (1,2))
478
479    indices, count = separate_points_by_polygon(points, polygon,
480                                                closed=closed,
481                                                verbose=verbose)
482
483    # Return indices of points outside polygon
484    if count == len(indices):
485        # No points are outside
486        return num.array([])
487    else:
488        return indices[count:][::-1]  #return reversed
489
490##
491# @brief Separate a list of points into two sets inside+outside a polygon.
492# @param points A set of points (tuple, list or array).
493# @param polygon A set of points defining a polygon (tuple, list or array).
494# @param closed True if points on boundary are considered 'inside' polygon.
495# @param verbose True if this function is to be verbose.
496# @return A tuple (in, out) of point indices for poinst inside amd outside.
497def in_and_outside_polygon(points, polygon, closed=True, verbose=False):
498    """Determine points inside and outside a polygon
499
500       See separate_points_by_polygon for documentation
501
502       Returns an array of points inside and array of points outside the polygon
503    """
504
505    try:
506        points = ensure_numeric(points, num.float)
507    except NameError, e:
508        raise NameError, e
509    except:
510        msg = 'Points could not be converted to numeric array'
511        raise Exception, msg
512
513    try:
514        polygon = ensure_numeric(polygon, num.float)
515    except NameError, e:
516        raise NameError, e
517    except:
518        msg = 'Polygon could not be converted to numeric array'
519        raise Exception, msg
520
521    if len(points.shape) == 1:
522        # Only one point was passed in. Convert to array of points
523        points = num.reshape(points, (1,2))
524
525    indices, count = separate_points_by_polygon(points, polygon,
526                                                closed=closed,
527                                                verbose=verbose)
528
529    # Returns indices of points inside and indices of points outside
530    # the polygon
531    if count == len(indices):
532        # No points are outside
533        return indices[:count],[]
534    else:
535        return  indices[:count], indices[count:][::-1]  #return reversed
536
537##
538# @brief Sort a list of points into contiguous points inside+outside a polygon.
539# @param points A set of points (tuple, list or array).
540# @param polygon A set of points defining a polygon (tuple, list or array).
541# @param closed True if points on boundary are considered 'inside' polygon.
542# @param verbose True if this function is to be verbose.
543# @return (indices, count) where indices are point indices and count is the
544#          delimiter index between point inside (on left) and others.
545def separate_points_by_polygon(points, polygon,
546                               closed=True, 
547                               check_input=True,
548                               verbose=False):
549    """Determine whether points are inside or outside a polygon
550
551    Input:
552       points - Tuple of (x, y) coordinates, or list of tuples
553       polygon - list of vertices of polygon
554       closed - (optional) determine whether points on boundary should be
555       regarded as belonging to the polygon (closed = True)
556       or not (closed = False)
557       check_input: Allows faster execution if set to False
558
559    Outputs:
560       indices: array of same length as points with indices of points falling
561       inside the polygon listed from the beginning and indices of points
562       falling outside listed from the end.
563
564       count: count of points falling inside the polygon
565
566       The indices of points inside are obtained as indices[:count]
567       The indices of points outside are obtained as indices[count:]
568
569    Examples:
570       U = [[0,0], [1,0], [1,1], [0,1]] #Unit square
571
572       separate_points_by_polygon( [[0.5, 0.5], [1, -0.5], [0.3, 0.2]], U)
573       will return the indices [0, 2, 1] and count == 2 as only the first
574       and the last point are inside the unit square
575
576    Remarks:
577       The vertices may be listed clockwise or counterclockwise and
578       the first point may optionally be repeated.
579       Polygons do not need to be convex.
580       Polygons can have holes in them and points inside a hole is
581       regarded as being outside the polygon.
582
583    Algorithm is based on work by Darel Finley,
584    http://www.alienryderflex.com/polygon/
585
586    Uses underlying C-implementation in polygon_ext.c
587    """
588
589    if check_input:
590        #Input checks
591        assert isinstance(closed, bool), 'Keyword argument "closed" must be boolean'
592        assert isinstance(verbose, bool), 'Keyword argument "verbose" must be boolean'
593
594        try:
595            points = ensure_numeric(points, num.float)
596        except NameError, e:
597            raise NameError, e
598        except:
599            msg = 'Points could not be converted to numeric array'
600            raise msg
601
602        try:
603            polygon = ensure_numeric(polygon, num.float)
604        except NameError, e:
605            raise NameError, e
606        except:
607            msg = 'Polygon could not be converted to numeric array'
608            raise msg
609
610        msg = 'Polygon array must be a 2d array of vertices'
611        assert len(polygon.shape) == 2, msg
612
613        msg = 'Polygon array must have two columns' 
614        assert polygon.shape[1]==2, msg
615
616        msg = ('Points array must be 1 or 2 dimensional. '
617               'I got %d dimensions' % len(points.shape))
618        assert 0 < len(points.shape) < 3, msg
619
620        if len(points.shape) == 1:
621            # Only one point was passed in.  Convert to array of points.
622            points = num.reshape(points, (1,2))
623   
624            msg = ('Point array must have two columns (x,y), '
625                   'I got points.shape[1]=%d' % points.shape[0])
626            assert points.shape[1]==2, msg
627
628       
629            msg = ('Points array must be a 2d array. I got %s.'
630                   % str(points[:30]))
631            assert len(points.shape)==2, msg
632
633            msg = 'Points array must have two columns'
634            assert points.shape[1]==2, msg
635
636    N = polygon.shape[0] # Number of vertices in polygon
637    M = points.shape[0]  # Number of points
638
639    indices = num.zeros(M, num.int)
640
641    count = _separate_points_by_polygon(points, polygon, indices,
642                                        int(closed), int(verbose))
643
644    if verbose:
645        log.critical('Found %d points (out of %d) inside polygon' % (count, M))
646
647    return indices, count
648
649##
650# @brief Determine area of a polygon.
651# @param input_polygon The polygon to get area of.
652# @return A scalar value for the polygon area.
653def polygon_area(input_polygon):
654    """ Determine area of arbitrary polygon.
655
656    Reference:     http://mathworld.wolfram.com/PolygonArea.html
657    """
658
659    # Move polygon to origin (0,0) to avoid rounding errors
660    # This makes a copy of the polygon to avoid destroying it
661    input_polygon = ensure_numeric(input_polygon)
662    min_x = min(input_polygon[:,0])
663    min_y = min(input_polygon[:,1])
664    polygon = input_polygon - [min_x, min_y]
665
666    # Compute area
667    n = len(polygon)
668    poly_area = 0.0
669
670    for i in range(n):
671        pti = polygon[i]
672        if i == n-1:
673            pt1 = polygon[0]
674        else:
675            pt1 = polygon[i+1]
676        xi = pti[0]
677        yi1 = pt1[1]
678        xi1 = pt1[0]
679        yi = pti[1]
680        poly_area += xi*yi1 - xi1*yi
681
682    return abs(poly_area/2)
683
684##
685# @brief Plot a set of polygons.
686# @param polygons_points List of polygons to plot.
687# @param style List of styles for each polygon.
688# @param figname Name to save figure to.
689# @param label Title for the plot.
690# @param verbose True if this function is to be verbose.
691# @return A list of min/max x and y values [minx, maxx, miny, maxy].
692# @note A style value is 'line' for polygons, 'outside' for points outside.
693def plot_polygons(polygons_points,
694                  style=None,
695                  figname=None,
696                  label=None,
697                  alpha=None,
698                  verbose=False):
699    """ Take list of polygons and plot.
700
701    Inputs:
702
703    polygons         - list of polygons
704
705    style            - style list corresponding to each polygon
706                     - for a polygon, use 'line'
707                     - for points falling outside a polygon, use 'outside'
708                     - style can also be user defined as in normal pylab plot.
709
710    figname          - name to save figure to
711
712    label            - title for plotA
713
714    alpha            - transparency of polygon fill, 0.0=none, 1.0=solid
715                       if not supplied, no fill.
716
717    Outputs:
718
719    - list of min and max of x and y coordinates
720    - plot of polygons
721    """
722
723    from pylab import ion, hold, plot, axis, figure, legend, savefig, xlabel, \
724                      ylabel, title, close, title, fill
725
726    assert type(polygons_points) == list, \
727                'input must be a list of polygons and/or points'
728
729    ion()
730    hold(True)
731
732    minx = 1e10
733    maxx = 0.0
734    miny = 1e10
735    maxy = 0.0
736
737    if label is None:
738        label = ''
739
740    # clamp alpha to sensible range
741    if alpha:
742        try:
743            alpha = float(alpha)
744        except ValueError:
745            alpha = None
746        else:
747            if alpha < 0.0:
748                alpha = 0.0
749            if alpha > 1.0:
750                alpha = 1.0
751
752    n = len(polygons_points)
753    colour = []
754    if style is None:
755        style_type = 'line'
756        style = []
757        for i in range(n):
758            style.append(style_type)
759            colour.append('b-')
760    else:
761        for s in style:
762            if s == 'line': colour.append('b-')
763            if s == 'outside': colour.append('r.')
764            if s <> 'line':
765                if s <> 'outside':
766                    colour.append(s)
767
768    for i, item in enumerate(polygons_points):
769        x, y = poly_xy(item)
770        if min(x) < minx: minx = min(x)
771        if max(x) > maxx: maxx = max(x)
772        if min(y) < miny: miny = min(y)
773        if max(y) > maxy: maxy = max(y)
774        plot(x,y,colour[i])
775        if alpha:
776            fill(x, y, colour[i], alpha=alpha)
777        xlabel('x')
778        ylabel('y')
779        title(label)
780
781    #raw_input('wait 1')
782    #FIXME(Ole): This makes for some strange scalings sometimes.
783    #if minx <> 0:
784    #    axis([minx*0.9,maxx*1.1,miny*0.9,maxy*1.1])
785    #else:
786    #    if miny == 0:
787    #        axis([-maxx*.01,maxx*1.1,-maxy*0.01,maxy*1.1])
788    #    else:
789    #        axis([-maxx*.01,maxx*1.1,miny*0.9,maxy*1.1])
790
791    if figname is not None:
792        savefig(figname)
793    else:
794        savefig('test_image')
795
796    close('all')
797
798    vec = [minx, maxx, miny, maxy]
799    return vec
800
801##
802# @brief
803# @param polygon A set of points defining a polygon.
804# @param verbose True if this function is to be verbose.
805# @return A tuple (x, y) of X and Y coordinates of the polygon.
806# @note We duplicate the first point so can have closed polygon in plot.
807def poly_xy(polygon, verbose=False):
808    """ this is used within plot_polygons so need to duplicate
809        the first point so can have closed polygon in plot
810    """
811
812    try:
813        polygon = ensure_numeric(polygon, num.float)
814    except NameError, e:
815        raise NameError, e
816    except:
817        msg = ('Polygon %s could not be converted to numeric array'
818               % (str(polygon)))
819        raise Exception, msg
820
821    x = polygon[:,0]
822    y = polygon[:,1]
823    x = num.concatenate((x, [polygon[0,0]]), axis = 0)
824    y = num.concatenate((y, [polygon[0,1]]), axis = 0)
825
826    return x, y
827
828
829##
830# @brief Define a class that defines a callable object for a polygon.
831# @note Object created is function: f: x,y -> z
832#       where x, y and z are vectors and z depends on whether x,y belongs
833#       to specified polygons.
834class Polygon_function:
835    """Create callable object f: x,y -> z, where a,y,z are vectors and
836    where f will return different values depending on whether x,y belongs
837    to specified polygons.
838
839    To instantiate:
840
841       Polygon_function(polygons)
842
843    where polygons is a list of tuples of the form
844
845      [ (P0, v0), (P1, v1), ...]
846
847      with Pi being lists of vertices defining polygons and vi either
848      constants or functions of x,y to be applied to points with the polygon.
849
850    The function takes an optional argument, default which is the value
851    (or function) to used for points not belonging to any polygon.
852    For example:
853
854       Polygon_function(polygons, default = 0.03)
855
856    If omitted the default value will be 0.0
857
858    Note: If two polygons overlap, the one last in the list takes precedence
859
860    Coordinates specified in the call are assumed to be relative to the
861    origin (georeference) e.g. used by domain.
862    By specifying the optional argument georeference,
863    all points are made relative.
864
865    FIXME: This should really work with geo_spatial point sets.
866    """
867
868    ##
869    # @brief Create instance of a polygon function.
870    # @param regions A list of (x,y) tuples defining a polygon.
871    # @param default Value or function returning value for points outside poly.
872    # @param geo_reference ??
873    def __init__(self, regions, default=0.0, geo_reference=None):
874        try:
875            len(regions)
876        except:
877            msg = ('Polygon_function takes a list of pairs (polygon, value).'
878                   'Got %s' % str(regions))
879            raise Exception, msg
880
881        T = regions[0]
882
883        if isinstance(T, basestring):
884            msg = ('You passed in a list of text values into polygon_function '
885                   'instead of a list of pairs (polygon, value): "%s"'
886                   % str(T))
887            raise Exception, msg
888
889        try:
890            a = len(T)
891        except:
892            msg = ('Polygon_function takes a list of pairs (polygon, value). '
893                   'Got %s' % str(T))
894            raise Exception, msg
895
896        msg = ('Each entry in regions have two components: (polygon, value). '
897               'I got %s' % str(T))
898        assert a == 2, msg
899
900        if geo_reference is None:
901            from anuga.coordinate_transforms.geo_reference import Geo_reference
902            geo_reference = Geo_reference()
903
904        self.default = default
905
906        # Make points in polygons relative to geo_reference
907        self.regions = []
908        for polygon, value in regions:
909            P = geo_reference.change_points_geo_ref(polygon)
910            self.regions.append((P, value))
911
912    ##
913    # @brief Implement the 'callable' property of Polygon_function.
914    # @param x List of x coordinates of points ot interest.
915    # @param y List of y coordinates of points ot interest.
916    def __call__(self, x, y):
917        x = num.array(x, num.float)
918        y = num.array(y, num.float)
919
920        # x and y must be one-dimensional and same length
921        assert len(x.shape) == 1 and len(y.shape) == 1
922        N = x.shape[0]
923        assert y.shape[0] == N
924
925        points = num.ascontiguousarray(num.concatenate((x[:,num.newaxis],
926                                                        y[:,num.newaxis]),
927                                                       axis=1 ))
928
929        if callable(self.default):
930            z = self.default(x, y)
931        else:
932            z = num.ones(N, num.float) * self.default
933
934        for polygon, value in self.regions:
935            indices = inside_polygon(points, polygon)
936
937            # FIXME: This needs to be vectorised
938            if callable(value):
939                for i in indices:
940                    xx = num.array([x[i]])
941                    yy = num.array([y[i]])
942                    z[i] = value(xx, yy)[0]
943            else:
944                for i in indices:
945                    z[i] = value
946
947        if len(z) == 0:
948            msg = ('Warning: points provided to Polygon function did not fall '
949                   'within its regions in [%.2f, %.2f], y in [%.2f, %.2f]'
950                   % (min(x), max(x), min(y), max(y)))
951            log.critical(msg)
952
953        return z
954
955################################################################################
956# Functions to read and write polygon information
957################################################################################
958
959##
960# @brief Read polygon data from a file.
961# @param filename Path to file containing polygon data.
962# @param delimiter Delimiter to split polygon data with.
963# @return A list of point data from the polygon file.
964def read_polygon(filename, delimiter=','):
965    """Read points assumed to form a polygon.
966
967    There must be exactly two numbers in each line separated by the delimiter.
968    No header.
969    """
970
971    fid = open(filename)
972    lines = fid.readlines()
973    fid.close()
974    polygon = []
975    for line in lines:
976        fields = line.split(delimiter)
977        polygon.append([float(fields[0]), float(fields[1])])
978
979    return polygon
980
981##
982# @brief Write polygon data to a file.
983# @param polygon Polygon points to write to file.
984# @param filename Path to file to write.
985# @note Delimiter is assumed to be a comma.
986def write_polygon(polygon, filename=None):
987    """Write polygon to csv file.
988
989    There will be exactly two numbers, easting and northing, in each line
990    separated by a comma.
991
992    No header.
993    """
994
995    fid = open(filename, 'w')
996    for point in polygon:
997        fid.write('%f, %f\n' % point)
998    fid.close()
999
1000##
1001# @brief Unimplemented.
1002def read_tagged_polygons(filename):
1003    """
1004    """
1005    pass
1006
1007##
1008# @brief Populate given polygon with uniformly distributed points.
1009# @param polygon Polygon to uniformly fill.
1010# @param number_of_points Number of points required in polygon.
1011# @param seed Seed for random number generator.
1012# @param exclude List of polygons inside main where points should be excluded.
1013# @return List of random points inside input polygon.
1014# @note Delimiter is assumed to be a comma.
1015def populate_polygon(polygon, number_of_points, seed=None, exclude=None):
1016    """Populate given polygon with uniformly distributed points.
1017
1018    Input:
1019       polygon - list of vertices of polygon
1020       number_of_points - (optional) number of points
1021       seed - seed for random number generator (default=None)
1022       exclude - list of polygons (inside main polygon) from where points
1023                 should be excluded
1024
1025    Output:
1026       points - list of points inside polygon
1027
1028    Examples:
1029       populate_polygon( [[0,0], [1,0], [1,1], [0,1]], 5 )
1030       will return five randomly selected points inside the unit square
1031    """
1032
1033    from random import uniform, seed as seed_function
1034
1035    seed_function(seed)
1036
1037    points = []
1038
1039    # Find outer extent of polygon
1040    max_x = min_x = polygon[0][0]
1041    max_y = min_y = polygon[0][1]
1042    for point in polygon[1:]:
1043        x = point[0]
1044        if x > max_x: max_x = x
1045        if x < min_x: min_x = x
1046        y = point[1]
1047        if y > max_y: max_y = y
1048        if y < min_y: min_y = y
1049
1050    while len(points) < number_of_points:
1051        x = uniform(min_x, max_x)
1052        y = uniform(min_y, max_y)
1053
1054        append = False
1055        if is_inside_polygon([x,y], polygon):
1056            append = True
1057
1058            #Check exclusions
1059            if exclude is not None:
1060                for ex_poly in exclude:
1061                    if is_inside_polygon([x,y], ex_poly):
1062                        append = False
1063
1064        if append is True:
1065            points.append([x,y])
1066
1067    return points
1068
1069##
1070# @brief Get a point inside a polygon that is close to an edge.
1071# @param polygon List  of vertices of polygon.
1072# @param delta Maximum distance from an edge is delta * sqrt(2).
1073# @return A point that is inside polgon and close to the polygon edge.
1074def point_in_polygon(polygon, delta=1e-8):
1075    """Return a point inside a given polygon which will be close to the
1076    polygon edge.
1077
1078    Input:
1079       polygon - list of vertices of polygon
1080       delta - the square root of 2 * delta is the maximum distance from the
1081       polygon points and the returned point.
1082    Output:
1083       points - a point inside polygon
1084
1085       searches in all diagonals and up and down (not left and right).
1086    """
1087
1088    import exceptions
1089
1090    class Found(exceptions.Exception): pass
1091
1092    polygon = ensure_numeric(polygon)
1093   
1094    point_in = False
1095    while not point_in:
1096        try:
1097            for poly_point in polygon:     # [1:]:
1098                for x_mult in range(-1, 2):
1099                    for y_mult in range(-1, 2):
1100                        x = poly_point[0]
1101                        y = poly_point[1]
1102
1103                        if x == 0:
1104                            x_delta = x_mult * delta
1105                        else:
1106                            x_delta = x + x_mult*x*delta
1107
1108                        if y == 0:
1109                            y_delta = y_mult * delta
1110                        else:
1111                            y_delta = y + y_mult*y*delta
1112
1113                        point = [x_delta, y_delta]
1114
1115                        if is_inside_polygon(point, polygon, closed=False):
1116                            raise Found
1117        except Found:
1118            point_in = True
1119        else:
1120            delta = delta * 0.1
1121
1122    return point
1123
1124##
1125# @brief Calculate approximate number of triangles inside a bounding polygon.
1126# @param interior_regions
1127# @param bounding_poly
1128# @param remainder_res
1129# @return The number of triangles.
1130def number_mesh_triangles(interior_regions, bounding_poly, remainder_res):
1131    """Calculate the approximate number of triangles inside the
1132    bounding polygon and the other interior regions
1133
1134    Polygon areas are converted to square Kms
1135
1136    FIXME: Add tests for this function
1137    """
1138
1139    from anuga.utilities.polygon import polygon_area
1140
1141    # TO DO check if any of the regions fall inside one another
1142
1143    log.critical('-' * 80)
1144    log.critical('Polygon   Max triangle area (m^2)   Total area (km^2) '
1145                 'Estimated #triangles')
1146    log.critical('-' * 80)
1147       
1148    no_triangles = 0.0
1149    area = polygon_area(bounding_poly)
1150
1151    for poly, resolution in interior_regions:
1152        this_area = polygon_area(poly)
1153        this_triangles = this_area/resolution
1154        no_triangles += this_triangles
1155        area -= this_area
1156
1157        log.critical('Interior %s%s%d'
1158                     % (('%.0f' % resolution).ljust(25),
1159                        ('%.2f' % (this_area/1000000)).ljust(19), 
1160                        this_triangles))
1161        #print 'Interior ',
1162        #print ('%.0f' % resolution).ljust(25),
1163        #print ('%.2f' % (this_area/1000000)).ljust(19),
1164        #print '%d' % (this_triangles)
1165
1166    bound_triangles = area/remainder_res
1167    no_triangles += bound_triangles
1168
1169    log.critical('Bounding %s%s%d'
1170                 % (('%.0f' % remainder_res).ljust(25),
1171                    ('%.2f' % (area/1000000)).ljust(19),
1172                    bound_triangles))
1173    #print 'Bounding ',
1174    #print ('%.0f' % remainder_res).ljust(25),
1175    #print ('%.2f' % (area/1000000)).ljust(19),
1176    #print '%d' % (bound_triangles)
1177
1178    total_number_of_triangles = no_triangles/0.7
1179
1180    log.critical('Estimated total number of triangles: %d'
1181                 % total_number_of_triangles)
1182    log.critical('Note: This is generally about 20%% '
1183                 'less than the final amount')
1184
1185    return int(total_number_of_triangles)
1186
1187##
1188# @brief Reduce number of points in polygon by the specified factor.
1189# @param polygon The polygon to reduce.
1190# @param factor The factor to reduce polygon points by (default 10).
1191# @return The reduced polygon points list.
1192# @note The extrema of both axes are preserved.
1193def decimate_polygon(polygon, factor=10):
1194    """Reduce number of points in polygon by the specified
1195    factor (default=10, hence the name of the function) such that
1196    the extrema in both axes are preserved.
1197
1198    Return reduced polygon
1199    """
1200
1201    # FIXME(Ole): This doesn't work at present,
1202    # but it isn't critical either
1203
1204    # Find outer extent of polygon
1205    num_polygon = ensure_numeric(polygon)
1206    max_x = max(num_polygon[:,0])
1207    max_y = max(num_polygon[:,1])
1208    min_x = min(num_polygon[:,0])
1209    min_y = min(num_polygon[:,1])
1210
1211    # Keep only some points making sure extrema are kept
1212    reduced_polygon = []
1213    for i, point in enumerate(polygon):
1214        x = point[0]
1215        y = point[1]
1216        if x in [min_x, max_x] and y in [min_y, max_y]:
1217            # Keep
1218            reduced_polygon.append(point)
1219        else:
1220            if len(reduced_polygon)*factor < i:
1221                reduced_polygon.append(point)
1222
1223    return reduced_polygon
1224
1225##
1226# @brief Interpolate linearly from polyline nodes to midpoints of triangles.
1227# @param data The data on the polyline nodes.
1228# @param polyline_nodes ??
1229# @param gauge_neighbour_id ??  FIXME(Ole): I want to get rid of this
1230# @param point_coordinates ??
1231# @param verbose True if this function is to be verbose.
1232def interpolate_polyline(data,
1233                         polyline_nodes,
1234                         gauge_neighbour_id,
1235                         interpolation_points=None,
1236                         rtol=1.0e-6,
1237                         atol=1.0e-8,
1238                         verbose=False):
1239    """Interpolate linearly between values data on polyline nodes
1240    of a polyline to list of interpolation points.
1241
1242    data is the data on the polyline nodes.
1243
1244    Inputs:
1245      data: Vector or array of data at the polyline nodes.
1246      polyline_nodes: Location of nodes where data is available.
1247      gauge_neighbour_id: ?
1248      interpolation_points: Interpolate polyline data to these positions.
1249          List of coordinate pairs [x, y] of
1250          data points or an nx2 numeric array or a Geospatial_data object
1251      rtol, atol: Used to determine whether a point is on the polyline or not.
1252                  See point_on_line.
1253
1254    Output:
1255      Interpolated values at interpolation points
1256    """
1257
1258    if isinstance(interpolation_points, Geospatial_data):
1259        interpolation_points = interpolation_points.\
1260                                    get_data_points(absolute=True)
1261
1262    interpolated_values = num.zeros(len(interpolation_points), num.float)
1263
1264    data = ensure_numeric(data, num.float)
1265    polyline_nodes = ensure_numeric(polyline_nodes, num.float)
1266    interpolation_points = ensure_numeric(interpolation_points, num.float)
1267    gauge_neighbour_id = ensure_numeric(gauge_neighbour_id, num.int)
1268
1269    n = polyline_nodes.shape[0]    # Number of nodes in polyline
1270
1271    # Input sanity check
1272    msg = 'interpolation_points are not given (interpolate.py)'
1273    assert interpolation_points is not None, msg
1274
1275    msg = 'function value must be specified at every interpolation node'
1276    assert data.shape[0] == polyline_nodes.shape[0], msg
1277
1278    msg = 'Must define function value at one or more nodes'
1279    assert data.shape[0] > 0, msg
1280
1281    if n == 1:
1282        msg = 'Polyline contained only one point. I need more. ' + str(data)
1283        raise Exception, msg
1284    elif n > 1:
1285        _interpolate_polyline(data,
1286                              polyline_nodes,
1287                              gauge_neighbour_id,
1288                              interpolation_points,
1289                              interpolated_values,
1290                              rtol,
1291                              atol)
1292
1293    return interpolated_values
1294
1295##
1296# @brief
1297# @param data
1298# @param polyline_nodes
1299# @param gauge_neighbour_id
1300# @param interpolation_points
1301# @param interpolated_values
1302# @param rtol
1303# @param atol
1304# @return
1305# @note OBSOLETED BY C-EXTENSION
1306def _interpolate_polyline(data,
1307                          polyline_nodes,
1308                          gauge_neighbour_id,
1309                          interpolation_points,
1310                          interpolated_values,
1311                          rtol=1.0e-6,
1312                          atol=1.0e-8):
1313    """Auxiliary function used by interpolate_polyline
1314
1315    NOTE: OBSOLETED BY C-EXTENSION
1316    """
1317
1318    number_of_nodes = len(polyline_nodes)
1319    number_of_points = len(interpolation_points)
1320
1321    for j in range(number_of_nodes):
1322        neighbour_id = gauge_neighbour_id[j]
1323
1324        # FIXME(Ole): I am convinced that gauge_neighbour_id can be discarded,
1325        #             but need to check with John J.
1326        # Keep it for now (17 Jan 2009)
1327        # When gone, we can simply interpolate between neighbouring nodes,
1328        # i.e. neighbour_id = j+1.
1329        # and the test below becomes something like: if j < number_of_nodes...
1330
1331        if neighbour_id >= 0:
1332            x0, y0 = polyline_nodes[j,:]
1333            x1, y1 = polyline_nodes[neighbour_id,:]
1334
1335            segment_len = sqrt((x1-x0)**2 + (y1-y0)**2)
1336            segment_delta = data[neighbour_id] - data[j]
1337            slope = segment_delta/segment_len
1338
1339            for i in range(number_of_points):
1340                x, y = interpolation_points[i,:]
1341                if point_on_line([x, y], [[x0, y0], [x1, y1]],
1342                                 rtol=rtol, atol=atol):
1343                    dist = sqrt((x-x0)**2 + (y-y0)**2)
1344                    interpolated_values[i] = slope*dist + data[j]
1345
1346
1347################################################################################
1348# Initialise module
1349################################################################################
1350
1351from anuga.utilities import compile
1352if compile.can_use_C_extension('polygon_ext.c'):
1353    # Underlying C implementations can be accessed
1354    from polygon_ext import _point_on_line
1355    from polygon_ext import _separate_points_by_polygon
1356    from polygon_ext import _interpolate_polyline   
1357    from polygon_ext import _is_inside_triangle       
1358    #from polygon_ext import _intersection
1359
1360else:
1361    msg = 'C implementations could not be accessed by %s.\n ' %__file__
1362    msg += 'Make sure compile_all.py has been run as described in '
1363    msg += 'the ANUGA installation guide.'
1364    raise Exception, msg
1365
1366
1367if __name__ == "__main__":
1368    pass
Note: See TracBrowser for help on using the repository browser.