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

Last change on this file since 5518 was 5403, checked in by ole, 17 years ago

Added tolerances for point_on_line and changed interpolate_polyline to use a general relative tolerance.

File size: 28.2 KB
Line 
1#!/usr/bin/env python
2"""Polygon manipulations
3
4"""
5
6
7#try:
8#    from scipy import Float, Int, zeros, ones, array, concatenate, reshape, dot
9#except:
10#    #print 'Could not find scipy - using Numeric'
11
12from Numeric import Float, Int, zeros, ones, array, concatenate, reshape, dot, allclose
13
14
15from math import sqrt
16from anuga.utilities.numerical_tools import ensure_numeric
17from anuga.geospatial_data.geospatial_data import ensure_absolute
18
19
20def point_on_line(point, line, rtol=0.0, atol=0.0):
21    """Determine whether a point is on a line segment
22
23    Input:
24        point is given by [x, y]
25        line is given by [x0, y0], [x1, y1]] or
26        the equivalent 2x2 Numeric array with each row corresponding to a point.
27
28    Output:
29
30    Note: Line can be degenerate and function still works to discern coinciding points from non-coinciding.
31    """
32
33    # FIXME(Ole): Perhaps make defaults as in allclose: rtol=1.0e-5, atol=1.0e-8
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
48
49def intersection(line0, line1):
50    """Returns intersecting point between two line segments or None
51    (if parallel or no intersection is found).
52
53    However, if parallel lines coincide partly (i.e. shara a common segment,
54    the line segment where lines coincide is returned
55   
56
57    Inputs:
58        line0, line1: Each defined by two end points as in: [[x0, y0], [x1, y1]]
59                      A line can also be a 2x2 numeric array with each row
60                      corresponding to a point.
61
62
63    Output:
64        status, value
65
66        where status is interpreted as follows
67       
68        status == 0: no intersection with value set to None
69        status == 1: One intersection point found and returned in value as [x,y]
70        status == 2: Coinciding line segment found. Value taks the form [[x0,y0], [x1,y1]]
71        status == 3: Lines would coincide but only if extended. Value set to None
72        status == 4: Lines are parallel with a fixed distance apart. Value set to None.
73   
74    """
75
76    # FIXME (Ole): Write this in C
77
78    line0 = ensure_numeric(line0, Float)
79    line1 = ensure_numeric(line1, Float)   
80
81    x0 = line0[0,0]; y0 = line0[0,1]
82    x1 = line0[1,0]; y1 = line0[1,1]
83
84    x2 = line1[0,0]; y2 = line1[0,1]
85    x3 = line1[1,0]; y3 = line1[1,1]
86
87    denom = (y3-y2)*(x1-x0) - (x3-x2)*(y1-y0)
88    u0 = (x3-x2)*(y0-y2) - (y3-y2)*(x0-x2)
89    u1 = (x2-x0)*(y1-y0) - (y2-y0)*(x1-x0)
90       
91    if allclose(denom, 0.0):
92        # Lines are parallel - check if they coincide on a shared a segment
93
94        if allclose( [u0, u1], 0.0 ):
95            # We now know that the lines if continued coincide
96            # The remaining check will establish if the finite lines share a segment
97
98            line0_starts_on_line1 = line0_ends_on_line1 =\
99            line1_starts_on_line0 = line1_ends_on_line0 = False
100               
101            if point_on_line([x0, y0], line1):
102                line0_starts_on_line1 = True
103
104            if point_on_line([x1, y1], line1):
105                line0_ends_on_line1 = True
106 
107            if point_on_line([x2, y2], line0):
108                line1_starts_on_line0 = True
109
110            if point_on_line([x3, y3], line0):
111                line1_ends_on_line0 = True                               
112
113            if not(line0_starts_on_line1 or line0_ends_on_line1\
114               or line1_starts_on_line0 or line1_ends_on_line0):
115                # Lines are parallel and would coincide if extended, but not as they are.
116                return 3, None
117
118
119            # One line fully included in the other. Use direction of included line
120            if line0_starts_on_line1 and line0_ends_on_line1:
121                # Shared segment is line0 fully included in line1
122                segment = array([[x0, y0], [x1, y1]])               
123
124            if line1_starts_on_line0 and line1_ends_on_line0:
125                # Shared segment is line1 fully included in line0
126                segment = array([[x2, y2], [x3, y3]])
127           
128
129            # Overlap with lines are oriented the same way
130            if line0_starts_on_line1 and line1_ends_on_line0:
131                # Shared segment from line0 start to line 1 end
132                segment = array([[x0, y0], [x3, y3]])
133
134            if line1_starts_on_line0 and line0_ends_on_line1:
135                # Shared segment from line1 start to line 0 end
136                segment = array([[x2, y2], [x1, y1]])                               
137
138
139            # Overlap in opposite directions - use direction of line0
140            if line0_starts_on_line1 and line1_starts_on_line0:
141                # Shared segment from line0 start to line 1 end
142                segment = array([[x0, y0], [x2, y2]])
143
144            if line0_ends_on_line1 and line1_ends_on_line0:
145                # Shared segment from line0 start to line 1 end
146                segment = array([[x3, y3], [x1, y1]])               
147
148               
149            return 2, segment
150        else:
151            # Lines are parallel but they don't coincide
152            return 4, None #FIXME (Ole): Add distance here instead of None
153           
154    else:
155        # Lines are not parallel or coinciding
156        u0 = u0/denom
157        u1 = u1/denom       
158
159        x = x0 + u0*(x1-x0)
160        y = y0 + u0*(y1-y0)
161
162        # Sanity check - can be removed to speed up if needed
163        assert allclose(x, x2 + u1*(x3-x2))
164        assert allclose(y, y2 + u1*(y3-y2))       
165
166        # Check if point found lies within given line segments
167        if 0.0 <= u0 <= 1.0 and 0.0 <= u1 <= 1.0: 
168            # We have intersection
169
170            return 1, array([x, y])
171        else:
172            # No intersection
173            return 0, None
174
175
176def NEW_C_intersection(line0, line1):
177    #FIXME(Ole): To write in C
178    """Returns intersecting point between two line segments or None
179    (if parallel or no intersection is found).
180
181    However, if parallel lines coincide partly (i.e. shara a common segment,
182    the line segment where lines coincide is returned
183   
184
185    Inputs:
186        line0, line1: Each defined by two end points as in: [[x0, y0], [x1, y1]]
187                      A line can also be a 2x2 numeric array with each row
188                      corresponding to a point.
189
190
191    Output:
192        status, value
193
194        where status is interpreted as follows
195       
196        status == 0: no intersection with value set to None
197        status == 1: One intersection point found and returned in value as [x,y]
198        status == 2: Coinciding line segment found. Value taks the form [[x0,y0], [x1,y1]]
199        status == 3: Lines would coincide but only if extended. Value set to None
200        status == 4: Lines are parallel with a fixed distance apart. Value set to None.
201   
202    """
203
204
205    line0 = ensure_numeric(line0, Float)
206    line1 = ensure_numeric(line1, Float)   
207
208    status, value = _intersection(line0[0,0], line0[0,1],
209                                  line0[1,0], line0[1,1],
210                                  line1[0,0], line1[0,1],
211                                  line1[1,0], line1[1,1])
212
213    return status, value
214
215
216
217
218def is_inside_polygon(point, polygon, closed=True, verbose=False):
219    """Determine if one point is inside a polygon
220
221    See inside_polygon for more details
222    """
223
224    indices = inside_polygon(point, polygon, closed, verbose)
225
226    if indices.shape[0] == 1:
227        return True
228    elif indices.shape[0] == 0:
229        return False
230    else:
231        msg = 'is_inside_polygon must be invoked with one point only'
232        raise msg
233   
234
235def inside_polygon(points, polygon, closed=True, verbose=False):
236    """Determine points inside a polygon
237
238       Functions inside_polygon and outside_polygon have been defined in
239       terms af separate_by_polygon which will put all inside indices in
240       the first part of the indices array and outside indices in the last
241
242       See separate_points_by_polygon for documentation
243
244       points and polygon can be a geospatial instance,
245       a list or a numeric array
246    """
247
248    #if verbose: print 'Checking input to inside_polygon'
249
250    try:
251        points = ensure_absolute(points)
252    except NameError, e:
253        raise NameError, e
254    except:
255        # If this fails it is going to be because the points can't be
256        # converted to a numeric array.
257        msg = 'Points could not be converted to Numeric array' 
258        raise msg
259
260    try:
261        polygon = ensure_absolute(polygon)
262    except NameError, e:
263        raise NameError, e
264    except:
265        # If this fails it is going to be because the points can't be
266        # converted to a numeric array.
267        msg = 'Polygon %s could not be converted to Numeric array' %(str(polygon))
268        raise msg
269
270    if len(points.shape) == 1:
271        # Only one point was passed in. Convert to array of points
272        points = reshape(points, (1,2))
273
274    indices, count = separate_points_by_polygon(points, polygon,
275                                                closed=closed,
276                                                verbose=verbose)
277
278    # Return indices of points inside polygon
279    return indices[:count]
280
281
282
283def is_outside_polygon(point, polygon, closed=True, verbose=False,
284                       points_geo_ref=None, polygon_geo_ref=None):
285    """Determine if one point is outside a polygon
286
287    See outside_polygon for more details
288    """
289
290    indices = outside_polygon(point, polygon, closed, verbose)
291                              #points_geo_ref, polygon_geo_ref)
292
293    if indices.shape[0] == 1:
294        return True
295    elif indices.shape[0] == 0:
296        return False
297    else:
298        msg = 'is_outside_polygon must be invoked with one point only'
299        raise msg
300   
301
302def outside_polygon(points, polygon, closed = True, verbose = False):
303    """Determine points outside a polygon
304
305       Functions inside_polygon and outside_polygon have been defined in
306       terms af separate_by_polygon which will put all inside indices in
307       the first part of the indices array and outside indices in the last
308
309       See separate_points_by_polygon for documentation
310    """
311
312    #if verbose: print 'Checking input to outside_polygon'
313    try:
314        points = ensure_numeric(points, Float)
315    except NameError, e:
316        raise NameError, e
317    except:
318        msg = 'Points could not be converted to Numeric array'
319        raise msg
320
321    try:
322        polygon = ensure_numeric(polygon, Float)
323    except NameError, e:
324        raise NameError, e
325    except:
326        msg = 'Polygon could not be converted to Numeric array'
327        raise msg
328
329
330    if len(points.shape) == 1:
331        # Only one point was passed in. Convert to array of points
332        points = reshape(points, (1,2))
333
334    indices, count = separate_points_by_polygon(points, polygon,
335                                                closed=closed,
336                                                verbose=verbose)
337
338    # Return indices of points outside polygon
339    if count == len(indices):
340        # No points are outside
341        return array([])
342    else:
343        return indices[count:][::-1]  #return reversed
344       
345
346def in_and_outside_polygon(points, polygon, closed = True, verbose = False):
347    """Determine points inside and outside a polygon
348
349       See separate_points_by_polygon for documentation
350
351       Returns an array of points inside and an array of points outside the polygon
352    """
353
354    #if verbose: print 'Checking input to outside_polygon'
355    try:
356        points = ensure_numeric(points, Float)
357    except NameError, e:
358        raise NameError, e
359    except:
360        msg = 'Points could not be converted to Numeric array'
361        raise msg
362
363    try:
364        polygon = ensure_numeric(polygon, Float)
365    except NameError, e:
366        raise NameError, e
367    except:
368        msg = 'Polygon could not be converted to Numeric array'
369        raise msg
370
371    if len(points.shape) == 1:
372        # Only one point was passed in. Convert to array of points
373        points = reshape(points, (1,2))
374
375
376    indices, count = separate_points_by_polygon(points, polygon,
377                                                closed=closed,
378                                                verbose=verbose)
379   
380    # Returns indices of points inside and indices of points outside
381    # the polygon
382
383    if count == len(indices):
384        # No points are outside
385        return indices[:count],[]
386    else:
387        return  indices[:count], indices[count:][::-1]  #return reversed
388
389
390def separate_points_by_polygon(points, polygon,
391                               closed = True, verbose = False):
392    """Determine whether points are inside or outside a polygon
393
394    Input:
395       points - Tuple of (x, y) coordinates, or list of tuples
396       polygon - list of vertices of polygon
397       closed - (optional) determine whether points on boundary should be
398       regarded as belonging to the polygon (closed = True)
399       or not (closed = False)
400
401    Outputs:
402       indices: array of same length as points with indices of points falling
403       inside the polygon listed from the beginning and indices of points
404       falling outside listed from the end.
405
406       count: count of points falling inside the polygon
407
408       The indices of points inside are obtained as indices[:count]
409       The indices of points outside are obtained as indices[count:]
410
411
412    Examples:
413       U = [[0,0], [1,0], [1,1], [0,1]] #Unit square
414
415       separate_points_by_polygon( [[0.5, 0.5], [1, -0.5], [0.3, 0.2]], U)
416       will return the indices [0, 2, 1] and count == 2 as only the first
417       and the last point are inside the unit square
418
419    Remarks:
420       The vertices may be listed clockwise or counterclockwise and
421       the first point may optionally be repeated.
422       Polygons do not need to be convex.
423       Polygons can have holes in them and points inside a hole is
424       regarded as being outside the polygon.
425
426    Algorithm is based on work by Darel Finley,
427    http://www.alienryderflex.com/polygon/
428
429    Uses underlying C-implementation in polygon_ext.c
430    """
431
432
433    #if verbose: print 'Checking input to separate_points_by_polygon'
434
435
436    #Input checks
437
438    assert isinstance(closed, bool), 'Keyword argument "closed" must be boolean'
439    assert isinstance(verbose, bool), 'Keyword argument "verbose" must be boolean'
440
441
442    try:
443        points = ensure_numeric(points, Float)
444    except NameError, e:
445        raise NameError, e
446    except:
447        msg = 'Points could not be converted to Numeric array'
448        raise msg
449
450    #if verbose: print 'Checking input to separate_points_by_polygon 2'
451    try:
452        polygon = ensure_numeric(polygon, Float)
453    except NameError, e:
454        raise NameError, e
455    except:
456        msg = 'Polygon could not be converted to Numeric array'
457        raise msg
458
459    #if verbose: print 'check'
460
461    assert len(polygon.shape) == 2,\
462       'Polygon array must be a 2d array of vertices'
463
464    assert polygon.shape[1] == 2,\
465       'Polygon array must have two columns'
466
467    assert len(points.shape) == 2,\
468       'Points array must be a 2d array'
469
470    assert points.shape[1] == 2,\
471       'Points array must have two columns'
472
473    N = polygon.shape[0] #Number of vertices in polygon
474    M = points.shape[0]  #Number of points
475
476
477    indices = zeros( M, Int )
478
479    count = _separate_points_by_polygon(points, polygon, indices,
480                                        int(closed), int(verbose))
481
482    if verbose: print 'Found %d points (out of %d) inside polygon'\
483       %(count, M)
484    return indices, count
485
486
487def polygon_area(polygon):
488    """ Determin area of arbitrary polygon
489    Reference
490    http://mathworld.wolfram.com/PolygonArea.html
491    """
492   
493    n = len(polygon)
494    poly_area = 0.0
495
496    for i in range(n):
497        pti = polygon[i]
498        if i == n-1:
499            pt1 = polygon[0]
500        else:
501            pt1 = polygon[i+1]
502        xi = pti[0]
503        yi1 = pt1[1]
504        xi1 = pt1[0]
505        yi = pti[1]
506        poly_area += xi*yi1 - xi1*yi
507       
508    return abs(poly_area/2)
509
510def plot_polygons(polygons_points, style=None, 
511                  figname=None, label=None, verbose=False):
512   
513    """ Take list of polygons and plot.
514
515    Inputs:
516
517    polygons         - list of polygons
518
519    style            - style list corresponding to each polygon
520                     - for a polygon, use 'line'
521                     - for points falling outside a polygon, use 'outside'
522                       
523    figname          - name to save figure to
524
525    label            - title for plot
526
527    Outputs:
528
529    - list of min and max of x and y coordinates
530    - plot of polygons
531    """ 
532
533    from pylab import ion, hold, plot, axis, figure, legend, savefig, xlabel, ylabel, title, close, title
534
535    assert type(polygons_points) == list,\
536               'input must be a list of polygons and/or points'
537               
538    ion()
539    hold(True)
540
541    minx = 1e10
542    maxx = 0.0
543    miny = 1e10
544    maxy = 0.0
545
546    if label is None: label = ''
547
548    n = len(polygons_points)
549    colour = []
550    if style is None:
551        style_type = 'line' 
552        style = []
553        for i in range(n):
554            style.append(style_type)
555            colour.append('b-')
556    else:
557        for s in style:
558            if s == 'line': colour.append('b-')           
559            if s == 'outside': colour.append('r.')
560            if s <> 'line':
561                if s <> 'outside':
562                    colour.append('g.')
563           
564    for i, item in enumerate(polygons_points):
565        x, y = poly_xy(item)
566        if min(x) < minx: minx = min(x)
567        if max(x) > maxx: maxx = max(x)
568        if min(y) < miny: miny = min(y)
569        if max(y) > maxy: maxy = max(y)
570        plot(x,y,colour[i])
571        xlabel('x')
572        ylabel('y')
573        title(label)
574
575    #raw_input('wait 1')
576    #FIXME(Ole): This makes for some strange scalings sometimes.
577    #if minx <> 0:
578    #    axis([minx*0.9,maxx*1.1,miny*0.9,maxy*1.1])
579    #else:
580    #    if miny == 0:
581    #        axis([-maxx*.01,maxx*1.1,-maxy*0.01,maxy*1.1])
582    #    else:
583    #        axis([-maxx*.01,maxx*1.1,miny*0.9,maxy*1.1])
584
585    if figname is not None:
586        savefig(figname)
587    else:
588        savefig('test_image')
589
590    close('all')
591
592    vec = [minx,maxx,miny,maxy]
593
594    return vec
595
596def poly_xy(polygon, verbose=False):
597    """ this is used within plot_polygons so need to duplicate
598        the first point so can have closed polygon in plot
599    """
600
601    #if verbose: print 'Checking input to poly_xy'
602
603    try:
604        polygon = ensure_numeric(polygon, Float)
605    except NameError, e:
606        raise NameError, e
607    except:
608        msg = 'Polygon %s could not be converted to Numeric array' %(str(polygon))
609        raise msg
610
611    x = polygon[:,0]
612    y = polygon[:,1]
613    x = concatenate((x, [polygon[0,0]]), axis = 0)
614    y = concatenate((y, [polygon[0,1]]), axis = 0)
615   
616    return x, y
617   
618#    x = []
619#    y = []
620#    n = len(poly)
621#    firstpt = poly[0]
622#    for i in range(n):
623#        thispt = poly[i]
624#        x.append(thispt[0])
625#        y.append(thispt[1])
626
627#    x.append(firstpt[0])
628#    y.append(firstpt[1])
629   
630#    return x, y
631
632class Polygon_function:
633    """Create callable object f: x,y -> z, where a,y,z are vectors and
634    where f will return different values depending on whether x,y belongs
635    to specified polygons.
636
637    To instantiate:
638
639       Polygon_function(polygons)
640
641    where polygons is a list of tuples of the form
642
643      [ (P0, v0), (P1, v1), ...]
644
645      with Pi being lists of vertices defining polygons and vi either
646      constants or functions of x,y to be applied to points with the polygon.
647
648    The function takes an optional argument, default which is the value
649    (or function) to used for points not belonging to any polygon.
650    For example:
651
652       Polygon_function(polygons, default = 0.03)
653
654    If omitted the default value will be 0.0
655
656    Note: If two polygons overlap, the one last in the list takes precedence
657
658    Coordinates specified in the call are assumed to be relative to the
659    origin (georeference) e.g. used by domain.
660    By specifying the optional argument georeference,
661    all points are made relative.
662
663    FIXME: This should really work with geo_spatial point sets.
664    """
665
666    def __init__(self, regions, default=0.0, geo_reference=None):
667
668        try:
669            len(regions)
670        except:
671            msg = 'Polygon_function takes a list of pairs (polygon, value).'
672            msg += 'Got %s' %polygons
673            raise msg
674
675
676        T = regions[0]
677        try:
678            a = len(T)
679        except:
680            msg = 'Polygon_function takes a list of pairs (polygon, value).'
681            msg += 'Got %s' %polygons
682            raise msg
683
684        assert a == 2, 'Must have two component each: %s' %T
685
686
687        if geo_reference is None:
688            from anuga.coordinate_transforms.geo_reference import Geo_reference
689            geo_reference = Geo_reference()
690
691
692        self.default = default
693
694        # Make points in polygons relative to geo_reference
695        self.regions = []
696        for polygon, value in regions:
697            P = geo_reference.change_points_geo_ref(polygon)
698            self.regions.append( (P, value) )
699
700
701
702
703    def __call__(self, x, y):
704        x = array(x).astype(Float)
705        y = array(y).astype(Float)
706
707        N = len(x)
708        assert len(y) == N
709
710        points = concatenate( (reshape(x, (N, 1)),
711                               reshape(y, (N, 1))), axis=1 )
712
713        if callable(self.default):
714            z = self.default(x,y)
715        else:
716            z = ones(N, Float) * self.default
717
718        for polygon, value in self.regions:
719            indices = inside_polygon(points, polygon)
720
721            # FIXME: This needs to be vectorised
722            if callable(value):
723                for i in indices:
724                    xx = array([x[i]])
725                    yy = array([y[i]])
726                    z[i] = value(xx, yy)[0]
727            else:
728                for i in indices:
729                    z[i] = value
730
731        return z
732
733
734def read_polygon(filename, split=','):
735    """Read points assumed to form a polygon.
736       There must be exactly two numbers in each line separated by a comma.
737       No header.
738    """
739
740    #Get polygon
741    fid = open(filename)
742    lines = fid.readlines()
743    fid.close()
744    polygon = []
745    for line in lines:
746        fields = line.split(split)
747        polygon.append( [float(fields[0]), float(fields[1])] )
748
749    return polygon
750
751
752def write_polygon(polygon, filename=None):
753    """Write polygon to csv file.
754       There will be exactly two numbers, easting and northing,
755       in each line separated by a comma.
756       
757       No header.   
758    """
759
760    fid = open(filename, 'w')
761    for point in polygon:
762        fid.write('%f, %f\n' %point)
763    fid.close()
764   
765
766def populate_polygon(polygon, number_of_points, seed=None, exclude=None):
767    """Populate given polygon with uniformly distributed points.
768
769    Input:
770       polygon - list of vertices of polygon
771       number_of_points - (optional) number of points
772       seed - seed for random number generator (default=None)
773       exclude - list of polygons (inside main polygon) from where points should be excluded
774
775    Output:
776       points - list of points inside polygon
777
778    Examples:
779       populate_polygon( [[0,0], [1,0], [1,1], [0,1]], 5 )
780       will return five randomly selected points inside the unit square
781    """
782
783    from random import uniform, seed as seed_function
784
785    seed_function(seed)
786
787    points = []
788
789    # Find outer extent of polygon
790    max_x = min_x = polygon[0][0]
791    max_y = min_y = polygon[0][1]
792    for point in polygon[1:]:
793        x = point[0]
794        if x > max_x: max_x = x
795        if x < min_x: min_x = x
796        y = point[1]
797        if y > max_y: max_y = y
798        if y < min_y: min_y = y
799
800
801    while len(points) < number_of_points:
802        x = uniform(min_x, max_x)
803        y = uniform(min_y, max_y)
804
805        append = False
806        if is_inside_polygon([x,y], polygon):
807
808            append = True
809
810            #Check exclusions
811            if exclude is not None:
812                for ex_poly in exclude:
813                    if is_inside_polygon([x,y], ex_poly):
814                        append = False
815
816
817        if append is True:
818            points.append([x,y])
819
820    return points
821
822
823def point_in_polygon(polygon, delta=1e-8):
824    """Return a point inside a given polygon which will be close to the
825    polygon edge.
826
827    Input:
828       polygon - list of vertices of polygon
829       delta - the square root of 2 * delta is the maximum distance from the
830       polygon points and the returned point.
831    Output:
832       points - a point inside polygon
833
834       searches in all diagonals and up and down (not left and right)
835    """
836    import exceptions
837    class Found(exceptions.Exception): pass
838
839    point_in = False
840    while not point_in:
841        try:
842            for poly_point in polygon: #[1:]:
843                for x_mult in range (-1,2):
844                    for y_mult in range (-1,2):
845                        x = poly_point[0]
846                        y = poly_point[1]
847                        if x == 0:
848                            x_delta = x_mult*delta
849                        else:
850                            x_delta = x+x_mult*x*delta
851
852                        if y == 0:
853                            y_delta = y_mult*delta
854                        else:
855                            y_delta = y+y_mult*y*delta
856
857                        point = [x_delta, y_delta]
858                        #print "point",point
859                        if is_inside_polygon(point, polygon, closed=False):
860                            raise Found
861        except Found:
862            point_in = True
863        else:
864            delta = delta*0.1
865    return point
866
867
868def number_mesh_triangles(interior_regions, bounding_poly, remainder_res):
869    """Calculate the approximate number of triangles inside the
870    bounding polygon and the other interior regions
871
872    Polygon areas are converted to square Kms
873
874    FIXME: Add tests for this function
875    """
876   
877    from anuga.utilities.polygon import polygon_area
878
879
880    # TO DO check if any of the regions fall inside one another
881
882    print '----------------------------------------------------------------------------'
883    print 'Polygon   Max triangle area (m^2)   Total area (km^2)   Estimated #triangles'
884    print '----------------------------------------------------------------------------'   
885       
886    no_triangles = 0.0
887    area = polygon_area(bounding_poly)
888   
889    for poly, resolution in interior_regions:
890        this_area = polygon_area(poly)
891        this_triangles = this_area/resolution
892        no_triangles += this_triangles
893        area -= this_area
894       
895        print 'Interior ',
896        print ('%.0f' %resolution).ljust(25),
897        print ('%.2f' %(this_area/1000000)).ljust(19),
898        print '%d' %(this_triangles)
899       
900    bound_triangles = area/remainder_res
901    no_triangles += bound_triangles
902
903    print 'Bounding ',
904    print ('%.0f' %remainder_res).ljust(25),
905    print ('%.2f' %(area/1000000)).ljust(19),
906    print '%d' %(bound_triangles)   
907
908    total_number_of_triangles = no_triangles/0.7
909
910    print 'Estimated total number of triangles: %d' %total_number_of_triangles
911    print 'Note: This is generally about 20% less than the final amount'   
912
913    return int(total_number_of_triangles)
914
915
916def decimate_polygon(polygon, factor=10):
917    """Reduce number of points in polygon by the specified
918    factor (default=10, hence the name of the function) such that
919    the extrema in both axes are preserved.
920
921    Return reduced polygon
922    """
923
924    # FIXME(Ole): This doesn't work at present,
925    # but it isn't critical either
926
927    # Find outer extent of polygon
928    num_polygon = ensure_numeric(polygon)
929    max_x = max(num_polygon[:,0])
930    max_y = max(num_polygon[:,1])
931    min_x = min(num_polygon[:,0])
932    min_y = min(num_polygon[:,1])       
933
934    # Keep only some points making sure extrema are kept
935    reduced_polygon = []   
936    for i, point in enumerate(polygon):
937        x = point[0]
938        y = point[1]       
939        if x in [min_x, max_x] and y in [min_y, max_y]:
940            # Keep
941            reduced_polygon.append(point)
942        else:
943            if len(reduced_polygon)*factor < i:
944                reduced_polygon.append(point)               
945
946    return reduced_polygon
947
948##############################################
949#Initialise module
950
951from anuga.utilities.compile import can_use_C_extension
952if can_use_C_extension('polygon_ext.c'):
953    # Underlying C implementations can be accessed
954    from polygon_ext import _point_on_line
955    from polygon_ext import _separate_points_by_polygon
956    #from polygon_ext import _intersection
957
958else:
959    msg = 'C implementations could not be accessed by %s.\n ' %__file__
960    msg += 'Make sure compile_all.py has been run as described in '
961    msg += 'the ANUGA installation guide.'
962    raise Exception, msg
963
964
965if __name__ == "__main__":
966    pass
Note: See TracBrowser for help on using the repository browser.