source: trunk/anuga_core/source/anuga/utilities/polygon_ext.c @ 8059

Last change on this file since 8059 was 8055, checked in by habili, 13 years ago

Changed function names from "polyline" to "line"

File size: 29.6 KB
Line 
1// Python - C extension for polygon module.
2//
3// To compile (Python2.3):
4//  gcc -c polygon_ext.c -I/usr/include/python2.3 -o polygon_ext.o -Wall -O
5//  gcc -shared polygon_ext.o  -o polygon_ext.so
6//
7// See the module polygon.py
8//
9//
10// Ole Nielsen, GA 2004
11//
12// NOTE: We use long* instead of int* for numeric arrays as this will work both
13//       for 64 as well as 32 bit systems
14
15
16#include "Python.h"
17#include "numpy/arrayobject.h"
18#include "math.h"
19
20#include "util_ext.h"
21
22#define YES 1
23#define NO 0
24
25
26double dist(double x,
27            double y) {
28 
29  return sqrt(x*x + y*y);
30}
31
32
33int __point_on_line(double x, double y,
34                    double x0, double y0,
35                    double x1, double y1,
36                    double rtol,
37                    double atol) {
38  /*Determine whether a point is on a line segment
39
40    Input: x, y, x0, x0, x1, y1: where
41        point is given by x, y
42        line is given by (x0, y0) and (x1, y1)
43
44  */
45
46  double a0, a1, a_normal0, a_normal1, b0, b1, len_a, len_b;
47  double nominator, denominator;
48  int is_parallel;
49
50  a0 = x - x0;
51  a1 = y - y0;
52
53  a_normal0 = a1;
54  a_normal1 = -a0;
55
56  b0 = x1 - x0;
57  b1 = y1 - y0;
58
59  nominator = fabs(a_normal0*b0 + a_normal1*b1);
60  denominator = b0*b0 + b1*b1;
61 
62  // Determine if line is parallel to point vector up to a tolerance
63  is_parallel = 0;
64  if (denominator == 0.0) {
65    // Use absolute tolerance
66    if (nominator <= atol) {
67      is_parallel = 1;
68    }
69  } else {
70    // Denominator is positive - use relative tolerance
71    if (nominator/denominator <= rtol) {
72      is_parallel = 1;
73    }   
74  }
75   
76  if (is_parallel) {
77    // Point is somewhere on the infinite extension of the line
78    // subject to specified absolute tolerance
79
80    len_a = dist(a0, a1); //sqrt(a0*a0 + a1*a1);
81    len_b = dist(b0, b1); //sqrt(b0*b0 + b1*b1);
82
83    if (a0*b0 + a1*b1 >= 0 && len_a <= len_b) {
84      return 1;
85    } else {
86      return 0;
87    }
88  } else {
89    return 0;
90  }
91}
92
93
94//  public domain function by Darel Rex Finley, 2006
95//  http://www.alienryderflex.com/intersect/
96//
97//  Determines the intersection point of the line segment defined by points A and B
98//  with the line segment defined by points C and D.
99//
100//  Returns YES if the intersection point was found, and stores that point in X,Y.
101//  Returns NO if there is no determinable intersection point, in which case X,Y will
102//  be unmodified.
103
104int __line_segment_intersection(
105        double Ax, double Ay,
106        double Bx, double By,
107        double Cx, double Cy,
108        double Dx, double Dy,
109        double *X, double *Y) {
110
111    double distAB, theCos, theSin, newX, ABpos;
112
113    //  Fail if either line segment is zero-length.
114    if ( (Ax == Bx && Ay == By) || (Cx == Dx && Cy == Dy) ) return NO ;
115
116    //  Fail if the segments share an end-point.
117    if ( (Ax == Cx && Ay == Cy) || (Bx == Cx && By == Cy)
118            || (Ax == Dx && Ay == Dy) || (Bx == Dx && By == Dy) ) {
119        return NO;
120    }
121
122    //  (1) Translate the system so that point A is on the origin.
123    Bx -= Ax;
124    By -= Ay;
125    Cx -= Ax;
126    Cy -= Ay;
127    Dx -= Ax;
128    Dy -= Ay;
129
130    //  Discover the length of segment A-B.
131    distAB = sqrt(Bx * Bx + By * By);
132
133    //  (2) Rotate the system so that point B is on the positive X axis.
134    theCos = Bx / distAB;
135    theSin = By / distAB;
136    newX = Cx * theCos + Cy*theSin;
137    Cy = Cy * theCos - Cx*theSin;
138    Cx = newX;
139    newX = Dx * theCos + Dy*theSin;
140    Dy = Dy * theCos - Dx*theSin;
141    Dx = newX;
142
143    //  Fail if segment C-D doesn't cross line A-B.
144    if ( (Cy < 0. && Dy < 0.) || (Cy >= 0. && Dy >= 0.) ) return NO;
145
146    //  (3) Discover the position of the intersection point along line A-B.
147    ABpos = Dx + (Cx - Dx) * Dy / (Dy - Cy);
148
149    //  Fail if segment C-D crosses line A-B outside of segment A-B.
150    if (ABpos < 0. || ABpos > distAB) return NO;
151
152    //  (4) Apply the discovered position to line A-B in the original coordinate system.
153    *X = Ax + ABpos*theCos;
154    *Y = Ay + ABpos*theSin;
155
156    //  Success.
157    return YES;
158}
159
160/*
161WORK IN PROGRESS TO OPTIMISE INTERSECTION
162int __intersection(double x0, double y0,
163                   double x1, double y1) {
164
165
166    x0 = line0[0,0]; y0 = line0[0,1]
167    x1 = line0[1,0]; y1 = line0[1,1]
168
169    x2 = line1[0,0]; y2 = line1[0,1]
170    x3 = line1[1,0]; y3 = line1[1,1]
171
172    denom = (y3-y2)*(x1-x0) - (x3-x2)*(y1-y0)
173    u0 = (x3-x2)*(y0-y2) - (y3-y2)*(x0-x2)
174    u1 = (x2-x0)*(y1-y0) - (y2-y0)*(x1-x0)
175       
176    if allclose(denom, 0.0):
177        # Lines are parallel - check if they coincide on a shared a segment
178
179        if allclose( [u0, u1], 0.0 ):
180            # We now know that the lines if continued coincide
181            # The remaining check will establish if the finite lines share a segment
182
183            line0_starts_on_line1 = line0_ends_on_line1 =\
184            line1_starts_on_line0 = line1_ends_on_line0 = False
185               
186            if point_on_line([x0, y0], line1):
187                line0_starts_on_line1 = True
188
189            if point_on_line([x1, y1], line1):
190                line0_ends_on_line1 = True
191 
192            if point_on_line([x2, y2], line0):
193                line1_starts_on_line0 = True
194
195            if point_on_line([x3, y3], line0):
196                line1_ends_on_line0 = True                               
197
198            if not(line0_starts_on_line1 or line0_ends_on_line1\
199               or line1_starts_on_line0 or line1_ends_on_line0):
200                # Lines are parallel and would coincide if extended, but not as they are.
201                return 3, None
202
203
204            # One line fully included in the other. Use direction of included line
205            if line0_starts_on_line1 and line0_ends_on_line1:
206                # Shared segment is line0 fully included in line1
207                segment = array([[x0, y0], [x1, y1]])               
208
209            if line1_starts_on_line0 and line1_ends_on_line0:
210                # Shared segment is line1 fully included in line0
211                segment = array([[x2, y2], [x3, y3]])
212           
213
214            # Overlap with lines are oriented the same way
215            if line0_starts_on_line1 and line1_ends_on_line0:
216                # Shared segment from line0 start to line 1 end
217                segment = array([[x0, y0], [x3, y3]])
218
219            if line1_starts_on_line0 and line0_ends_on_line1:
220                # Shared segment from line1 start to line 0 end
221                segment = array([[x2, y2], [x1, y1]])                               
222
223
224            # Overlap in opposite directions - use direction of line0
225            if line0_starts_on_line1 and line1_starts_on_line0:
226                # Shared segment from line0 start to line 1 end
227                segment = array([[x0, y0], [x2, y2]])
228
229            if line0_ends_on_line1 and line1_ends_on_line0:
230                # Shared segment from line0 start to line 1 end
231                segment = array([[x3, y3], [x1, y1]])               
232
233               
234            return 2, segment
235        else:
236            # Lines are parallel but they do not coincide
237            return 4, None #FIXME (Ole): Add distance here instead of None
238           
239    else:
240        # Lines are not parallel or coinciding
241        u0 = u0/denom
242        u1 = u1/denom       
243
244        x = x0 + u0*(x1-x0)
245        y = y0 + u0*(y1-y0)
246
247        # Sanity check - can be removed to speed up if needed
248        assert allclose(x, x2 + u1*(x3-x2))
249        assert allclose(y, y2 + u1*(y3-y2))       
250
251        # Check if point found lies within given line segments
252        if 0.0 <= u0 <= 1.0 and 0.0 <= u1 <= 1.0:
253            # We have intersection
254
255            return 1, array([x, y])
256        else:
257            # No intersection
258            return 0, None
259
260
261}
262*/
263
264
265
266int __interpolate_polyline(int number_of_nodes,
267                           int number_of_points,
268                           double* data,
269                           double* polyline_nodes,
270                           long* gauge_neighbour_id,
271                           double* interpolation_points,                               
272                           double* interpolated_values,
273                           double rtol,
274                           double atol) {
275                           
276  int j, i, neighbour_id;
277  double x0, y0, x1, y1, x, y;
278  double segment_len, segment_delta, slope, alpha;
279
280  for (j=0; j<number_of_nodes; j++) { 
281
282    neighbour_id = gauge_neighbour_id[j];
283       
284    // FIXME(Ole): I am convinced that gauge_neighbour_id can be discarded, but need to check with John J.
285    // Keep it for now (17 Jan 2009)
286    // When gone, we can simply interpolate between neighbouring nodes, i.e. neighbour_id = j+1.
287    // and the test below becomes something like: if j < number_of_nodes... 
288       
289    if (neighbour_id >= 0) {
290      x0 = polyline_nodes[2*j];
291      y0 = polyline_nodes[2*j+1];
292     
293      x1 = polyline_nodes[2*neighbour_id];
294      y1 = polyline_nodes[2*neighbour_id+1];     
295     
296           
297      segment_len = dist(x1-x0, y1-y0);
298      segment_delta = data[neighbour_id] - data[j];           
299      slope = segment_delta/segment_len;
300           
301      for (i=0; i<number_of_points; i++) {               
302        x = interpolation_points[2*i];
303        y = interpolation_points[2*i+1];       
304       
305        if (__point_on_line(x, y, x0, y0, x1, y1, rtol, atol)) {
306          alpha = dist(x-x0, y-y0);
307          interpolated_values[i] = slope*alpha + data[j];
308        }
309      }
310    }
311  }
312                           
313  return 0;                         
314}                                                     
315
316
317int __triangle_polygon_overlap(double* polygon,
318                               double* triangle,
319                               int polygon_number_of_vertices)
320{
321    int i, ii, j, jj, A, B;
322    double p0_x, p0_y, p1_x, p1_y, pp_x, pp_y;
323    double t0_x, t0_y, t1_x, t1_y, tp_x, tp_y;
324    double u_x, u_y, v_x, v_y, w_x, w_y;
325    double u_dot_tp, v_dot_tp, v_dot_pp, w_dot_pp;
326    double a, b;
327   
328    p0_x = polygon[0];
329    p0_y = polygon[1];
330   
331    A = 0;
332    B = 0;
333   
334    for (i = 1; i < polygon_number_of_vertices + 1; i++)
335    {
336        ii = i%polygon_number_of_vertices;
337       
338        p1_x = polygon[2*ii];
339        p1_y = polygon[2*ii + 1];
340       
341        pp_x = -(p1_y - p0_y);
342        pp_y = p1_x - p0_x;
343 
344        t0_x = triangle[0];
345        t0_y = triangle[1];
346 
347        for (j = 1; j < 4; j++)
348        {
349            jj = j%3;
350                     
351            t1_x = triangle[2*jj];
352            t1_y = triangle[2*jj + 1];
353           
354            tp_x = -(t1_y - t0_y); //perpendicular to triangle vector
355            tp_y = t1_x - t0_x; //perpendicular to polygon vector
356       
357            u_x = p1_x - p0_x;
358            u_y = p1_y - p0_y;
359            v_x = t0_x - p0_x;
360            v_y = t0_y - p0_y;
361            w_x = t1_x - t0_x;
362            w_y = t1_y - t0_y;
363
364            u_dot_tp = (u_x*tp_x) + (u_y*tp_y);
365           
366            if (u_dot_tp != 0.0f) //Vectors are not parallel
367            {
368                v_dot_tp = (v_x*tp_x) + (v_y*tp_y);
369                v_dot_pp = (v_x*pp_x) + (v_y*pp_y);
370                w_dot_pp = (w_x*pp_x) + (w_y*pp_y);
371               
372                a = v_dot_tp/u_dot_tp;
373                b = -v_dot_pp/w_dot_pp;
374                             
375                if (a >= 0.0f && a <= 1.0f && b >=0.0f && b <=1.0f)
376                {
377                    return 1; //overlap
378                }
379               
380                if (b >= 0.0f && b <= 1.0f && a > 1.0f)
381                {
382                    A++; 
383                }
384               
385                if (a >= 0.0f && a <= 1.0f && b > 1.0f)
386                {
387                    B++; 
388                }
389               
390                if (A == 4 || B == 3)
391                {
392                    return 1; //overlap
393                }
394            }
395           
396            t0_x = t1_x;
397            t0_y = t1_y;
398        }
399       
400        p0_x = p1_x;
401        p0_y = p1_y;
402    }
403
404    return 0; //no overlap
405}
406                 
407
408int __polygon_overlap(double* polygon,
409                      double* triangles,
410                      long* indices,
411                      int M, //number of triangles
412                      int polygon_number_of_vertices)
413{
414    double* triangle;
415    int i, inside_index, outside_index;
416   
417    inside_index = 0;    // Keep track of triangles that overlap
418    outside_index = M - 1; // Keep track of triangles that don't overlap (starting from end)
419   
420    for (i = 0; i < M; i++)
421    {
422        triangle = triangles + 6*i;
423       
424        if (__triangle_polygon_overlap(polygon, 
425                                      triangle, 
426                                      polygon_number_of_vertices))
427        {
428            indices[inside_index] = i;
429            inside_index++;
430        }
431        else
432        {
433            indices[outside_index] = i;
434            outside_index -= 1;           
435        }
436    }
437   
438    return inside_index;
439}             
440
441
442int __triangle_line_intersect(double* line,
443                              double* triangle)
444{
445    int j, jj, A, B;
446    double p0_x, p0_y, p1_x, p1_y, pp_x, pp_y;
447    double t0_x, t0_y, t1_x, t1_y, tp_x, tp_y;
448    double u_x, u_y, v_x, v_y, w_x, w_y;
449    double u_dot_tp, v_dot_tp, v_dot_pp, w_dot_pp;
450    double a, b;
451   
452    p0_x = line[0];
453    p0_y = line[1];
454    p1_x = line[2];
455    p1_y = line[3];
456   
457    pp_x = -(p1_y - p0_y);
458    pp_y = p1_x - p0_x;
459   
460    A = 0;
461    B = 0;
462   
463    t0_x = triangle[0];
464    t0_y = triangle[1];
465
466    for (j = 1; j < 4; j++)
467    {
468        jj = j%3;
469                 
470        t1_x = triangle[2*jj];
471        t1_y = triangle[2*jj + 1];
472       
473        tp_x = -(t1_y - t0_y); //perpendicular to triangle vector
474        tp_y = t1_x - t0_x; 
475   
476        u_x = p1_x - p0_x;
477        u_y = p1_y - p0_y;
478        v_x = t0_x - p0_x;
479        v_y = t0_y - p0_y;
480        w_x = t1_x - t0_x;
481        w_y = t1_y - t0_y;
482
483        u_dot_tp = (u_x*tp_x) + (u_y*tp_y);
484       
485        if (u_dot_tp != 0.0f) //If vectors are not parallel, continue
486        {
487            v_dot_tp = (v_x*tp_x) + (v_y*tp_y);
488            v_dot_pp = (v_x*pp_x) + (v_y*pp_y);
489            w_dot_pp = (w_x*pp_x) + (w_y*pp_y);
490           
491            a = v_dot_tp/u_dot_tp;
492            b = -v_dot_pp/w_dot_pp;
493                         
494            if (a >= 0.0f && a <= 1.0f && b >=0.0f && b <=1.0f)
495            {
496                return 1; //intersect
497            }
498           
499            if (a > 1.0f && b >= 0.0f && b <= 1.0f)
500            {
501                A++; 
502            }
503           
504            if (a < 0.0f && b >= 0.0f && b <= 1.0f)
505            {
506                B++; 
507            }
508        }
509       
510        t0_x = t1_x;
511        t0_y = t1_y;
512    }
513   
514    if (A >= 1 && B >= 1)
515    {
516        return 1; //polygon sits completely inside a triangle
517    }
518   
519    return 0; //no intersection
520}
521                 
522
523int __line_intersect(double* line,
524                     double* triangles,
525                     long* indices,
526                     int M) //number of triangles
527{
528    double* triangle;
529    int i, inside_index, outside_index;
530   
531    inside_index = 0;    // Keep track of triangles that intersect
532    outside_index = M - 1; // Keep track of triangles that don't intersect (starting from end)
533   
534    for (i = 0; i < M; i++)
535    {
536        triangle = triangles + 6*i;
537       
538        if (__triangle_line_intersect(line, 
539                                      triangle))
540        {
541            indices[inside_index] = i;
542            inside_index++;
543        }
544        else
545        {
546            indices[outside_index] = i;
547            outside_index -= 1;           
548        }
549    }
550   
551    return inside_index;
552}             
553
554
555
556int __is_inside_triangle(double* point,
557                         double* triangle,
558                         int closed,
559                         double rtol,
560                         double atol) {
561                         
562  double vx, vy, v0x, v0y, v1x, v1y;
563  double a00, a10, a01, a11, b0, b1;
564  double denom, alpha, beta;
565 
566  double x, y; // Point coordinates
567  int i, j, res;
568
569  x = point[0];
570  y = point[1];
571 
572  // Quickly reject points that are clearly outside
573  if ((x < triangle[0]) && 
574      (x < triangle[2]) && 
575      (x < triangle[4])) return 0;       
576     
577  if ((x > triangle[0]) && 
578      (x > triangle[2]) && 
579      (x > triangle[4])) return 0;             
580 
581  if ((y < triangle[1]) && 
582      (y < triangle[3]) && 
583      (y < triangle[5])) return 0;       
584     
585  if ((y > triangle[1]) && 
586      (y > triangle[3]) && 
587      (y > triangle[5])) return 0;             
588 
589 
590  // v0 = C-A
591  v0x = triangle[4]-triangle[0]; 
592  v0y = triangle[5]-triangle[1];
593 
594  // v1 = B-A   
595  v1x = triangle[2]-triangle[0]; 
596  v1y = triangle[3]-triangle[1];
597
598  // First check if point lies wholly inside triangle
599  a00 = v0x*v0x + v0y*v0y; // innerproduct(v0, v0)
600  a01 = v0x*v1x + v0y*v1y; // innerproduct(v0, v1)
601  a10 = a01;               // innerproduct(v1, v0)
602  a11 = v1x*v1x + v1y*v1y; // innerproduct(v1, v1)
603   
604  denom = a11*a00 - a01*a10;
605
606  if (fabs(denom) > 0.0) {
607    // v = point-A 
608    vx = x - triangle[0]; 
609    vy = y - triangle[1];     
610   
611    b0 = v0x*vx + v0y*vy; // innerproduct(v0, v)       
612    b1 = v1x*vx + v1y*vy; // innerproduct(v1, v)           
613   
614    alpha = (b0*a11 - b1*a01)/denom;
615    beta = (b1*a00 - b0*a10)/denom;       
616   
617    if ((alpha > 0.0) && (beta > 0.0) && (alpha+beta < 1.0)) return 1;
618  }
619
620  if (closed) {
621    // Check if point lies on one of the edges
622       
623    for (i=0; i<3; i++) {
624      j = (i+1) % 3; // Circular index into triangle array
625      res = __point_on_line(x, y,
626                            triangle[2*i], triangle[2*i+1], 
627                            triangle[2*j], triangle[2*j+1],                         
628                            rtol, atol);
629      if (res) return 1;
630    }
631  }
632               
633  // Default return if point is outside triangle                         
634  return 0;                                             
635}                                                                             
636
637
638int __separate_points_by_polygon(int M,     // Number of points
639                                 int N,     // Number of polygon vertices
640                                 double* points,
641                                 double* polygon,
642                                 long* indices,  // M-Array for storage indices
643                                 int closed,
644                                 int verbose) {
645
646  double minpx, maxpx, minpy, maxpy, x, y, px_i, py_i, px_j, py_j, rtol=0.0, atol=0.0;
647  int i, j, k, outside_index, inside_index, inside;
648
649  // Find min and max of poly used for optimisation when points
650  // are far away from polygon
651 
652  // FIXME(Ole): Pass in rtol and atol from Python
653
654  minpx = polygon[0]; maxpx = minpx;
655  minpy = polygon[1]; maxpy = minpy;
656
657  for (i=0; i<N; i++) {
658    px_i = polygon[2*i];
659    py_i = polygon[2*i + 1];
660
661    if (px_i < minpx) minpx = px_i;
662    if (px_i > maxpx) maxpx = px_i;
663    if (py_i < minpy) minpy = py_i;
664    if (py_i > maxpy) maxpy = py_i;
665  }
666
667  // Begin main loop (for each point)
668  inside_index = 0;    // Keep track of points inside
669  outside_index = M-1; // Keep track of points outside (starting from end)   
670  if (verbose){
671     printf("Separating %d points\n", M);
672  } 
673  for (k=0; k<M; k++) {
674    if (verbose){
675      if (k %((M+10)/10)==0) printf("Doing %d of %d\n", k, M);
676    }
677   
678    x = points[2*k];
679    y = points[2*k + 1];
680
681    inside = 0;
682
683    // Optimisation
684    if ((x > maxpx) || (x < minpx) || (y > maxpy) || (y < minpy)) {
685      // Nothing
686    } else {   
687      // Check polygon
688      for (i=0; i<N; i++) {
689        j = (i+1)%N;
690
691        px_i = polygon[2*i];
692        py_i = polygon[2*i+1];
693        px_j = polygon[2*j];
694        py_j = polygon[2*j+1];
695
696        // Check for case where point is contained in line segment
697        if (__point_on_line(x, y, px_i, py_i, px_j, py_j, rtol, atol)) {
698          if (closed == 1) {
699            inside = 1;
700          } else {
701            inside = 0;
702          }
703          break;
704        } else {
705          //Check if truly inside polygon
706          if ( ((py_i < y) && (py_j >= y)) ||
707               ((py_j < y) && (py_i >= y)) ) {
708            if (px_i + (y-py_i)/(py_j-py_i)*(px_j-px_i) < x)
709              inside = 1-inside;
710          }
711        }
712      }
713    } 
714    if (inside == 1) {
715      indices[inside_index] = k;
716      inside_index += 1;
717    } else {
718      indices[outside_index] = k;
719      outside_index -= 1;   
720    }
721  } // End k
722
723  return inside_index;
724}
725
726
727
728// Gateways to Python
729PyObject *_point_on_line(PyObject *self, PyObject *args) {
730  //
731  // point_on_line(x, y, x0, y0, x1, y1)
732  //
733
734  double x, y, x0, y0, x1, y1, rtol, atol;
735  int res;
736  PyObject *result;
737
738  // Convert Python arguments to C
739  if (!PyArg_ParseTuple(args, "dddddddd", &x, &y, &x0, &y0, &x1, &y1, &rtol, &atol)) {
740    PyErr_SetString(PyExc_RuntimeError, 
741                    "point_on_line could not parse input");   
742    return NULL;
743  }
744
745
746  // Call underlying routine
747  res = __point_on_line(x, y, x0, y0, x1, y1, rtol, atol);
748
749  // Return values a and b
750  result = Py_BuildValue("i", res);
751  return result;
752}
753
754
755
756// Gateways to Python
757PyObject *_interpolate_polyline(PyObject *self, PyObject *args) {
758  //
759  // _interpolate_polyline(data, polyline_nodes, gauge_neighbour_id, interpolation_points
760  //                       interpolated_values):
761  //
762
763 
764  PyArrayObject
765    *data,
766    *polyline_nodes,
767    *gauge_neighbour_id,
768    *interpolation_points,
769    *interpolated_values;
770
771  double rtol, atol; 
772  int number_of_nodes, number_of_points, res;
773 
774  // Convert Python arguments to C
775  if (!PyArg_ParseTuple(args, "OOOOOdd",
776                        &data,
777                        &polyline_nodes,
778                        &gauge_neighbour_id,
779                        &interpolation_points,
780                        &interpolated_values,
781                        &rtol,
782                        &atol)) {
783   
784    PyErr_SetString(PyExc_RuntimeError, 
785                    "_interpolate_polyline could not parse input");
786    return NULL;
787  }
788
789  // check that numpy array objects arrays are C contiguous memory
790  CHECK_C_CONTIG(data);
791  CHECK_C_CONTIG(polyline_nodes);
792  CHECK_C_CONTIG(gauge_neighbour_id);
793  CHECK_C_CONTIG(interpolation_points);
794  CHECK_C_CONTIG(interpolated_values);
795
796  number_of_nodes = polyline_nodes -> dimensions[0];  // Number of nodes in polyline
797  number_of_points = interpolation_points -> dimensions[0];   //Number of points
798 
799
800  // Call underlying routine
801  res = __interpolate_polyline(number_of_nodes,
802                               number_of_points,
803                               (double*) data -> data,
804                               (double*) polyline_nodes -> data,
805                               (long*) gauge_neighbour_id -> data,
806                               (double*) interpolation_points -> data,                         
807                               (double*) interpolated_values -> data,
808                               rtol,
809                               atol);                                                 
810
811  // Return
812  return Py_BuildValue(""); 
813}
814
815
816PyObject *_polygon_overlap(PyObject *self, PyObject *args) 
817{
818  //
819  // _polygon_triangle_overlap(polygon, triangle)
820  //
821
822 
823  PyArrayObject
824    *polygon,
825    *triangles,
826    *indices;
827   
828    int res;
829
830  PyObject *result;
831     
832  // Convert Python arguments to C
833  if (!PyArg_ParseTuple(args, "OOO",
834                        &polygon,
835                        &triangles,
836            &indices)) {
837   
838    PyErr_SetString(PyExc_RuntimeError, 
839                    "_polygon_triangle_overlap could not parse input");
840    return NULL;
841  }
842
843  // Call underlying routine
844  res = __polygon_overlap((double*) polygon->data,
845                             (double*) triangles->data,
846                 (long*) indices->data,
847                 (int) triangles->dimensions[0]/3,
848                 (int) polygon->dimensions[0]);                                               
849
850
851  // Return result
852  result = Py_BuildValue("i", res);
853  return result; 
854}
855
856PyObject *_line_intersect(PyObject *self, PyObject *args) 
857{
858  //
859  // _polygon_triangle_overlap(polygon, triangle)
860  //
861 
862    //def separate_points_by_polygon(points, polygon, closed, verbose, one_point):
863  //  """Determine whether points are inside or outside a polygon
864  //
865  //  Input:
866  //     point - Tuple of (x, y) coordinates, or list of tuples
867  //     polygon - list of vertices of polygon
868  //     closed - (optional) determine whether points on boundary should be
869  //     regarded as belonging to the polygon (closed = True)
870  //     or not (closed = False)
871
872  //
873  //  Output:
874  //     indices: array of same length as points with indices of points falling
875  //     inside the polygon listed from the beginning and indices of points
876  //     falling outside listed from the end.
877  //     
878  //     count: count of points falling inside the polygon
879  //     
880  //     The indices of points inside are obtained as indices[:count]
881  //     The indices of points outside are obtained as indices[count:]       
882  //
883  //  Examples:
884  //     separate_polygon( [[0.5, 0.5], [1, -0.5], [0.3, 0.2]] )
885  //     will return the indices [0, 2, 1] as only the first and the last point
886  //     is inside the unit square
887  //
888  //  Remarks:
889  //     The vertices may be listed clockwise or counterclockwise and
890  //     the first point may optionally be repeated.
891  //     Polygons do not need to be convex.
892  //     Polygons can have holes in them and points inside a hole is
893  //     regarded as being outside the polygon.
894  //
895  //
896  //  Algorithm is based on work by Darel Finley,
897  //  http://www.alienryderflex.com/polygon/
898  //
899 
900 
901  PyArrayObject
902    *line,
903    *triangles,
904    *indices;
905   
906    int res;
907
908  PyObject *result;
909     
910  // Convert Python arguments to C
911  if (!PyArg_ParseTuple(args, "OOO",
912                        &line,
913                        &triangles,
914            &indices)) {
915   
916    PyErr_SetString(PyExc_RuntimeError, 
917                    "_polygon_triangle_overlap could not parse input");
918    return NULL;
919  }
920
921  // Call underlying routine
922  res = __line_intersect((double*) line->data,
923                                     (double*) triangles->data,
924                         (long*) indices->data,
925                         (int) triangles->dimensions[0]/3);                                                   
926
927
928  // Return result
929  result = Py_BuildValue("i", res);
930  return result; 
931}
932     
933PyObject *_is_inside_triangle(PyObject *self, PyObject *args) {
934  //
935  // _is_inside_triangle(point, triangle, int(closed), rtol, atol)
936  //
937
938 
939  PyArrayObject
940    *point,
941    *triangle;
942
943  double rtol, atol; 
944  int closed, res;
945
946  PyObject *result;
947     
948  // Convert Python arguments to C
949  if (!PyArg_ParseTuple(args, "OOidd",
950                        &point,
951                        &triangle,
952                        &closed,
953                        &rtol,
954                        &atol)) {
955   
956    PyErr_SetString(PyExc_RuntimeError, 
957                    "_is_inside_triangle could not parse input");
958    return NULL;
959  }
960
961  // Call underlying routine
962  res = __is_inside_triangle((double*) point -> data,
963                             (double*) triangle -> data,
964                             closed,
965                             rtol,
966                             atol);                                                   
967
968
969  // Return result
970  result = Py_BuildValue("i", res);
971  return result; 
972}
973     
974
975
976/*
977PyObject *_intersection(PyObject *self, PyObject *args) {
978  //
979  // intersection(x0, y0, x1, y1)
980  //
981
982  double x, y, x0, y0, x1, y1;
983  int res;
984  PyObject *result;
985
986  // Convert Python arguments to C
987  if (!PyArg_ParseTuple(args, "dddddd", &x, &y, &x0, &y0, &x1, &y1)) {
988    PyErr_SetString(PyExc_RuntimeError,
989                    "point_on_line could not parse input");   
990    return NULL;
991  }
992
993
994  // Call underlying routine
995  res = __intersection(x, y, x0, y0, x1, y1);
996
997  // Return values a and b
998  result = Py_BuildValue("i", res);
999  return result;
1000}
1001*/
1002
1003
1004PyObject *_separate_points_by_polygon(PyObject *self, PyObject *args) {
1005  //def separate_points_by_polygon(points, polygon, closed, verbose, one_point):
1006  //  """Determine whether points are inside or outside a polygon
1007  //
1008  //  Input:
1009  //     point - Tuple of (x, y) coordinates, or list of tuples
1010  //     polygon - list of vertices of polygon
1011  //     closed - (optional) determine whether points on boundary should be
1012  //     regarded as belonging to the polygon (closed = True)
1013  //     or not (closed = False)
1014
1015  //
1016  //  Output:
1017  //     indices: array of same length as points with indices of points falling
1018  //     inside the polygon listed from the beginning and indices of points
1019  //     falling outside listed from the end.
1020  //     
1021  //     count: count of points falling inside the polygon
1022  //     
1023  //     The indices of points inside are obtained as indices[:count]
1024  //     The indices of points outside are obtained as indices[count:]       
1025  //
1026  //  Examples:
1027  //     separate_polygon( [[0.5, 0.5], [1, -0.5], [0.3, 0.2]] )
1028  //     will return the indices [0, 2, 1] as only the first and the last point
1029  //     is inside the unit square
1030  //
1031  //  Remarks:
1032  //     The vertices may be listed clockwise or counterclockwise and
1033  //     the first point may optionally be repeated.
1034  //     Polygons do not need to be convex.
1035  //     Polygons can have holes in them and points inside a hole is
1036  //     regarded as being outside the polygon.
1037  //
1038  //
1039  //  Algorithm is based on work by Darel Finley,
1040  //  http://www.alienryderflex.com/polygon/
1041  //
1042  //
1043
1044  PyArrayObject
1045    *points,
1046    *polygon,
1047    *indices;
1048
1049  int closed, verbose; //Flags
1050  int count, M, N;
1051
1052  // Convert Python arguments to C
1053  if (!PyArg_ParseTuple(args, "OOOii",
1054                        &points,
1055                        &polygon,
1056                        &indices,
1057                        &closed,
1058                        &verbose)) {
1059   
1060
1061    PyErr_SetString(PyExc_RuntimeError, 
1062                    "separate_points_by_polygon could not parse input");
1063    return NULL;
1064  }
1065
1066  // check that points, polygon and indices arrays are C contiguous
1067  CHECK_C_CONTIG(points);
1068  CHECK_C_CONTIG(polygon);
1069  CHECK_C_CONTIG(indices);
1070
1071  M = points -> dimensions[0];   //Number of points
1072  N = polygon -> dimensions[0];  //Number of vertices in polygon
1073
1074  //FIXME (Ole): This isn't send to Python's sys.stdout
1075  if (verbose) printf("Got %d points and %d polygon vertices\n", M, N);
1076 
1077  //Call underlying routine
1078  count = __separate_points_by_polygon(M, N,
1079                                       (double*) points -> data,
1080                                       (double*) polygon -> data,
1081                                       (long*) indices -> data,
1082                                       closed, verbose);
1083 
1084  //NOTE: return number of points inside..
1085  return Py_BuildValue("i", count);
1086}
1087
1088
1089
1090// Method table for python module
1091static struct PyMethodDef MethodTable[] = {
1092  /* The cast of the function is necessary since PyCFunction values
1093   * only take two PyObject* parameters, and rotate() takes
1094   * three.
1095   */
1096
1097  {"_point_on_line", _point_on_line, METH_VARARGS, "Print out"},
1098  //{"_intersection", _intersection, METH_VARARGS, "Print out"}, 
1099  {"_separate_points_by_polygon", _separate_points_by_polygon, 
1100                                 METH_VARARGS, "Print out"},
1101  {"_interpolate_polyline", _interpolate_polyline, 
1102                                 METH_VARARGS, "Print out"},                             
1103  {"_polygon_overlap", _polygon_overlap, 
1104                                 METH_VARARGS, "Print out"},
1105  {"_line_intersect", _line_intersect, 
1106                                 METH_VARARGS, "Print out"},                               
1107  {"_is_inside_triangle", _is_inside_triangle, 
1108                                 METH_VARARGS, "Print out"},
1109  {NULL, NULL, 0, NULL}   /* sentinel */
1110};
1111
1112
1113// Module initialisation
1114void initpolygon_ext(void){
1115  Py_InitModule("polygon_ext", MethodTable);
1116
1117  import_array();     //Necessary for handling of NumPY structures
1118}
1119
1120
1121
Note: See TracBrowser for help on using the repository browser.