PostGIS 3.7.0dev-r@@SVN_REVISION@@
Loading...
Searching...
No Matches
measures.c
Go to the documentation of this file.
1/**********************************************************************
2 *
3 * PostGIS - Spatial Types for PostgreSQL
4 * http://postgis.net
5 *
6 * PostGIS is free software: you can redistribute it and/or modify
7 * it under the terms of the GNU General Public License as published by
8 * the Free Software Foundation, either version 2 of the License, or
9 * (at your option) any later version.
10 *
11 * PostGIS is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 * GNU General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public License
17 * along with PostGIS. If not, see <http://www.gnu.org/licenses/>.
18 *
19 **********************************************************************
20 *
21 * Copyright 2001-2006 Refractions Research Inc.
22 * Copyright 2010 Nicklas Avén
23 * Copyright 2012 Paul Ramsey
24 * Copyright 2019 Darafei Praliaskouski
25 *
26 **********************************************************************/
27
28#include <string.h>
29#include <stdlib.h>
30
31#include "../postgis_config.h"
32//#define POSTGIS_DEBUG_LEVEL 5
33
34#include "measures.h"
35#include "lwgeom_log.h"
36
37/*------------------------------------------------------------------------------------------------------------
38Initializing functions
39The functions starting the distance-calculation processes
40--------------------------------------------------------------------------------------------------------------*/
41
42LWGEOM *
43lwgeom_closest_line(const LWGEOM *lw1, const LWGEOM *lw2)
44{
45 return lw_dist2d_distanceline(lw1, lw2, lw1->srid, DIST_MIN);
46}
47
48LWGEOM *
49lwgeom_furthest_line(const LWGEOM *lw1, const LWGEOM *lw2)
50{
51 return lw_dist2d_distanceline(lw1, lw2, lw1->srid, DIST_MAX);
52}
53
54LWGEOM *
55lwgeom_closest_point(const LWGEOM *lw1, const LWGEOM *lw2)
56{
57 return lw_dist2d_distancepoint(lw1, lw2, lw1->srid, DIST_MIN);
58}
59
60LWGEOM *
61lwgeom_furthest_point(const LWGEOM *lw1, const LWGEOM *lw2)
62{
63 return lw_dist2d_distancepoint(lw1, lw2, lw1->srid, DIST_MAX);
64}
65
66void
68{
69 dl->twisted = -1;
70 dl->p1.x = dl->p1.y = 0.0;
71 dl->p2.x = dl->p2.y = 0.0;
72 dl->mode = mode;
73 dl->tolerance = 0.0;
74 if (mode == DIST_MIN)
75 dl->distance = FLT_MAX;
76 else
77 dl->distance = -1 * FLT_MAX;
78}
79
80static void
81lw_dist2d_distpts_set(DISTPTS *dl, double distance, const POINT2D *p1, const POINT2D *p2)
82{
83 int update = (dl->mode == DIST_MIN) ? (distance < dl->distance) : (distance > dl->distance);
84 if (update)
85 {
86 dl->distance = distance;
87 dl->p1 = *p1;
88 dl->p2 = *p2;
89 }
90}
91
95LWGEOM *
96lw_dist2d_distanceline(const LWGEOM *lw1, const LWGEOM *lw2, int32_t srid, int mode)
97{
98 double x1, x2, y1, y2;
99
100 double initdistance = (mode == DIST_MIN ? FLT_MAX : -1.0);
101 DISTPTS thedl;
102 LWPOINT *lwpoints[2];
103 LWGEOM *result;
104
105 thedl.mode = mode;
106 thedl.distance = initdistance;
107 thedl.tolerance = 0.0;
108
109 LWDEBUG(2, "lw_dist2d_distanceline is called");
110
111 if (!lw_dist2d_comp(lw1, lw2, &thedl))
112 {
113 /*should never get here. all cases ought to be error handled earlier*/
114 lwerror("Some unspecified error.");
116 }
117
118 /*if thedl.distance is unchanged there where only empty geometries input*/
119 if (thedl.distance == initdistance)
120 {
121 LWDEBUG(3, "didn't find geometries to measure between, returning null");
123 }
124 else
125 {
126 x1 = thedl.p1.x;
127 y1 = thedl.p1.y;
128 x2 = thedl.p2.x;
129 y2 = thedl.p2.y;
130
131 lwpoints[0] = lwpoint_make2d(srid, x1, y1);
132 lwpoints[1] = lwpoint_make2d(srid, x2, y2);
133
134 result = (LWGEOM *)lwline_from_ptarray(srid, 2, lwpoints);
135 }
136 return result;
137}
138
142LWGEOM *
143lw_dist2d_distancepoint(const LWGEOM *lw1, const LWGEOM *lw2, int32_t srid, int mode)
144{
145 double x, y;
146 DISTPTS thedl;
147 double initdistance = FLT_MAX;
148 LWGEOM *result;
149
150 thedl.mode = mode;
151 thedl.distance = initdistance;
152 thedl.tolerance = 0;
153
154 LWDEBUG(2, "lw_dist2d_distancepoint is called");
155
156 if (!lw_dist2d_comp(lw1, lw2, &thedl))
157 {
158 /*should never get here. all cases ought to be error handled earlier*/
159 lwerror("Some unspecified error.");
161 }
162 if (thedl.distance == initdistance)
163 {
164 LWDEBUG(3, "didn't find geometries to measure between, returning null");
166 }
167 else
168 {
169 x = thedl.p1.x;
170 y = thedl.p1.y;
171 result = (LWGEOM *)lwpoint_make2d(srid, x, y);
172 }
173 return result;
174}
175
179double
180lwgeom_maxdistance2d(const LWGEOM *lw1, const LWGEOM *lw2)
181{
182 LWDEBUG(2, "lwgeom_maxdistance2d is called");
183
184 return lwgeom_maxdistance2d_tolerance(lw1, lw2, 0.0);
185}
186
191double
192lwgeom_maxdistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
193{
194 /*double thedist;*/
195 DISTPTS thedl;
196 LWDEBUG(2, "lwgeom_maxdistance2d_tolerance is called");
197 thedl.mode = DIST_MAX;
198 thedl.distance = -1;
199 thedl.tolerance = tolerance;
200 if (lw_dist2d_comp(lw1, lw2, &thedl))
201 return thedl.distance;
202
203 /*should never get here. all cases ought to be error handled earlier*/
204 lwerror("Some unspecified error.");
205 return -1;
206}
207
211double
212lwgeom_mindistance2d(const LWGEOM *lw1, const LWGEOM *lw2)
213{
214 return lwgeom_mindistance2d_tolerance(lw1, lw2, 0.0);
215}
216
221double
222lwgeom_mindistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
223{
224 DISTPTS thedl;
225 LWDEBUG(2, "lwgeom_mindistance2d_tolerance is called");
226 thedl.mode = DIST_MIN;
227 thedl.distance = FLT_MAX;
228 thedl.tolerance = tolerance;
229 if (lw_dist2d_comp(lw1, lw2, &thedl))
230 return thedl.distance;
231 /*should never get here. all cases ought to be error handled earlier*/
232 lwerror("Some unspecified error.");
233 return FLT_MAX;
234}
235
236/*------------------------------------------------------------------------------------------------------------
237End of Initializing functions
238--------------------------------------------------------------------------------------------------------------*/
239
240/*------------------------------------------------------------------------------------------------------------
241Preprocessing functions
242Functions preparing geometries for distance-calculations
243--------------------------------------------------------------------------------------------------------------*/
244
250int
251lw_dist2d_comp(const LWGEOM *lw1, const LWGEOM *lw2, DISTPTS *dl)
252{
253 return lw_dist2d_recursive(lw1, lw2, dl);
254}
255
256static int
258{
259 /* Differs from lwgeom_is_collection by not treating CURVEPOLYGON as collection */
260 switch (g->type)
261 {
262 case TINTYPE:
263 case MULTIPOINTTYPE:
264 case MULTILINETYPE:
265 case MULTIPOLYGONTYPE:
266 case COLLECTIONTYPE:
267 case MULTICURVETYPE:
268 case MULTISURFACETYPE:
269 case COMPOUNDTYPE:
271 return LW_TRUE;
272 break;
273
274 default:
275 return LW_FALSE;
276 }
277}
278
279
280/* Check for overlapping bboxes */
281static int
282lw_dist2d_check_overlap(const LWGEOM *lwg1, const LWGEOM *lwg2)
283{
284 assert(lwg1 && lwg2 && lwg1->bbox && lwg2->bbox);
285
286 /* Check if the geometries intersect. */
287 if ((lwg1->bbox->xmax < lwg2->bbox->xmin || lwg1->bbox->xmin > lwg2->bbox->xmax ||
288 lwg1->bbox->ymax < lwg2->bbox->ymin || lwg1->bbox->ymin > lwg2->bbox->ymax))
289 {
290 return LW_FALSE;
291 }
292 return LW_TRUE;
293}
294
298int
299lw_dist2d_recursive(const LWGEOM *lwg1, const LWGEOM *lwg2, DISTPTS *dl)
300{
301 int i, j;
302 int n1 = 1;
303 int n2 = 1;
304 LWGEOM *g1 = NULL;
305 LWGEOM *g2 = NULL;
306 LWCOLLECTION *c1 = NULL;
307 LWCOLLECTION *c2 = NULL;
308
309 LWDEBUGF(2, "lw_dist2d_comp is called with type1=%d, type2=%d", lwg1->type, lwg2->type);
310
311 if (lw_dist2d_is_collection(lwg1))
312 {
313 LWDEBUG(3, "First geometry is collection");
314 c1 = lwgeom_as_lwcollection(lwg1);
315 n1 = c1->ngeoms;
316 }
317 if (lw_dist2d_is_collection(lwg2))
318 {
319 LWDEBUG(3, "Second geometry is collection");
320 c2 = lwgeom_as_lwcollection(lwg2);
321 n2 = c2->ngeoms;
322 }
323
324 for (i = 0; i < n1; i++)
325 {
326
327 if (lw_dist2d_is_collection(lwg1))
328 g1 = c1->geoms[i];
329 else
330 g1 = (LWGEOM *)lwg1;
331
332 if (!g1) continue;
333
334 if (lwgeom_is_empty(g1))
335 continue;
336
338 {
339 LWDEBUG(3, "Found collection inside first geometry collection, recursing");
340 if (!lw_dist2d_recursive(g1, lwg2, dl))
341 return LW_FALSE;
342 continue;
343 }
344 for (j = 0; j < n2; j++)
345 {
346 if (lw_dist2d_is_collection(lwg2))
347 g2 = c2->geoms[j];
348 else
349 g2 = (LWGEOM *)lwg2;
350
351 if (!g2) continue;
352
354 {
355 LWDEBUG(3, "Found collection inside second geometry collection, recursing");
356 if (!lw_dist2d_recursive(g1, g2, dl))
357 return LW_FALSE;
358 continue;
359 }
360
361 if (!g1->bbox)
362 lwgeom_add_bbox(g1);
363
364 if (!g2->bbox)
365 lwgeom_add_bbox(g2);
366
367 /* If one of geometries is empty, skip */
368 if (lwgeom_is_empty(g1) || lwgeom_is_empty(g2))
369 continue;
370
371 if ((dl->mode != DIST_MAX) && (!lw_dist2d_check_overlap(g1, g2)) &&
372 (g1->type == LINETYPE || g1->type == POLYGONTYPE || g1->type == TRIANGLETYPE) &&
373 (g2->type == LINETYPE || g2->type == POLYGONTYPE || g2->type == TRIANGLETYPE))
374 {
375 if (!lw_dist2d_distribute_fast(g1, g2, dl))
376 return LW_FALSE;
377 }
378 else
379 {
380 if (!lw_dist2d_distribute_bruteforce(g1, g2, dl))
381 return LW_FALSE;
382 LWDEBUGF(2, "Distance so far: %.15g (%.15g tolerated)", dl->distance, dl->tolerance);
383 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
384 return LW_TRUE; /*just a check if the answer is already given*/
385 LWDEBUG(2, "Not below tolerance yet");
386 }
387 }
388 }
389 return LW_TRUE;
390}
391
392int
394{
395
396 int t1 = lwg1->type;
397 int t2 = lwg2->type;
398
399 switch (t1)
400 {
401 case POINTTYPE:
402 {
403 dl->twisted = 1;
404 switch (t2)
405 {
406 case POINTTYPE:
407 return lw_dist2d_point_point((LWPOINT *)lwg1, (LWPOINT *)lwg2, dl);
408 case LINETYPE:
409 return lw_dist2d_point_line((LWPOINT *)lwg1, (LWLINE *)lwg2, dl);
410 case TRIANGLETYPE:
411 return lw_dist2d_point_tri((LWPOINT *)lwg1, (LWTRIANGLE *)lwg2, dl);
412 case POLYGONTYPE:
413 return lw_dist2d_point_poly((LWPOINT *)lwg1, (LWPOLY *)lwg2, dl);
414 case CIRCSTRINGTYPE:
415 return lw_dist2d_point_circstring((LWPOINT *)lwg1, (LWCIRCSTRING *)lwg2, dl);
416 case CURVEPOLYTYPE:
417 return lw_dist2d_point_curvepoly((LWPOINT *)lwg1, (LWCURVEPOLY *)lwg2, dl);
418 default:
419 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
420 return LW_FALSE;
421 }
422 }
423 case LINETYPE:
424 {
425 dl->twisted = 1;
426 switch (t2)
427 {
428 case POINTTYPE:
429 dl->twisted = -1;
430 return lw_dist2d_point_line((LWPOINT *)lwg2, (LWLINE *)lwg1, dl);
431 case LINETYPE:
432 return lw_dist2d_line_line((LWLINE *)lwg1, (LWLINE *)lwg2, dl);
433 case TRIANGLETYPE:
434 return lw_dist2d_line_tri((LWLINE *)lwg1, (LWTRIANGLE *)lwg2, dl);
435 case POLYGONTYPE:
436 return lw_dist2d_line_poly((LWLINE *)lwg1, (LWPOLY *)lwg2, dl);
437 case CIRCSTRINGTYPE:
438 return lw_dist2d_line_circstring((LWLINE *)lwg1, (LWCIRCSTRING *)lwg2, dl);
439 case CURVEPOLYTYPE:
440 return lw_dist2d_line_curvepoly((LWLINE *)lwg1, (LWCURVEPOLY *)lwg2, dl);
441 default:
442 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
443 return LW_FALSE;
444 }
445 }
446 case TRIANGLETYPE:
447 {
448 dl->twisted = 1;
449 switch (t2)
450 {
451 case POINTTYPE:
452 dl->twisted = -1;
453 return lw_dist2d_point_tri((LWPOINT *)lwg2, (LWTRIANGLE *)lwg1, dl);
454 case LINETYPE:
455 dl->twisted = -1;
456 return lw_dist2d_line_tri((LWLINE *)lwg2, (LWTRIANGLE *)lwg1, dl);
457 case TRIANGLETYPE:
458 return lw_dist2d_tri_tri((LWTRIANGLE *)lwg1, (LWTRIANGLE *)lwg2, dl);
459 case POLYGONTYPE:
460 return lw_dist2d_tri_poly((LWTRIANGLE *)lwg1, (LWPOLY *)lwg2, dl);
461 case CIRCSTRINGTYPE:
462 return lw_dist2d_tri_circstring((LWTRIANGLE *)lwg1, (LWCIRCSTRING *)lwg2, dl);
463 case CURVEPOLYTYPE:
464 return lw_dist2d_tri_curvepoly((LWTRIANGLE *)lwg1, (LWCURVEPOLY *)lwg2, dl);
465 default:
466 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
467 return LW_FALSE;
468 }
469 }
470 case CIRCSTRINGTYPE:
471 {
472 dl->twisted = 1;
473 switch (t2)
474 {
475 case POINTTYPE:
476 dl->twisted = -1;
477 return lw_dist2d_point_circstring((LWPOINT *)lwg2, (LWCIRCSTRING *)lwg1, dl);
478 case LINETYPE:
479 dl->twisted = -1;
480 return lw_dist2d_line_circstring((LWLINE *)lwg2, (LWCIRCSTRING *)lwg1, dl);
481 case TRIANGLETYPE:
482 dl->twisted = -1;
483 return lw_dist2d_tri_circstring((LWTRIANGLE *)lwg2, (LWCIRCSTRING *)lwg1, dl);
484 case POLYGONTYPE:
485 return lw_dist2d_circstring_poly((LWCIRCSTRING *)lwg1, (LWPOLY *)lwg2, dl);
486 case CIRCSTRINGTYPE:
487 return lw_dist2d_circstring_circstring((LWCIRCSTRING *)lwg1, (LWCIRCSTRING *)lwg2, dl);
488 case CURVEPOLYTYPE:
489 return lw_dist2d_circstring_curvepoly((LWCIRCSTRING *)lwg1, (LWCURVEPOLY *)lwg2, dl);
490 default:
491 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
492 return LW_FALSE;
493 }
494 }
495 case POLYGONTYPE:
496 {
497 dl->twisted = -1;
498 switch (t2)
499 {
500 case POINTTYPE:
501 return lw_dist2d_point_poly((LWPOINT *)lwg2, (LWPOLY *)lwg1, dl);
502 case LINETYPE:
503 return lw_dist2d_line_poly((LWLINE *)lwg2, (LWPOLY *)lwg1, dl);
504 case TRIANGLETYPE:
505 return lw_dist2d_tri_poly((LWTRIANGLE *)lwg2, (LWPOLY *)lwg1, dl);
506 case CIRCSTRINGTYPE:
507 return lw_dist2d_circstring_poly((LWCIRCSTRING *)lwg2, (LWPOLY *)lwg1, dl);
508 case POLYGONTYPE:
509 dl->twisted = 1;
510 return lw_dist2d_poly_poly((LWPOLY *)lwg1, (LWPOLY *)lwg2, dl);
511 case CURVEPOLYTYPE:
512 dl->twisted = 1;
513 return lw_dist2d_poly_curvepoly((LWPOLY *)lwg1, (LWCURVEPOLY *)lwg2, dl);
514 default:
515 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
516 return LW_FALSE;
517 }
518 }
519 case CURVEPOLYTYPE:
520 {
521 dl->twisted = -1;
522 switch (t2)
523 {
524 case POINTTYPE:
525 return lw_dist2d_point_curvepoly((LWPOINT *)lwg2, (LWCURVEPOLY *)lwg1, dl);
526 case LINETYPE:
527 return lw_dist2d_line_curvepoly((LWLINE *)lwg2, (LWCURVEPOLY *)lwg1, dl);
528 case TRIANGLETYPE:
529 return lw_dist2d_tri_curvepoly((LWTRIANGLE *)lwg2, (LWCURVEPOLY *)lwg1, dl);
530 case POLYGONTYPE:
531 return lw_dist2d_poly_curvepoly((LWPOLY *)lwg2, (LWCURVEPOLY *)lwg1, dl);
532 case CIRCSTRINGTYPE:
533 return lw_dist2d_circstring_curvepoly((LWCIRCSTRING *)lwg2, (LWCURVEPOLY *)lwg1, dl);
534 case CURVEPOLYTYPE:
535 dl->twisted = 1;
536 return lw_dist2d_curvepoly_curvepoly((LWCURVEPOLY *)lwg1, (LWCURVEPOLY *)lwg2, dl);
537 default:
538 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t2));
539 return LW_FALSE;
540 }
541 }
542 default:
543 {
544 lwerror("%s: Unsupported geometry type: %s", __func__, lwtype_name(t1));
545 return LW_FALSE;
546 }
547 }
548
549 return LW_FALSE;
550}
551
552
554int
556{
557 POINTARRAY *pa1, *pa2;
558 int type1 = lwg1->type;
559 int type2 = lwg2->type;
560
561 LWDEBUGF(2, "lw_dist2d_distribute_fast is called with typ1=%d, type2=%d", lwg1->type, lwg2->type);
562
563 switch (type1)
564 {
565 case LINETYPE:
566 pa1 = ((LWLINE *)lwg1)->points;
567 break;
568 case POLYGONTYPE:
569 pa1 = ((LWPOLY *)lwg1)->rings[0];
570 break;
571 case TRIANGLETYPE:
572 pa1 = ((LWTRIANGLE *)lwg1)->points;
573 break;
574 default:
575 lwerror("Unsupported geometry1 type: %s", lwtype_name(type1));
576 return LW_FALSE;
577 }
578 switch (type2)
579 {
580 case LINETYPE:
581 pa2 = ((LWLINE *)lwg2)->points;
582 break;
583 case POLYGONTYPE:
584 pa2 = ((LWPOLY *)lwg2)->rings[0];
585 break;
586 case TRIANGLETYPE:
587 pa2 = ((LWTRIANGLE *)lwg2)->points;
588 break;
589 default:
590 lwerror("Unsupported geometry2 type: %s", lwtype_name(type1));
591 return LW_FALSE;
592 }
593 dl->twisted = 1;
594 return lw_dist2d_fast_ptarray_ptarray(pa1, pa2, dl, lwg1->bbox, lwg2->bbox);
595}
596
597/*------------------------------------------------------------------------------------------------------------
598End of Preprocessing functions
599--------------------------------------------------------------------------------------------------------------*/
600
601/*------------------------------------------------------------------------------------------------------------
602Brute force functions
603The old way of calculating distances, now used for:
6041) distances to points (because there shouldn't be anything to gain by the new way of doing it)
6052) distances when subgeometries geometries bboxes overlaps
606--------------------------------------------------------------------------------------------------------------*/
607
611int
613{
614 const POINT2D *p1 = getPoint2d_cp(point1->point, 0);
615 const POINT2D *p2 = getPoint2d_cp(point2->point, 0);
616 return lw_dist2d_pt_pt(p1, p2, dl);
617}
621int
623{
624 const POINT2D *p = getPoint2d_cp(point->point, 0);
625 return lw_dist2d_pt_ptarray(p, line->points, dl);
626}
627
628int
630{
631 const POINT2D *pt = getPoint2d_cp(point->point, 0);
632 /* Is point inside triangle? */
633 if (dl->mode == DIST_MIN && ptarray_contains_point(tri->points, pt) != LW_OUTSIDE)
634 {
635 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
636 return LW_TRUE;
637 }
638
639 return lw_dist2d_pt_ptarray(pt, tri->points, dl);
640}
641
642int
644{
645 const POINT2D *p;
646 p = getPoint2d_cp(point->point, 0);
647 return lw_dist2d_pt_ptarrayarc(p, circ->points, dl);
648}
649
655int
657{
658 const POINT2D *p = getPoint2d_cp(point->point, 0);
659
660 /* Max distance? Check only outer ring.*/
661 if (dl->mode == DIST_MAX)
662 return lw_dist2d_pt_ptarray(p, poly->rings[0], dl);
663
664 /* Return distance to outer ring if not inside it */
665 if (ptarray_contains_point(poly->rings[0], p) == LW_OUTSIDE)
666 return lw_dist2d_pt_ptarray(p, poly->rings[0], dl);
667
668 /*
669 * Inside the outer ring.
670 * Scan though each of the inner rings looking to
671 * see if its inside. If not, distance==0.
672 * Otherwise, distance = pt to ring distance
673 */
674 for (uint32_t i = 1; i < poly->nrings; i++)
675 if (ptarray_contains_point(poly->rings[i], p) != LW_OUTSIDE)
676 return lw_dist2d_pt_ptarray(p, poly->rings[i], dl);
677
678 /* Is inside the polygon */
679 lw_dist2d_distpts_set(dl, 0.0, p, p);
680 return LW_TRUE;
681}
682
683int
685{
686 const POINT2D *p = getPoint2d_cp(point->point, 0);
687
688 if (dl->mode == DIST_MAX)
689 lwerror("lw_dist2d_point_curvepoly cannot calculate max distance");
690
691 /* Return distance to outer ring if not inside it */
692 if (lwgeom_contains_point(poly->rings[0], p) == LW_OUTSIDE)
693 return lw_dist2d_recursive((LWGEOM *)point, poly->rings[0], dl);
694
695 /* Inside the outer ring.
696 * Scan though each of the inner rings looking to see if its inside. If not, distance==0.
697 * Otherwise, distance = pt to ring distance.
698 */
699 for (uint32_t i = 1; i < poly->nrings; i++)
700 if (lwgeom_contains_point(poly->rings[i], p) == LW_INSIDE)
701 return lw_dist2d_recursive((LWGEOM *)point, poly->rings[i], dl);
702
703 /* Is inside the polygon */
704 lw_dist2d_distpts_set(dl, 0.0, p, p);
705 return LW_TRUE;
706}
707
708int
710{
711 POINTARRAY *pa1 = line1->points;
712 POINTARRAY *pa2 = line2->points;
713 return lw_dist2d_ptarray_ptarray(pa1, pa2, dl);
714}
715
716int
718{
719 const POINT2D *pt = getPoint2d_cp(line->points, 0);
720 /* Is there a point inside triangle? */
721 if (dl->mode == DIST_MIN && ptarray_contains_point(tri->points, pt) != LW_OUTSIDE)
722 {
723 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
724 return LW_TRUE;
725 }
726
727 return lw_dist2d_ptarray_ptarray(line->points, tri->points, dl);
728}
729
730int
732{
733 return lw_dist2d_ptarray_ptarrayarc(line1->points, line2->points, dl);
734}
735
747int
749{
750 POINTARRAY *pa = line->points;
751 const POINT2D *pt = getPoint2d_cp(pa, 0);
752
753 /* Line has a point outside poly. Check distance to outer ring only. */
754 if (ptarray_contains_point(poly->rings[0], pt) == LW_OUTSIDE || dl->mode == DIST_MAX)
755 return lw_dist2d_ptarray_ptarray(pa, poly->rings[0], dl);
756
757 for (uint32_t i = 1; i < poly->nrings; i++)
758 {
759 if (!lw_dist2d_ptarray_ptarray(pa, poly->rings[i], dl))
760 return LW_FALSE;
761
762 /* just a check if the answer is already given */
763 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
764 return LW_TRUE;
765 }
766
767 /* It's inside a hole, then the actual distance is the min ring distance */
768 for (uint32_t i = 1; i < poly->nrings; i++)
769 if (ptarray_contains_point(poly->rings[i], pt) != LW_OUTSIDE)
770 return LW_TRUE;
771
772 /* Not in hole, so inside polygon */
773 if (dl->mode == DIST_MIN)
774 {
775 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
776 }
777 return LW_TRUE;
778}
779
780int
782{
783 const POINT2D *pt = getPoint2d_cp(line->points, 0);
784
785 /* Line has a point outside curvepoly. Check distance to outer ring only. */
786 if (lwgeom_contains_point(poly->rings[0], pt) == LW_OUTSIDE)
787 return lw_dist2d_recursive((LWGEOM *)line, poly->rings[0], dl);
788
789 for (uint32_t i = 1; i < poly->nrings; i++)
790 {
791 if (!lw_dist2d_recursive((LWGEOM *)line, poly->rings[i], dl))
792 return LW_FALSE;
793
794 /* just a check if the answer is already given */
795 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
796 return LW_TRUE;
797 }
798
799 /* It's inside a hole, then distance is actual min distance */
800 for (uint32_t i = 1; i < poly->nrings; i++)
801 if (lwgeom_contains_point(poly->rings[i], pt) != LW_OUTSIDE)
802 return LW_TRUE;
803
804 /* Not in hole, so inside polygon */
805 if (dl->mode == DIST_MIN)
806 {
807 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
808 }
809 return LW_TRUE;
810}
811
812int
814{
815 POINTARRAY *pa1 = tri1->points;
816 POINTARRAY *pa2 = tri2->points;
817 const POINT2D *pt = getPoint2d_cp(pa2, 0);
818 if (dl->mode == DIST_MIN && ptarray_contains_point(pa1, pt) != LW_OUTSIDE)
819 {
820 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
821 return LW_TRUE;
822 }
823
824 pt = getPoint2d_cp(pa1, 0);
825 if (dl->mode == DIST_MIN && ptarray_contains_point(pa2, pt) != LW_OUTSIDE)
826 {
827 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
828 return LW_TRUE;
829 }
830
831 return lw_dist2d_ptarray_ptarray(pa1, pa2, dl);
832}
833
834int
836{
837 POINTARRAY *pa = tri->points;
838 const POINT2D *pt = getPoint2d_cp(pa, 0);
839
840 /* If we are looking for maxdistance, just check the outer rings.*/
841 if (dl->mode == DIST_MAX)
842 return lw_dist2d_ptarray_ptarray(pa, poly->rings[0], dl);
843
844 /* Triangle has a point outside poly. Check distance to outer ring only. */
845 if (ptarray_contains_point(poly->rings[0], pt) == LW_OUTSIDE)
846 {
847 if (!lw_dist2d_ptarray_ptarray(pa, poly->rings[0], dl))
848 return LW_FALSE;
849
850 /* just a check if the answer is already given */
851 if (dl->distance <= dl->tolerance)
852 return LW_TRUE;
853
854 /* Maybe poly is inside triangle? */
855 const POINT2D *pt2 = getPoint2d_cp(poly->rings[0], 0);
856 if (ptarray_contains_point(pa, pt2) != LW_OUTSIDE)
857 {
858 lw_dist2d_distpts_set(dl, 0.0, pt2, pt2);
859 return LW_TRUE;
860 }
861 }
862
863 for (uint32_t i = 1; i < poly->nrings; i++)
864 {
865 if (!lw_dist2d_ptarray_ptarray(pa, poly->rings[i], dl))
866 return LW_FALSE;
867
868 /* just a check if the answer is already given */
869 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
870 return LW_TRUE;
871 }
872
873 /* It's inside a hole, then the actual distance is the min ring distance */
874 for (uint32_t i = 1; i < poly->nrings; i++)
875 if (ptarray_contains_point(poly->rings[i], pt) != LW_OUTSIDE)
876 return LW_TRUE;
877
878 /* Not in hole, so inside polygon */
879 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
880 return LW_TRUE;
881}
882static const POINT2D *
884{
885 switch (geom->type)
886 {
887 case LINETYPE:
888 return getPoint2d_cp(((LWLINE *)geom)->points, 0);
889 case CIRCSTRINGTYPE:
890 return getPoint2d_cp(((LWCIRCSTRING *)geom)->points, 0);
891 case COMPOUNDTYPE:
892 {
893 LWCOMPOUND *comp = (LWCOMPOUND *)geom;
894 LWLINE *line = (LWLINE *)(comp->geoms[0]);
895 return getPoint2d_cp(line->points, 0);
896 }
897 default:
898 lwerror("lw_curvering_getfirstpoint2d_cp: unknown type");
899 }
900 return NULL;
901}
902
903int
905{
906 const POINT2D *pt = getPoint2d_cp(tri->points, 0);
907
908 /* If we are looking for maxdistance, just check the outer rings.*/
909 if (dl->mode == DIST_MAX)
910 return lw_dist2d_recursive((LWGEOM *)tri, poly->rings[0], dl);
911
912 /* Line has a point outside curvepoly. Check distance to outer ring only. */
913 if (lwgeom_contains_point(poly->rings[0], pt) == LW_OUTSIDE)
914 {
915 if (lw_dist2d_recursive((LWGEOM *)tri, poly->rings[0], dl))
916 return LW_TRUE;
917 /* Maybe poly is inside triangle? */
919 {
920 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
921 return LW_TRUE;
922 }
923 }
924
925 for (uint32_t i = 1; i < poly->nrings; i++)
926 {
927 if (!lw_dist2d_recursive((LWGEOM *)tri, poly->rings[i], dl))
928 return LW_FALSE;
929
930 /* just a check if the answer is already given */
931 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
932 return LW_TRUE;
933 }
934
935 /* It's inside a hole, then distance is actual min distance */
936 for (uint32_t i = 1; i < poly->nrings; i++)
937 if (lwgeom_contains_point(poly->rings[i], pt) != LW_OUTSIDE)
938 return LW_TRUE;
939
940 /* Not in hole, so inside polygon */
941 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
942 return LW_TRUE;
943}
944
945int
947{
949 if (ptarray_contains_point(tri->points, pt) != LW_OUTSIDE && dl->mode == DIST_MIN)
950 {
951 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
952 return LW_TRUE;
953 }
954
955 return lw_dist2d_ptarray_ptarrayarc(tri->points, line->points, dl);
956}
957
967int
969{
970 const POINT2D *pt;
971
972 LWDEBUG(2, "lw_dist2d_poly_poly called");
973
974 /*1 if we are looking for maxdistance, just check the outer rings.*/
975 if (dl->mode == DIST_MAX)
976 return lw_dist2d_ptarray_ptarray(poly1->rings[0], poly2->rings[0], dl);
977
978 /* 2 check if poly1 has first point outside poly2 and vice versa, if so, just check outer rings
979 here it would be possible to handle the information about which one is inside which one and only search for the
980 smaller ones in the bigger ones holes.*/
981 pt = getPoint2d_cp(poly1->rings[0], 0);
982 if (ptarray_contains_point(poly2->rings[0], pt) == LW_OUTSIDE)
983 {
984 pt = getPoint2d_cp(poly2->rings[0], 0);
985 if (ptarray_contains_point(poly1->rings[0], pt) == LW_OUTSIDE)
986 return lw_dist2d_ptarray_ptarray(poly1->rings[0], poly2->rings[0], dl);
987 }
988
989 /*3 check if first point of poly2 is in a hole of poly1. If so check outer ring of poly2 against that hole
990 * of poly1*/
991 pt = getPoint2d_cp(poly2->rings[0], 0);
992 for (uint32_t i = 1; i < poly1->nrings; i++)
993 if (ptarray_contains_point(poly1->rings[i], pt) != LW_OUTSIDE)
994 return lw_dist2d_ptarray_ptarray(poly1->rings[i], poly2->rings[0], dl);
995
996 /*4 check if first point of poly1 is in a hole of poly2. If so check outer ring of poly1 against that hole
997 * of poly2*/
998 pt = getPoint2d_cp(poly1->rings[0], 0);
999 for (uint32_t i = 1; i < poly2->nrings; i++)
1000 if (ptarray_contains_point(poly2->rings[i], pt) != LW_OUTSIDE)
1001 return lw_dist2d_ptarray_ptarray(poly1->rings[0], poly2->rings[i], dl);
1002
1003 /*5 If we have come all the way here we know that the first point of one of them is inside the other ones
1004 * outer ring and not in holes so we check which one is inside.*/
1005 pt = getPoint2d_cp(poly1->rings[0], 0);
1006 if (ptarray_contains_point(poly2->rings[0], pt) != LW_OUTSIDE)
1007 {
1008 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
1009 return LW_TRUE;
1010 }
1011
1012 pt = getPoint2d_cp(poly2->rings[0], 0);
1013 if (ptarray_contains_point(poly1->rings[0], pt) != LW_OUTSIDE)
1014 {
1015 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
1016 return LW_TRUE;
1017 }
1018
1019 lwerror("Unspecified error in function lw_dist2d_poly_poly");
1020 return LW_FALSE;
1021}
1022
1023int
1025{
1027 int rv = lw_dist2d_curvepoly_curvepoly(curvepoly1, curvepoly2, dl);
1028 lwgeom_free((LWGEOM *)curvepoly1);
1029 return rv;
1030}
1031
1032int
1034{
1036 int rv = lw_dist2d_line_curvepoly((LWLINE *)circ, curvepoly, dl);
1037 lwgeom_free((LWGEOM *)curvepoly);
1038 return rv;
1039}
1040
1041int
1043{
1044 return lw_dist2d_line_curvepoly((LWLINE *)circ, poly, dl);
1045}
1046
1047int
1052
1053int
1055{
1056 const POINT2D *pt;
1057
1058 /*1 if we are looking for maxdistance, just check the outer rings.*/
1059 if (dl->mode == DIST_MAX)
1060 return lw_dist2d_recursive(poly1->rings[0], poly2->rings[0], dl);
1061
1062 /* 2 check if poly1 has first point outside poly2 and vice versa, if so, just check outer rings
1063 here it would be possible to handle the information about which one is inside which one and only search for the
1064 smaller ones in the bigger ones holes.*/
1065 pt = lw_curvering_getfirstpoint2d_cp(poly1->rings[0]);
1066 if (lwgeom_contains_point(poly2->rings[0], pt) == LW_OUTSIDE)
1067 {
1068 pt = lw_curvering_getfirstpoint2d_cp(poly2->rings[0]);
1069 if (lwgeom_contains_point(poly1->rings[0], pt) == LW_OUTSIDE)
1070 return lw_dist2d_recursive(poly1->rings[0], poly2->rings[0], dl);
1071 }
1072
1073 /*3 check if first point of poly2 is in a hole of poly1. If so check outer ring of poly2 against that hole
1074 * of poly1*/
1075 pt = lw_curvering_getfirstpoint2d_cp(poly2->rings[0]);
1076 for (uint32_t i = 1; i < poly1->nrings; i++)
1077 if (lwgeom_contains_point(poly1->rings[i], pt) != LW_OUTSIDE)
1078 return lw_dist2d_recursive(poly1->rings[i], poly2->rings[0], dl);
1079
1080 /*4 check if first point of poly1 is in a hole of poly2. If so check outer ring of poly1 against that hole
1081 * of poly2*/
1082 pt = lw_curvering_getfirstpoint2d_cp(poly1->rings[0]);
1083 for (uint32_t i = 1; i < poly2->nrings; i++)
1084 if (lwgeom_contains_point(poly2->rings[i], pt) != LW_OUTSIDE)
1085 return lw_dist2d_recursive(poly1->rings[0], poly2->rings[i], dl);
1086
1087 /*5 If we have come all the way here we know that the first point of one of them is inside the other ones
1088 * outer ring and not in holes so we check which one is inside.*/
1089 pt = lw_curvering_getfirstpoint2d_cp(poly1->rings[0]);
1090 if (lwgeom_contains_point(poly2->rings[0], pt) != LW_OUTSIDE)
1091 {
1092 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
1093 return LW_TRUE;
1094 }
1095
1096 pt = lw_curvering_getfirstpoint2d_cp(poly2->rings[0]);
1097 if (lwgeom_contains_point(poly1->rings[0], pt) != LW_OUTSIDE)
1098 {
1099 lw_dist2d_distpts_set(dl, 0.0, pt, pt);
1100 return LW_TRUE;
1101 }
1102
1103 lwerror("Unspecified error in function lw_dist2d_curvepoly_curvepoly");
1104 return LW_FALSE;
1105}
1106
1111int
1113{
1114 const POINT2D *start, *end;
1115 int twist = dl->twisted;
1116
1117 start = getPoint2d_cp(pa, 0);
1118
1119 LWDEBUG(2, "lw_dist2d_pt_ptarray enter");
1120
1121 if (!lw_dist2d_pt_pt(p, start, dl))
1122 return LW_FALSE;
1123
1124 LWDEBUGF(2, "lw_dist2d_pt_ptarray: distance from first point ? : %.15g", dl->distance);
1125
1126 for (uint32_t t = 1; t < pa->npoints; t++)
1127 {
1128 dl->twisted = twist;
1129 end = getPoint2d_cp(pa, t);
1130 if (!lw_dist2d_pt_seg(p, start, end, dl))
1131 return LW_FALSE;
1132
1133 LWDEBUGF(2, "lw_dist2d_pt_ptarray: distance from seg %u ? : %.15g", t, dl->distance);
1134
1135 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
1136 return LW_TRUE; /*just a check if the answer is already given*/
1137 start = end;
1138 }
1139
1140 return LW_TRUE;
1141}
1142
1147int
1149{
1150 uint32_t t;
1151 const POINT2D *A1;
1152 const POINT2D *A2;
1153 const POINT2D *A3;
1154 int twist = dl->twisted;
1155
1156 LWDEBUG(2, "lw_dist2d_pt_ptarrayarc is called");
1157
1158 if (pa->npoints % 2 == 0 || pa->npoints < 3)
1159 {
1160 lwerror("lw_dist2d_pt_ptarrayarc called with non-arc input");
1161 return LW_FALSE;
1162 }
1163
1164 if (dl->mode == DIST_MAX)
1165 {
1166 lwerror("lw_dist2d_pt_ptarrayarc does not currently support DIST_MAX mode");
1167 return LW_FALSE;
1168 }
1169
1170 A1 = getPoint2d_cp(pa, 0);
1171
1172 if (!lw_dist2d_pt_pt(p, A1, dl))
1173 return LW_FALSE;
1174
1175 for (t = 1; t < pa->npoints; t += 2)
1176 {
1177 dl->twisted = twist;
1178 A2 = getPoint2d_cp(pa, t);
1179 A3 = getPoint2d_cp(pa, t + 1);
1180
1181 if (lw_dist2d_pt_arc(p, A1, A2, A3, dl) == LW_FALSE)
1182 return LW_FALSE;
1183
1184 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
1185 return LW_TRUE; /*just a check if the answer is already given*/
1186
1187 A1 = A3;
1188 }
1189
1190 return LW_TRUE;
1191}
1192
1196int
1198{
1199 uint32_t t, u;
1200 const POINT2D *start, *end;
1201 const POINT2D *start2, *end2;
1202 int twist = dl->twisted;
1203
1204 LWDEBUGF(2, "lw_dist2d_ptarray_ptarray called (points: %d-%d)", l1->npoints, l2->npoints);
1205
1206 /* If we are searching for maxdistance we go straight to point-point calculation since the maxdistance have
1207 * to be between two vertices*/
1208 if (dl->mode == DIST_MAX)
1209 {
1210 for (t = 0; t < l1->npoints; t++) /*for each segment in L1 */
1211 {
1212 start = getPoint2d_cp(l1, t);
1213 for (u = 0; u < l2->npoints; u++) /*for each segment in L2 */
1214 {
1215 start2 = getPoint2d_cp(l2, u);
1216 lw_dist2d_pt_pt(start, start2, dl);
1217 }
1218 }
1219 }
1220 else
1221 {
1222 start = getPoint2d_cp(l1, 0);
1223 for (t = 1; t < l1->npoints; t++) /*for each segment in L1 */
1224 {
1225 end = getPoint2d_cp(l1, t);
1226 start2 = getPoint2d_cp(l2, 0);
1227 for (u = 1; u < l2->npoints; u++) /*for each segment in L2 */
1228 {
1229 end2 = getPoint2d_cp(l2, u);
1230 dl->twisted = twist;
1231 lw_dist2d_seg_seg(start, end, start2, end2, dl);
1232 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
1233 return LW_TRUE; /*just a check if the answer is already given*/
1234 start2 = end2;
1235 }
1236 start = end;
1237 }
1238 }
1239 return LW_TRUE;
1240}
1241
1245int
1247{
1248 uint32_t t, u;
1249 const POINT2D *A1;
1250 const POINT2D *A2;
1251 const POINT2D *B1;
1252 const POINT2D *B2;
1253 const POINT2D *B3;
1254 int twist = dl->twisted;
1255
1256 LWDEBUGF(2, "lw_dist2d_ptarray_ptarrayarc called (points: %d-%d)", pa->npoints, pb->npoints);
1257
1258 if (pb->npoints % 2 == 0 || pb->npoints < 3)
1259 {
1260 lwerror("lw_dist2d_ptarray_ptarrayarc called with non-arc input");
1261 return LW_FALSE;
1262 }
1263
1264 if (dl->mode == DIST_MAX)
1265 {
1266 lwerror("lw_dist2d_ptarray_ptarrayarc does not currently support DIST_MAX mode");
1267 return LW_FALSE;
1268 }
1269 else
1270 {
1271 A1 = getPoint2d_cp(pa, 0);
1272 for (t = 1; t < pa->npoints; t++) /* For each segment in pa */
1273 {
1274 A2 = getPoint2d_cp(pa, t);
1275 B1 = getPoint2d_cp(pb, 0);
1276 for (u = 1; u < pb->npoints; u += 2) /* For each arc in pb */
1277 {
1278 B2 = getPoint2d_cp(pb, u);
1279 B3 = getPoint2d_cp(pb, u + 1);
1280 dl->twisted = twist;
1281
1282 lw_dist2d_seg_arc(A1, A2, B1, B2, B3, dl);
1283
1284 /* If we've found a distance within tolerance, we're done */
1285 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
1286 return LW_TRUE;
1287
1288 B1 = B3;
1289 }
1290 A1 = A2;
1291 }
1292 }
1293 return LW_TRUE;
1294}
1295
1299int
1301{
1302 uint32_t t, u;
1303 const POINT2D *A1;
1304 const POINT2D *A2;
1305 const POINT2D *A3;
1306 const POINT2D *B1;
1307 const POINT2D *B2;
1308 const POINT2D *B3;
1309 int twist = dl->twisted;
1310
1311 LWDEBUGF(2, "lw_dist2d_ptarrayarc_ptarrayarc called (points: %d-%d)", pa->npoints, pb->npoints);
1312
1313 if (dl->mode == DIST_MAX)
1314 {
1315 lwerror("lw_dist2d_ptarrayarc_ptarrayarc does not currently support DIST_MAX mode");
1316 return LW_FALSE;
1317 }
1318 else
1319 {
1320 A1 = getPoint2d_cp(pa, 0);
1321 for (t = 1; t < pa->npoints; t += 2) /* For each segment in pa */
1322 {
1323 A2 = getPoint2d_cp(pa, t);
1324 A3 = getPoint2d_cp(pa, t + 1);
1325 B1 = getPoint2d_cp(pb, 0);
1326 for (u = 1; u < pb->npoints; u += 2) /* For each arc in pb */
1327 {
1328 B2 = getPoint2d_cp(pb, u);
1329 B3 = getPoint2d_cp(pb, u + 1);
1330 dl->twisted = twist;
1331
1332 lw_dist2d_arc_arc(A1, A2, A3, B1, B2, B3, dl);
1333
1334 /* If we've found a distance within tolerance, we're done */
1335 if (dl->distance <= dl->tolerance && dl->mode == DIST_MIN)
1336 return LW_TRUE;
1337
1338 B1 = B3;
1339 }
1340 A1 = A3;
1341 }
1342 }
1343 return LW_TRUE;
1344}
1345
1350int
1352 const POINT2D *A2,
1353 const POINT2D *B1,
1354 const POINT2D *B2,
1355 const POINT2D *B3,
1356 DISTPTS *dl)
1357{
1358 POINT2D C; /* center of arc circle */
1359 double radius_C; /* radius of arc circle */
1360 POINT2D D; /* point on A closest to C */
1361 double dist_C_D; /* distance from C to D */
1362 int pt_in_arc, pt_in_seg;
1363 DISTPTS dltmp;
1364
1365 /* Bail out on crazy modes */
1366 if (dl->mode < 0)
1367 lwerror("lw_dist2d_seg_arc does not support maxdistance mode");
1368
1369 /* What if the "arc" is a point? */
1370 if (lw_arc_is_pt(B1, B2, B3))
1371 return lw_dist2d_pt_seg(B1, A1, A2, dl);
1372
1373 /* Calculate center and radius of the circle. */
1374 radius_C = lw_arc_center(B1, B2, B3, &C);
1375
1376 /* This "arc" is actually a line (B2 is collinear with B1,B3) */
1377 if (radius_C < 0.0)
1378 return lw_dist2d_seg_seg(A1, A2, B1, B3, dl);
1379
1380 /* Calculate distance between the line and circle center */
1382 if (lw_dist2d_pt_seg(&C, A1, A2, &dltmp) == LW_FALSE)
1383 lwerror("lw_dist2d_pt_seg failed in lw_dist2d_seg_arc");
1384
1385 D = dltmp.p1;
1386 dist_C_D = dltmp.distance;
1387
1388 /* Line intersects circle, maybe arc intersects edge? */
1389 /* If so, that's the closest point. */
1390 /* If not, the closest point is one of the end points of A */
1391 if (dist_C_D < radius_C)
1392 {
1393 double length_A; /* length of the segment A */
1394 POINT2D E, F; /* points of intersection of edge A and circle(B) */
1395 double dist_D_EF; /* distance from D to E or F (same distance both ways) */
1396
1397 dist_D_EF = sqrt(radius_C * radius_C - dist_C_D * dist_C_D);
1398 length_A = sqrt((A2->x - A1->x) * (A2->x - A1->x) + (A2->y - A1->y) * (A2->y - A1->y));
1399
1400 /* Point of intersection E */
1401 E.x = D.x - (A2->x - A1->x) * dist_D_EF / length_A;
1402 E.y = D.y - (A2->y - A1->y) * dist_D_EF / length_A;
1403 /* Point of intersection F */
1404 F.x = D.x + (A2->x - A1->x) * dist_D_EF / length_A;
1405 F.y = D.y + (A2->y - A1->y) * dist_D_EF / length_A;
1406
1407 /* If E is within A and within B then it's an intersection point */
1408 pt_in_arc = lw_pt_in_arc(&E, B1, B2, B3);
1409 pt_in_seg = lw_pt_in_seg(&E, A1, A2);
1410
1411 if (pt_in_arc && pt_in_seg)
1412 {
1413 lw_dist2d_distpts_set(dl, 0.0, &E, &E);
1414 return LW_TRUE;
1415 }
1416
1417 /* If F is within A and within B then it's an intersection point */
1418 pt_in_arc = lw_pt_in_arc(&F, B1, B2, B3);
1419 pt_in_seg = lw_pt_in_seg(&F, A1, A2);
1420
1421 if (pt_in_arc && pt_in_seg)
1422 {
1423 lw_dist2d_distpts_set(dl, 0.0, &F, &F);
1424 return LW_TRUE;
1425 }
1426 }
1427
1428 /* Line grazes circle, maybe arc intersects edge? */
1429 /* If so, grazing point is the closest point. */
1430 /* If not, the closest point is one of the end points of A */
1431 else if (dist_C_D == radius_C)
1432 {
1433 /* Closest point D is also the point of grazing */
1434 pt_in_arc = lw_pt_in_arc(&D, B1, B2, B3);
1435 pt_in_seg = lw_pt_in_seg(&D, A1, A2);
1436
1437 /* Is D contained in both A and B? */
1438 if (pt_in_arc && pt_in_seg)
1439 {
1440 lw_dist2d_distpts_set(dl, 0.0, &D, &D);
1441 return LW_TRUE;
1442 }
1443 }
1444 /* Line misses circle. */
1445 /* If closest point to A on circle is within B, then that's the closest */
1446 /* Otherwise, the closest point will be an end point of A */
1447 else
1448 {
1449 POINT2D G; /* Point on circle closest to A */
1450 G.x = C.x + (D.x - C.x) * radius_C / dist_C_D;
1451 G.y = C.y + (D.y - C.y) * radius_C / dist_C_D;
1452
1453 pt_in_arc = lw_pt_in_arc(&G, B1, B2, B3);
1454 pt_in_seg = lw_pt_in_seg(&D, A1, A2);
1455
1456 /* Closest point is on the interior of A and B */
1457 if (pt_in_arc && pt_in_seg)
1458 return lw_dist2d_pt_pt(&D, &G, dl);
1459 }
1460
1461 /* Now we test the many combinations of end points with either */
1462 /* arcs or edges. Each previous check determined if the closest */
1463 /* potential point was within the arc/segment inscribed on the */
1464 /* line/circle holding the arc/segment. */
1465
1466 /* Closest point is in the arc, but not in the segment, so */
1467 /* one of the segment end points must be the closest. */
1468 if (pt_in_arc && !pt_in_seg)
1469 {
1470 lw_dist2d_pt_arc(A1, B1, B2, B3, dl);
1471 lw_dist2d_pt_arc(A2, B1, B2, B3, dl);
1472 return LW_TRUE;
1473 }
1474 /* or, one of the arc end points is the closest */
1475 else if (pt_in_seg && !pt_in_arc)
1476 {
1477 lw_dist2d_pt_seg(B1, A1, A2, dl);
1478 lw_dist2d_pt_seg(B3, A1, A2, dl);
1479 return LW_TRUE;
1480 }
1481 /* Finally, one of the end-point to end-point combos is the closest. */
1482 else
1483 {
1484 lw_dist2d_pt_pt(A1, B1, dl);
1485 lw_dist2d_pt_pt(A1, B3, dl);
1486 lw_dist2d_pt_pt(A2, B1, dl);
1487 lw_dist2d_pt_pt(A2, B3, dl);
1488 return LW_TRUE;
1489 }
1490
1491 return LW_FALSE;
1492}
1493
1494int
1495lw_dist2d_pt_arc(const POINT2D *P, const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, DISTPTS *dl)
1496{
1497 double radius_A, d;
1498 POINT2D C; /* center of circle defined by arc A */
1499 POINT2D X; /* point circle(A) where line from C to P crosses */
1500
1501 if (dl->mode < 0)
1502 lwerror("lw_dist2d_pt_arc does not support maxdistance mode");
1503
1504 /* What if the arc is a point? */
1505 if (lw_arc_is_pt(A1, A2, A3))
1506 return lw_dist2d_pt_pt(P, A1, dl);
1507
1508 /* Calculate centers and radii of circles. */
1509 radius_A = lw_arc_center(A1, A2, A3, &C);
1510
1511 /* This "arc" is actually a line (A2 is colinear with A1,A3) */
1512 if (radius_A < 0.0)
1513 return lw_dist2d_pt_seg(P, A1, A3, dl);
1514
1515 /* Distance from point to center */
1516 d = distance2d_pt_pt(&C, P);
1517
1518 /* P is the center of the circle */
1519 if (FP_EQUALS(d, 0.0))
1520 {
1521 lw_dist2d_distpts_set(dl, radius_A, A1, P);
1522 return LW_TRUE;
1523 }
1524
1525 /* X is the point on the circle where the line from P to C crosses */
1526 X.x = C.x + (P->x - C.x) * radius_A / d;
1527 X.y = C.y + (P->y - C.y) * radius_A / d;
1528
1529 /* Is crossing point inside the arc? Or arc is actually circle? */
1530 if (p2d_same(A1, A3) || lw_pt_in_arc(&X, A1, A2, A3))
1531 {
1532 lw_dist2d_pt_pt(P, &X, dl);
1533 }
1534 else
1535 {
1536 /* Distance is the minimum of the distances to the arc end points */
1537 lw_dist2d_pt_pt(A1, P, dl);
1538 lw_dist2d_pt_pt(A3, P, dl);
1539 }
1540 return LW_TRUE;
1541}
1542
1543
1544
1565static int
1567 const POINT2D *A1, const POINT2D *A2, const POINT2D *A3,
1568 const POINT2D *center_A,
1569 double radius_A,
1570 const POINT2D *P,
1571 POINT2D *I,
1572 uint32_t *ni)
1573{
1574 POINT2D R;
1575
1576 // If the test point is on the center of the other
1577 // arc, some other point has to be closer, by definition.
1578 if (p2d_same(center_A, P))
1579 return 0;
1580
1581 // Calculate vector from the center to the pt
1582 double dir_x = center_A->x - P->x;
1583 double dir_y = center_A->y - P->y;
1584 double dist = sqrt(dir_x * dir_x + dir_y * dir_y);
1585
1586 // Normalize the direction vector to get a unit vector (length = 1)
1587 double unit_x = dir_x / dist;
1588 double unit_y = dir_y / dist;
1589
1590 // Calculate the two intersection points on the circle
1591 // Point 1: Move from center_A along the unit vector by distance radius_A
1592 R.x = center_A->x + unit_x * radius_A;
1593 R.y = center_A->y + unit_y * radius_A;
1594 if (lw_pt_in_arc(&R, A1, A2, A3))
1595 I[(*ni)++] = R;
1596
1597 // Point 2: Move from center_A against the unit vector by distance radius_A
1598 R.x = center_A->x - unit_x * radius_A;
1599 R.y = center_A->y - unit_y * radius_A;
1600 if (lw_pt_in_arc(&R, A1, A2, A3))
1601 I[(*ni)++] = R;
1602
1603 return 0;
1604}
1605
1606
1625static uint32_t
1627 const POINT2D *cA, double rA,
1628 const POINT2D *cB, double rB,
1629 POINT2D *I)
1630{
1631 // Vector from center A to center B
1632 double dx = cB->x - cA->x;
1633 double dy = cB->y - cA->y;
1634
1635 // Distance between the centers
1636 double d = sqrt(dx * dx + dy * dy);
1637
1638 // 'a' is the distance from the center of circle A to the point P,
1639 // which is the base of the right triangle formed by cA, P, and an intersection point.
1640 double a = (rA * rA - rB * rB + d * d) / (2.0 * d);
1641
1642 // 'h' is the height of that right triangle.
1643 double h_squared = rA * rA - a * a;
1644
1645 // Due to floating point errors, h_squared can be slightly negative.
1646 // This happens when the circles are perfectly tangent. Clamp to 0.
1647 if (h_squared < 0.0)
1648 h_squared = 0.0;
1649
1650 double h = sqrt(h_squared);
1651
1652 // Find the coordinates of point P
1653 double Px = cA->x + a * (dx / d);
1654 double Py = cA->y + a * (dy / d);
1655
1656 // The two intersection points are found by moving from P by a distance 'h'
1657 // in directions perpendicular to the line connecting the centers.
1658 // The perpendicular vector to (dx, dy) is (-dy, dx).
1659
1660 // Intersection point 1
1661 I[0].x = Px - h * (dy / d);
1662 I[0].y = Py + h * (dx / d);
1663
1664 // If h is very close to 0, the circles are tangent and there's only one intersection point.
1665 if (FP_IS_ZERO(h))
1666 return 1;
1667
1668 // Intersection point 2
1669 I[1].x = Px + h * (dy / d);
1670 I[1].y = Py - h * (dx / d);
1671
1672 return 2;
1673}
1674
1675
1676int
1678 const POINT2D *A1, const POINT2D *A2, const POINT2D *A3,
1679 const POINT2D *B1, const POINT2D *B2, const POINT2D *B3,
1680 DISTPTS *dl)
1681{
1682 POINT2D pts_A[16]; /* Points on A that might be the nearest */
1683 POINT2D pts_B[16]; /* Points on B that might be the nearest */
1684 uint32_t ai = 0, bi = 0; /* Number of points in pts_A and pts_B */
1685 POINT2D center_A, center_B; /* Center points of arcs A and B */
1686 double radius_A, radius_B, d; /* Radii of arcs A and B */
1687 int is_disjoint, is_overlapping, is_contained, is_same_center;
1688 POINT2D intersectionPts[2];
1689
1690 if (dl->mode != DIST_MIN)
1691 lwerror("lw_dist2d_arc_arc only supports mindistance");
1692
1693 /* What if one or both of our "arcs" is actually a point? */
1694 if (lw_arc_is_pt(B1, B2, B3) && lw_arc_is_pt(A1, A2, A3))
1695 return lw_dist2d_pt_pt(B1, A1, dl);
1696 else if (lw_arc_is_pt(B1, B2, B3))
1697 return lw_dist2d_pt_arc(B1, A1, A2, A3, dl);
1698 else if (lw_arc_is_pt(A1, A2, A3))
1699 return lw_dist2d_pt_arc(A1, B1, B2, B3, dl);
1700
1701 /* Calculate centers and radii of circles. */
1702 radius_A = lw_arc_center(A1, A2, A3, &center_A);
1703 radius_B = lw_arc_center(B1, B2, B3, &center_B);
1704
1705 /* Two co-linear arcs?!? That's two segments. */
1706 if (radius_A < 0 && radius_B < 0)
1707 return lw_dist2d_seg_seg(A1, A3, B1, B3, dl);
1708
1709 /* A is co-linear, delegate to lw_dist_seg_arc here. */
1710 if (radius_A < 0)
1711 return lw_dist2d_seg_arc(A1, A3, B1, B2, B3, dl);
1712
1713 /* B is co-linear, delegate to lw_dist_seg_arc here. */
1714 if (radius_B < 0)
1715 return lw_dist2d_seg_arc(B1, B3, A1, A2, A3, dl);
1716
1717 /* Circle relationships */
1718 d = distance2d_pt_pt(&center_A, &center_B);
1719 is_disjoint = (d > (radius_A + radius_B));
1720 is_contained = (d < fabs(radius_A - radius_B));
1721 is_same_center = p2d_same(&center_A, &center_B);
1722 is_overlapping = ! (is_disjoint || is_contained || is_same_center);
1723
1724 /*
1725 * Prime the array of potential closest points with the
1726 * arc end points, which frequently participate in closest
1727 * points.
1728 */
1729 pts_A[ai++] = *A1;
1730 pts_A[ai++] = *A3;
1731 pts_B[bi++] = *B1;
1732 pts_B[bi++] = *B3;
1733
1734 /*
1735 * Overlapping circles might have a zero distance
1736 * case if the circle intersection points are inside both
1737 * arcs.
1738 */
1739 if (is_overlapping)
1740 {
1741 /*
1742 * Find the two points the circles intersect at.
1743 */
1744 uint32_t npoints = lw_dist2d_circle_circle_intersections(
1745 &center_A, radius_A,
1746 &center_B, radius_B,
1747 intersectionPts);
1748 for (uint32_t i = 0; i < npoints; i++)
1749 {
1750 /*
1751 * If an intersection point is contained in both
1752 * arcs, that is a location of zero distance, so
1753 * we are done calculating.
1754 */
1755 if (lw_pt_in_arc(&intersectionPts[i], A1, A2, A3) &&
1756 lw_pt_in_arc(&intersectionPts[i], B1, B2, B3))
1757 {
1758 lw_dist2d_distpts_set(dl, 0.0, &intersectionPts[i], &intersectionPts[i]);
1759 return LW_TRUE;
1760 }
1761 }
1762 }
1763
1764 /*
1765 * Join the circle centers and find the places that
1766 * line intersects the circles. Where those places
1767 * are in the arcs, they are potential sites of the
1768 * closest points.
1769 */
1770 if (is_disjoint || is_contained || is_overlapping)
1771 {
1772 if (!is_same_center)
1773 {
1774 /* Add points on A that intersect line from center_A to center_B */
1776 A1, A2, A3,
1777 &center_A, radius_A, &center_B,
1778 pts_A, &ai);
1779
1780 /* Add points on B that intersect line from center_B to center_A */
1782 B1, B2, B3,
1783 &center_B, radius_B, &center_A,
1784 pts_B, &bi);
1785 }
1786
1787 /* Add points on A that intersect line to B1 */
1789 A1, A2, A3,
1790 &center_A, radius_A, B1,
1791 pts_A, &ai);
1792
1793 /* Add points on A that intersect line to B3 */
1795 A1, A2, A3,
1796 &center_A, radius_A, B3,
1797 pts_A, &ai);
1798
1799 /* Add points on B that intersect line to A1 */
1801 B1, B2, B3,
1802 &center_B, radius_B, A1,
1803 pts_B, &bi);
1804
1805 /* Add points on B that intersect line to A3 */
1807 B1, B2, B3,
1808 &center_B, radius_B, A3,
1809 pts_B, &bi);
1810 }
1811
1812 /*
1813 * Now just brute force check all pairs of participating
1814 * points, to find the pair that is closest together.
1815 */
1816 for (uint32_t i = 0; i < ai; i++)
1817 for (uint32_t j = 0; j < bi; j++)
1818 lw_dist2d_pt_pt(&pts_A[i], &pts_B[j], dl);
1819
1820 return LW_TRUE;
1821}
1822
1823
1829int
1830lw_dist2d_seg_seg(const POINT2D *A, const POINT2D *B, const POINT2D *C, const POINT2D *D, DISTPTS *dl)
1831{
1832 double s_top, s_bot, s;
1833 double r_top, r_bot, r;
1834
1835 /*A and B are the same point */
1836 if ((A->x == B->x) && (A->y == B->y))
1837 {
1838 return lw_dist2d_pt_seg(A, C, D, dl);
1839 }
1840
1841 /*U and V are the same point */
1842 if ((C->x == D->x) && (C->y == D->y))
1843 {
1844 dl->twisted = ((dl->twisted) * (-1));
1845 return lw_dist2d_pt_seg(D, A, B, dl);
1846 }
1847
1848 /* AB and CD are line segments */
1849 /* from comp.graphics.algo
1850
1851 Solving the above for r and s yields
1852 (Ay-Cy)(Dx-Cx)-(Ax-Cx)(Dy-Cy)
1853 r = ----------------------------- (eqn 1)
1854 (Bx-Ax)(Dy-Cy)-(By-Ay)(Dx-Cx)
1855
1856 (Ay-Cy)(Bx-Ax)-(Ax-Cx)(By-Ay)
1857 s = ----------------------------- (eqn 2)
1858 (Bx-Ax)(Dy-Cy)-(By-Ay)(Dx-Cx)
1859 Let P be the position vector of the intersection point, then
1860 P=A+r(B-A) or
1861 Px=Ax+r(Bx-Ax)
1862 Py=Ay+r(By-Ay)
1863 By examining the values of r & s, you can also determine some other limiting conditions:
1864 If 0<=r<=1 & 0<=s<=1, intersection exists
1865 r<0 or r>1 or s<0 or s>1 line segments do not intersect
1866 If the denominator in eqn 1 is zero, AB & CD are parallel
1867 If the numerator in eqn 1 is also zero, AB & CD are collinear.
1868
1869 */
1870 r_top = (A->y - C->y) * (D->x - C->x) - (A->x - C->x) * (D->y - C->y);
1871 r_bot = (B->x - A->x) * (D->y - C->y) - (B->y - A->y) * (D->x - C->x);
1872
1873 s_top = (A->y - C->y) * (B->x - A->x) - (A->x - C->x) * (B->y - A->y);
1874 s_bot = (B->x - A->x) * (D->y - C->y) - (B->y - A->y) * (D->x - C->x);
1875
1876 if ((r_bot == 0) || (s_bot == 0))
1877 {
1878 if ((lw_dist2d_pt_seg(A, C, D, dl)) && (lw_dist2d_pt_seg(B, C, D, dl)))
1879 {
1880 /* change the order of inputted geometries and that we notice by changing sign on dl->twisted*/
1881 dl->twisted *= -1;
1882 return ((lw_dist2d_pt_seg(C, A, B, dl)) && (lw_dist2d_pt_seg(D, A, B, dl)));
1883 }
1884 else
1885 return LW_FALSE; /* if any of the calls to lw_dist2d_pt_seg goes wrong we return false*/
1886 }
1887
1888 s = s_top / s_bot;
1889 r = r_top / r_bot;
1890
1891 if (((r < 0) || (r > 1) || (s < 0) || (s > 1)) || (dl->mode == DIST_MAX))
1892 {
1893 if ((lw_dist2d_pt_seg(A, C, D, dl)) && (lw_dist2d_pt_seg(B, C, D, dl)))
1894 {
1895 /* change the order of inputted geometries and that we notice by changing sign on dl->twisted*/
1896 dl->twisted *= -1;
1897 return ((lw_dist2d_pt_seg(C, A, B, dl)) && (lw_dist2d_pt_seg(D, A, B, dl)));
1898 }
1899 else
1900 return LW_FALSE; /* if any of the calls to lw_dist2d_pt_seg goes wrong we return false*/
1901 }
1902 else
1903 {
1904 /* If there is intersection we identify the intersection point and return it but only if we are looking
1905 * for mindistance */
1906 if (dl->mode == DIST_MIN)
1907 {
1908 POINT2D theP;
1909
1910 if (((A->x == C->x) && (A->y == C->y)) || ((A->x == D->x) && (A->y == D->y)))
1911 {
1912 theP.x = A->x;
1913 theP.y = A->y;
1914 }
1915 else if (((B->x == C->x) && (B->y == C->y)) || ((B->x == D->x) && (B->y == D->y)))
1916 {
1917 theP.x = B->x;
1918 theP.y = B->y;
1919 }
1920 else
1921 {
1922 theP.x = A->x + r * (B->x - A->x);
1923 theP.y = A->y + r * (B->y - A->y);
1924 }
1925 lw_dist2d_distpts_set(dl, 0, &theP, &theP);
1926 }
1927 return LW_TRUE;
1928 }
1929}
1930
1931/*------------------------------------------------------------------------------------------------------------
1932End of Brute force functions
1933--------------------------------------------------------------------------------------------------------------*/
1934
1935/*------------------------------------------------------------------------------------------------------------
1936New faster distance calculations
1937--------------------------------------------------------------------------------------------------------------*/
1938
1946int
1948{
1949 /*here we define two lists to hold our calculated "z"-values and the order number in the geometry*/
1950
1951 double k, thevalue;
1952 float deltaX, deltaY, c1m, c2m;
1953 POINT2D c1, c2;
1954 const POINT2D *theP;
1955 float min1X, max1X, max1Y, min1Y, min2X, max2X, max2Y, min2Y;
1956 int t;
1957 int n1 = l1->npoints;
1958 int n2 = l2->npoints;
1959
1960 LISTSTRUCT *list1, *list2;
1961 list1 = (LISTSTRUCT *)lwalloc(sizeof(LISTSTRUCT) * n1);
1962 list2 = (LISTSTRUCT *)lwalloc(sizeof(LISTSTRUCT) * n2);
1963
1964 LWDEBUG(2, "lw_dist2d_fast_ptarray_ptarray is called");
1965
1966 max1X = box1->xmax;
1967 min1X = box1->xmin;
1968 max1Y = box1->ymax;
1969 min1Y = box1->ymin;
1970 max2X = box2->xmax;
1971 min2X = box2->xmin;
1972 max2Y = box2->ymax;
1973 min2Y = box2->ymin;
1974 /*we want the center of the bboxes, and calculate the slope between the centerpoints*/
1975 c1.x = min1X + (max1X - min1X) / 2;
1976 c1.y = min1Y + (max1Y - min1Y) / 2;
1977 c2.x = min2X + (max2X - min2X) / 2;
1978 c2.y = min2Y + (max2Y - min2Y) / 2;
1979
1980 deltaX = (c2.x - c1.x);
1981 deltaY = (c2.y - c1.y);
1982
1983 /*Here we calculate where the line perpendicular to the center-center line crosses the axes for each vertex
1984 if the center-center line is vertical the perpendicular line will be horizontal and we find it's crossing the
1985 Y-axes with z = y-kx */
1986 if ((deltaX * deltaX) < (deltaY * deltaY)) /*North or South*/
1987 {
1988 k = -deltaX / deltaY;
1989 for (t = 0; t < n1; t++) /*for each segment in L1 */
1990 {
1991 theP = getPoint2d_cp(l1, t);
1992 thevalue = theP->y - (k * theP->x);
1993 list1[t].themeasure = thevalue;
1994 list1[t].pnr = t;
1995 }
1996 for (t = 0; t < n2; t++) /*for each segment in L2*/
1997 {
1998 theP = getPoint2d_cp(l2, t);
1999 thevalue = theP->y - (k * theP->x);
2000 list2[t].themeasure = thevalue;
2001 list2[t].pnr = t;
2002 }
2003 c1m = c1.y - (k * c1.x);
2004 c2m = c2.y - (k * c2.x);
2005 }
2006
2007 /*if the center-center line is horizontal the perpendicular line will be vertical. To eliminate problems with
2008 dividing by zero we are here mirroring the coordinate-system and we find it's crossing the X-axes with z =
2009 x-(1/k)y */
2010 else /*West or East*/
2011 {
2012 k = -deltaY / deltaX;
2013 for (t = 0; t < n1; t++) /*for each segment in L1 */
2014 {
2015 theP = getPoint2d_cp(l1, t);
2016 thevalue = theP->x - (k * theP->y);
2017 list1[t].themeasure = thevalue;
2018 list1[t].pnr = t;
2019 /* lwnotice("l1 %d, measure=%f",t,thevalue ); */
2020 }
2021 for (t = 0; t < n2; t++) /*for each segment in L2*/
2022 {
2023 theP = getPoint2d_cp(l2, t);
2024 thevalue = theP->x - (k * theP->y);
2025 list2[t].themeasure = thevalue;
2026 list2[t].pnr = t;
2027 /* lwnotice("l2 %d, measure=%f",t,thevalue ); */
2028 }
2029 c1m = c1.x - (k * c1.y);
2030 c2m = c2.x - (k * c2.y);
2031 }
2032
2033 /*we sort our lists by the calculated values*/
2034 qsort(list1, n1, sizeof(LISTSTRUCT), struct_cmp_by_measure);
2035 qsort(list2, n2, sizeof(LISTSTRUCT), struct_cmp_by_measure);
2036
2037 if (c1m < c2m)
2038 {
2039 if (!lw_dist2d_pre_seg_seg(l1, l2, list1, list2, k, dl))
2040 {
2041 lwfree(list1);
2042 lwfree(list2);
2043 return LW_FALSE;
2044 }
2045 }
2046 else
2047 {
2048 dl->twisted = ((dl->twisted) * (-1));
2049 if (!lw_dist2d_pre_seg_seg(l2, l1, list2, list1, k, dl))
2050 {
2051 lwfree(list1);
2052 lwfree(list2);
2053 return LW_FALSE;
2054 }
2055 }
2056 lwfree(list1);
2057 lwfree(list2);
2058 return LW_TRUE;
2059}
2060
2061int
2062struct_cmp_by_measure(const void *a, const void *b)
2063{
2064 LISTSTRUCT *ia = (LISTSTRUCT *)a;
2065 LISTSTRUCT *ib = (LISTSTRUCT *)b;
2066 return (ia->themeasure > ib->themeasure) ? 1 : ((ia->themeasure < ib->themeasure) ? -1 : 0);
2067}
2068
2070int
2072{
2073 const POINT2D *p1, *p2, *p3, *p4, *p01, *p02;
2074 int pnr1, pnr2, pnr3, pnr4, n1, n2, i, u, r, twist;
2075 double maxmeasure;
2076 n1 = l1->npoints;
2077 n2 = l2->npoints;
2078
2079 LWDEBUG(2, "lw_dist2d_pre_seg_seg is called");
2080
2081 p1 = getPoint2d_cp(l1, list1[0].pnr);
2082 p3 = getPoint2d_cp(l2, list2[0].pnr);
2083 lw_dist2d_pt_pt(p1, p3, dl);
2084 maxmeasure = sqrt(dl->distance * dl->distance + (dl->distance * dl->distance * k * k));
2085 twist = dl->twisted; /*to keep the incoming order between iterations*/
2086 for (i = (n1 - 1); i >= 0; --i)
2087 {
2088 /*we break this iteration when we have checked every point closer to our perpendicular "checkline" than
2089 * our shortest found distance*/
2090 if (((list2[0].themeasure - list1[i].themeasure)) > maxmeasure)
2091 break;
2092 /*because we are not iterating in the original point order we have to check the segment before and after
2093 * every point*/
2094 for (r = -1; r <= 1; r += 2)
2095 {
2096 pnr1 = list1[i].pnr;
2097 p1 = getPoint2d_cp(l1, pnr1);
2098 if (pnr1 + r < 0)
2099 {
2100 p01 = getPoint2d_cp(l1, (n1 - 1));
2101 if ((p1->x == p01->x) && (p1->y == p01->y))
2102 pnr2 = (n1 - 1);
2103 else
2104 pnr2 = pnr1; /* if it is a line and the last and first point is not the same we
2105 avoid the edge between start and end this way*/
2106 }
2107
2108 else if (pnr1 + r > (n1 - 1))
2109 {
2110 p01 = getPoint2d_cp(l1, 0);
2111 if ((p1->x == p01->x) && (p1->y == p01->y))
2112 pnr2 = 0;
2113 else
2114 pnr2 = pnr1; /* if it is a line and the last and first point is not the same we
2115 avoid the edge between start and end this way*/
2116 }
2117 else
2118 pnr2 = pnr1 + r;
2119
2120 p2 = getPoint2d_cp(l1, pnr2);
2121 for (u = 0; u < n2; ++u)
2122 {
2123 if (((list2[u].themeasure - list1[i].themeasure)) >= maxmeasure)
2124 break;
2125 pnr3 = list2[u].pnr;
2126 p3 = getPoint2d_cp(l2, pnr3);
2127 if (pnr3 == 0)
2128 {
2129 p02 = getPoint2d_cp(l2, (n2 - 1));
2130 if ((p3->x == p02->x) && (p3->y == p02->y))
2131 pnr4 = (n2 - 1);
2132 else
2133 pnr4 = pnr3; /* if it is a line and the last and first point is not the
2134 same we avoid the edge between start and end this way*/
2135 }
2136 else
2137 pnr4 = pnr3 - 1;
2138
2139 p4 = getPoint2d_cp(l2, pnr4);
2140 dl->twisted = twist;
2141 if (!lw_dist2d_selected_seg_seg(p1, p2, p3, p4, dl))
2142 return LW_FALSE;
2143
2144 if (pnr3 >= (n2 - 1))
2145 {
2146 p02 = getPoint2d_cp(l2, 0);
2147 if ((p3->x == p02->x) && (p3->y == p02->y))
2148 pnr4 = 0;
2149 else
2150 pnr4 = pnr3; /* if it is a line and the last and first point is not the
2151 same we avoid the edge between start and end this way*/
2152 }
2153
2154 else
2155 pnr4 = pnr3 + 1;
2156
2157 p4 = getPoint2d_cp(l2, pnr4);
2158 dl->twisted = twist; /*we reset the "twist" for each iteration*/
2159 if (!lw_dist2d_selected_seg_seg(p1, p2, p3, p4, dl))
2160 return LW_FALSE;
2161 /*here we "translate" the found mindistance so it can be compared to our "z"-values*/
2162 maxmeasure = sqrt(dl->distance * dl->distance + (dl->distance * dl->distance * k * k));
2163 }
2164 }
2165 }
2166
2167 return LW_TRUE;
2168}
2169
2175int
2176lw_dist2d_selected_seg_seg(const POINT2D *A, const POINT2D *B, const POINT2D *C, const POINT2D *D, DISTPTS *dl)
2177{
2178 /*A and B are the same point */
2179 if ((A->x == B->x) && (A->y == B->y))
2180 {
2181 return lw_dist2d_pt_seg(A, C, D, dl);
2182 }
2183 /*U and V are the same point */
2184
2185 if ((C->x == D->x) && (C->y == D->y))
2186 {
2187 dl->twisted *= -1;
2188 return lw_dist2d_pt_seg(D, A, B, dl);
2189 }
2190
2191 if ((lw_dist2d_pt_seg(A, C, D, dl)) && (lw_dist2d_pt_seg(B, C, D, dl)))
2192 {
2193 /* change the order of inputted geometries and that we notice by changing sign on dl->twisted */
2194 dl->twisted *= -1;
2195 return ((lw_dist2d_pt_seg(C, A, B, dl)) && (lw_dist2d_pt_seg(D, A, B, dl)));
2196 }
2197 else
2198 return LW_FALSE; /* if any of the calls to lw_dist2d_pt_seg goes wrong we return false*/
2199}
2200
2201/*------------------------------------------------------------------------------------------------------------
2202End of New faster distance calculations
2203--------------------------------------------------------------------------------------------------------------*/
2204
2205/*------------------------------------------------------------------------------------------------------------
2206Functions in common for Brute force and new calculation
2207--------------------------------------------------------------------------------------------------------------*/
2208
2216int
2217lw_dist2d_pt_seg(const POINT2D *p, const POINT2D *A, const POINT2D *B, DISTPTS *dl)
2218{
2219 double r;
2220
2221 LWDEBUG(2, "lw_dist2d_pt_seg called");
2222
2223 /*if start==end, then use pt distance */
2224 if ((A->x == B->x) && (A->y == B->y))
2225 {
2226 LWDEBUG(2, "lw_dist2d_pt_seg found first and last segment points being the same");
2227 return lw_dist2d_pt_pt(p, A, dl);
2228 }
2229
2230
2231 /*
2232 * otherwise, we use comp.graphics.algorithms
2233 * Frequently Asked Questions method
2234 *
2235 * (1) AC dot AB
2236 * r = ---------
2237 * ||AB||^2
2238 * r has the following meaning:
2239 * r=0 P = A
2240 * r=1 P = B
2241 * r<0 P is on the backward extension of AB
2242 * r>1 P is on the forward extension of AB
2243 * 0<r<1 P is interior to AB
2244 */
2245
2246 r = ((p->x - A->x) * (B->x - A->x) + (p->y - A->y) * (B->y - A->y)) /
2247 ((B->x - A->x) * (B->x - A->x) + (B->y - A->y) * (B->y - A->y));
2248
2249 LWDEBUGF(2, "lw_dist2d_pt_seg found r = %.15g", r);
2250
2251 /*This is for finding the maxdistance.
2252 the maxdistance have to be between two vertices, compared to mindistance which can be between two vertices.*/
2253 if (dl->mode == DIST_MAX)
2254 {
2255 if (r >= 0.5)
2256 return lw_dist2d_pt_pt(p, A, dl);
2257 else /* (r < 0.5) */
2258 return lw_dist2d_pt_pt(p, B, dl);
2259 }
2260
2261 if (r < 0) /*If p projected on the line is outside point A*/
2262 return lw_dist2d_pt_pt(p, A, dl);
2263 if (r >= 1) /*If p projected on the line is outside point B or on point B*/
2264 return lw_dist2d_pt_pt(p, B, dl);
2265
2266 /*If the point p is on the segment this is a more robust way to find out that*/
2267 if ((((A->y - p->y) * (B->x - A->x) == (A->x - p->x) * (B->y - A->y))) && (dl->mode == DIST_MIN))
2268 {
2269 lw_dist2d_distpts_set(dl, 0, p, p);
2270 }
2271
2272
2273 /*
2274 (2)
2275 (Ay-Cy)(Bx-Ax)-(Ax-Cx)(By-Ay)
2276 s = -----------------------------
2277 L^2
2278
2279 Then the distance from C to P = |s|*L.
2280 */
2281
2282 double s = ((A->y - p->y) * (B->x - A->x) - (A->x - p->x) * (B->y - A->y)) /
2283 ((B->x - A->x) * (B->x - A->x) + (B->y - A->y) * (B->y - A->y));
2284
2285 double dist = fabs(s) * sqrt(((B->x - A->x) * (B->x - A->x) + (B->y - A->y) * (B->y - A->y)));
2286 if ( dist < dl->distance )
2287 {
2288 dl->distance = dist;
2289 {
2290 POINT2D c;
2291 c.x = A->x + r * (B->x - A->x);
2292 c.y = A->y + r * (B->y - A->y);
2293 if (dl->twisted > 0)
2294 {
2295 dl->p1 = *p;
2296 dl->p2 = c;
2297 }
2298 else
2299 {
2300 dl->p1 = c;
2301 dl->p2 = *p;
2302 }
2303 }
2304 }
2305
2306 return LW_TRUE;
2307}
2308
2311int
2312lw_dist2d_pt_pt(const POINT2D *thep1, const POINT2D *thep2, DISTPTS *dl)
2313{
2314 double hside = thep2->x - thep1->x;
2315 double vside = thep2->y - thep1->y;
2316 double dist = sqrt(hside * hside + vside * vside);
2317
2318 /*multiplication with mode to handle mindistance (mode=1) and maxdistance (mode = (-1)*/
2319 if (((dl->distance - dist) * (dl->mode)) > 0)
2320 {
2321 dl->distance = dist;
2322
2323 /* To get the points in right order. twisted is updated between 1 and (-1) every time the order is
2324 * changed earlier in the chain*/
2325 if (dl->twisted > 0)
2326 {
2327 dl->p1 = *thep1;
2328 dl->p2 = *thep2;
2329 }
2330 else
2331 {
2332 dl->p1 = *thep2;
2333 dl->p2 = *thep1;
2334 }
2335 }
2336 return LW_TRUE;
2337}
2338
2339/*------------------------------------------------------------------------------------------------------------
2340End of Functions in common for Brute force and new calculation
2341--------------------------------------------------------------------------------------------------------------*/
2342
2343inline double
2344distance2d_pt_pt(const POINT2D *p1, const POINT2D *p2)
2345{
2346 double hside = p2->x - p1->x;
2347 double vside = p2->y - p1->y;
2348
2349 return hypot(hside, vside);
2350}
2351
2352/* return distance squared, useful to avoid sqrt calculations */
2353double
2354distance2d_sqr_pt_seg(const POINT2D *C, const POINT2D *A, const POINT2D *B)
2355{
2356 /*if start==end, then use pt distance */
2357 if ((A->x == B->x) && (A->y == B->y))
2358 return distance2d_sqr_pt_pt(C, A);
2359
2360 /*
2361 * otherwise, we use comp.graphics.algorithms
2362 * Frequently Asked Questions method
2363 *
2364 * (1) AC dot AB
2365 * r = ---------
2366 * ||AB||^2
2367 * r has the following meaning:
2368 * r=0 P = A
2369 * r=1 P = B
2370 * r<0 P is on the backward extension of AB
2371 * r>1 P is on the forward extension of AB
2372 * 0<r<1 P is interior to AB
2373 */
2374
2375 double ba_x = (B->x - A->x);
2376 double ba_y = (B->y - A->y);
2377 double ab_length_sqr = (ba_x * ba_x + ba_y * ba_y);
2378 double ca_x = (C->x - A->x);
2379 double ca_y = (C->y - A->y);
2380 double dot_ac_ab = (ca_x * ba_x + ca_y * ba_y);
2381
2382 if (dot_ac_ab <= 0)
2383 return distance2d_sqr_pt_pt(C, A);
2384 if (dot_ac_ab >= ab_length_sqr)
2385 return distance2d_sqr_pt_pt(C, B);
2386
2387 /*
2388 * (2)
2389 * (Ay-Cy)(Bx-Ax)-(Ax-Cx)(By-Ay)
2390 * s = -----------------------------
2391 * L^2
2392 *
2393 * Then the distance from C to P = |s|*L.
2394 *
2395 */
2396
2397 double s_numerator = ca_x * ba_y - ca_y * ba_x;
2398
2399 /* Distance = (s_num / ab) * (s_num / ab) * ab == s_num * s_num / ab) */
2400 return s_numerator * s_numerator / ab_length_sqr;
2401}
2402
2407int
2408azimuth_pt_pt(const POINT2D *A, const POINT2D *B, double *d)
2409{
2410 if (A->x == B->x && A->y == B->y)
2411 return LW_FALSE;
2412 *d = fmod(2 * M_PI + M_PI / 2 - atan2(B->y - A->y, B->x - A->x), 2 * M_PI);
2413 return LW_TRUE;
2414}
2415
2420int
2421project_pt(const POINT2D *P, double distance, double azimuth, POINT2D *R)
2422{
2423 const double TWOPI = 2.0 * M_PI;
2424 double slope;
2425 /* Deal with azimuth out of (-360,360) range */
2426 int orbits = floor(azimuth / TWOPI);
2427 azimuth -= TWOPI * orbits;
2428 /* Convert from azimuth to conventional slope */
2429 slope = TWOPI - azimuth + M_PI_2;
2430 if (slope > 0 && slope > TWOPI) slope -= TWOPI;
2431 if (slope < 0 && slope < -TWOPI) slope += TWOPI;
2432
2433 double dx = cos(slope) * distance;
2434 double dy = sin(slope) * distance;
2435 R->x = P->x + dx;
2436 R->y = P->y + dy;
2437 return LW_TRUE;
2438}
2439
2444int
2445project_pt_pt(const POINT4D *A, const POINT4D *B, double distance, POINT4D *R)
2446{
2447 /* Convert from azimuth to conventional slope */
2448 double len = distance2d_pt_pt((const POINT2D *)A, (const POINT2D *)B);
2449 double prop = distance / len;
2450 double dx = (B->x - A->x) * prop;
2451 double dy = (B->y - A->y) * prop;
2452 double dz = (B->z - A->z) * prop;
2453 double dm = (B->m - A->m) * prop;
2454 R->x = B->x + dx;
2455 R->y = B->y + dy;
2456 if (isfinite(dz)) R->z = B->z + dz;
2457 if (isfinite(dm)) R->m = B->m + dm;
2458 return LW_TRUE;
2459}
2460
char * s
Definition cu_in_wkt.c:23
char * r
Definition cu_in_wkt.c:24
char result[OUT_DOUBLE_BUFFER_SIZE]
Definition cu_print.c:267
const char * lwtype_name(uint8_t type)
Return the type name string associated with a type number (e.g.
Definition lwutil.c:216
#define LW_FALSE
Definition liblwgeom.h:94
#define COLLECTIONTYPE
Definition liblwgeom.h:108
#define COMPOUNDTYPE
Definition liblwgeom.h:110
void lwgeom_free(LWGEOM *geom)
Definition lwgeom.c:1246
#define CURVEPOLYTYPE
Definition liblwgeom.h:111
#define MULTILINETYPE
Definition liblwgeom.h:106
#define MULTISURFACETYPE
Definition liblwgeom.h:113
#define LINETYPE
Definition liblwgeom.h:103
LWCOLLECTION * lwgeom_as_lwcollection(const LWGEOM *lwgeom)
Definition lwgeom.c:261
LWCURVEPOLY * lwcurvepoly_construct_from_lwpoly(LWPOLY *lwpoly)
Construct an equivalent curve polygon from a polygon.
Definition lwcurvepoly.c:52
#define MULTIPOINTTYPE
Definition liblwgeom.h:105
#define POINTTYPE
LWTYPE numbers, used internally by PostGIS.
Definition liblwgeom.h:102
void * lwalloc(size_t size)
Definition lwutil.c:227
#define TINTYPE
Definition liblwgeom.h:116
#define MULTIPOLYGONTYPE
Definition liblwgeom.h:107
void lwfree(void *mem)
Definition lwutil.c:248
LWPOINT * lwpoint_make2d(int32_t srid, double x, double y)
Definition lwpoint.c:163
#define POLYGONTYPE
Definition liblwgeom.h:104
#define POLYHEDRALSURFACETYPE
Definition liblwgeom.h:114
#define CIRCSTRINGTYPE
Definition liblwgeom.h:109
LWCOLLECTION * lwcollection_construct_empty(uint8_t type, int32_t srid, char hasz, char hasm)
#define MULTICURVETYPE
Definition liblwgeom.h:112
#define TRIANGLETYPE
Definition liblwgeom.h:115
#define LW_TRUE
Return types for functions with status returns.
Definition liblwgeom.h:93
void lwgeom_add_bbox(LWGEOM *lwgeom)
Compute a bbox if not already computed.
Definition lwgeom.c:723
LWLINE * lwline_from_ptarray(int32_t srid, uint32_t npoints, LWPOINT **points)
Definition lwline.c:238
#define LW_INSIDE
Constants for point-in-polygon return values.
double lw_arc_center(const POINT2D *p1, const POINT2D *p2, const POINT2D *p3, POINT2D *result)
Determines the center of the circle defined by the three given points.
#define FP_EQUALS(A, B)
int ptarray_contains_point(const POINTARRAY *pa, const POINT2D *pt)
The following is based on the "Fast Winding Number Inclusion of a Point in a Polygon" algorithm by Da...
Definition ptarray.c:755
int lw_pt_in_seg(const POINT2D *P, const POINT2D *A1, const POINT2D *A2)
Returns true if P is between A1/A2.
#define LW_OUTSIDE
int lwgeom_contains_point(const LWGEOM *geom, const POINT2D *pt)
Definition lwcompound.c:172
int lw_arc_is_pt(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3)
Returns true if arc A is actually a point (all vertices are the same) .
int lw_pt_in_arc(const POINT2D *P, const POINT2D *A1, const POINT2D *A2, const POINT2D *A3)
Returns true if P is on the same side of the plane partition defined by A1/A3 as A2 is.
Definition lwalgorithm.c:91
#define FP_IS_ZERO(A)
int p2d_same(const POINT2D *p1, const POINT2D *p2)
Definition lwalgorithm.c:57
#define LWDEBUG(level, msg)
Definition lwgeom_log.h:101
#define LWDEBUGF(level, msg,...)
Definition lwgeom_log.h:106
void void lwerror(const char *fmt,...) __attribute__((format(printf
Write a notice out to the error handler.
static double distance2d_sqr_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition lwinline.h:33
static int lwgeom_is_empty(const LWGEOM *geom)
Return true or false depending on whether a geometry is an "empty" geometry (no vertices members)
Definition lwinline.h:199
static const POINT2D * getPoint2d_cp(const POINTARRAY *pa, uint32_t n)
Returns a POINT2D pointer into the POINTARRAY serialized_ptlist, suitable for reading from.
Definition lwinline.h:97
static double distance(double x1, double y1, double x2, double y2)
Definition lwtree.c:1032
static uint32_t lw_dist2d_circle_circle_intersections(const POINT2D *cA, double rA, const POINT2D *cB, double rB, POINT2D *I)
Calculates the intersection points of two overlapping circles.
Definition measures.c:1626
int lw_dist2d_circstring_curvepoly(LWCIRCSTRING *circ, LWCURVEPOLY *poly, DISTPTS *dl)
Definition measures.c:1042
int struct_cmp_by_measure(const void *a, const void *b)
Definition measures.c:2062
int lw_dist2d_line_curvepoly(LWLINE *line, LWCURVEPOLY *poly, DISTPTS *dl)
Definition measures.c:781
double distance2d_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition measures.c:2344
double lwgeom_maxdistance2d(const LWGEOM *lw1, const LWGEOM *lw2)
Function initializing max distance calculation.
Definition measures.c:180
LWGEOM * lwgeom_closest_point(const LWGEOM *lw1, const LWGEOM *lw2)
Definition measures.c:55
int lw_dist2d_line_circstring(LWLINE *line1, LWCIRCSTRING *line2, DISTPTS *dl)
Definition measures.c:731
int lw_dist2d_poly_poly(LWPOLY *poly1, LWPOLY *poly2, DISTPTS *dl)
Function handling polygon to polygon calculation 1 if we are looking for maxdistance,...
Definition measures.c:968
int lw_dist2d_comp(const LWGEOM *lw1, const LWGEOM *lw2, DISTPTS *dl)
This function just deserializes geometries Bboxes is not checked here since it is the subgeometries b...
Definition measures.c:251
LWGEOM * lwgeom_closest_line(const LWGEOM *lw1, const LWGEOM *lw2)
Definition measures.c:43
int lw_dist2d_pt_arc(const POINT2D *P, const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, DISTPTS *dl)
Definition measures.c:1495
static void lw_dist2d_distpts_set(DISTPTS *dl, double distance, const POINT2D *p1, const POINT2D *p2)
Definition measures.c:81
int lw_dist2d_pre_seg_seg(POINTARRAY *l1, POINTARRAY *l2, LISTSTRUCT *list1, LISTSTRUCT *list2, double k, DISTPTS *dl)
preparation before lw_dist2d_seg_seg.
Definition measures.c:2071
int azimuth_pt_pt(const POINT2D *A, const POINT2D *B, double *d)
Compute the azimuth of segment AB in radians.
Definition measures.c:2408
int lw_dist2d_point_curvepoly(LWPOINT *point, LWCURVEPOLY *poly, DISTPTS *dl)
Definition measures.c:684
int lw_dist2d_point_circstring(LWPOINT *point, LWCIRCSTRING *circ, DISTPTS *dl)
Definition measures.c:643
int lw_dist2d_pt_seg(const POINT2D *p, const POINT2D *A, const POINT2D *B, DISTPTS *dl)
lw_dist2d_comp from p to line A->B This one is now sending every occasion to lw_dist2d_pt_pt Before i...
Definition measures.c:2217
int lw_dist2d_circstring_poly(LWCIRCSTRING *circ, LWPOLY *poly, DISTPTS *dl)
Definition measures.c:1033
int lw_dist2d_tri_circstring(LWTRIANGLE *tri, LWCIRCSTRING *line, DISTPTS *dl)
Definition measures.c:946
int lw_dist2d_line_tri(LWLINE *line, LWTRIANGLE *tri, DISTPTS *dl)
Definition measures.c:717
static int lw_dist2d_check_overlap(const LWGEOM *lwg1, const LWGEOM *lwg2)
Definition measures.c:282
int lw_dist2d_pt_ptarray(const POINT2D *p, POINTARRAY *pa, DISTPTS *dl)
search all the segments of pointarray to see which one is closest to p1 Returns minimum distance betw...
Definition measures.c:1112
int lw_dist2d_recursive(const LWGEOM *lwg1, const LWGEOM *lwg2, DISTPTS *dl)
This is a recursive function delivering every possible combination of subgeometries.
Definition measures.c:299
LWGEOM * lw_dist2d_distanceline(const LWGEOM *lw1, const LWGEOM *lw2, int32_t srid, int mode)
Function initializing shortestline and longestline calculations.
Definition measures.c:96
int lw_dist2d_seg_seg(const POINT2D *A, const POINT2D *B, const POINT2D *C, const POINT2D *D, DISTPTS *dl)
Finds the shortest distance between two segments.
Definition measures.c:1830
int lw_dist2d_curvepoly_curvepoly(LWCURVEPOLY *poly1, LWCURVEPOLY *poly2, DISTPTS *dl)
Definition measures.c:1054
int lw_dist2d_poly_curvepoly(LWPOLY *poly1, LWCURVEPOLY *curvepoly2, DISTPTS *dl)
Definition measures.c:1024
int lw_dist2d_seg_arc(const POINT2D *A1, const POINT2D *A2, const POINT2D *B1, const POINT2D *B2, const POINT2D *B3, DISTPTS *dl)
Calculate the shortest distance between an arc and an edge.
Definition measures.c:1351
int lw_dist2d_point_line(LWPOINT *point, LWLINE *line, DISTPTS *dl)
point to line calculation
Definition measures.c:622
int lw_dist2d_tri_curvepoly(LWTRIANGLE *tri, LWCURVEPOLY *poly, DISTPTS *dl)
Definition measures.c:904
int project_pt(const POINT2D *P, double distance, double azimuth, POINT2D *R)
Azimuth is angle in radians from vertical axis.
Definition measures.c:2421
double lwgeom_mindistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
Function handling min distance calculations and dwithin calculations.
Definition measures.c:222
void lw_dist2d_distpts_init(DISTPTS *dl, int mode)
Definition measures.c:67
int lw_dist2d_fast_ptarray_ptarray(POINTARRAY *l1, POINTARRAY *l2, DISTPTS *dl, GBOX *box1, GBOX *box2)
The new faster calculation comparing pointarray to another pointarray the arrays can come from both p...
Definition measures.c:1947
static int lw_dist2d_is_collection(const LWGEOM *g)
Definition measures.c:257
int lw_dist2d_arc_arc(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, const POINT2D *B1, const POINT2D *B2, const POINT2D *B3, DISTPTS *dl)
Definition measures.c:1677
int lw_dist2d_line_line(LWLINE *line1, LWLINE *line2, DISTPTS *dl)
Definition measures.c:709
int lw_dist2d_tri_tri(LWTRIANGLE *tri1, LWTRIANGLE *tri2, DISTPTS *dl)
Definition measures.c:813
int lw_dist2d_point_tri(LWPOINT *point, LWTRIANGLE *tri, DISTPTS *dl)
Definition measures.c:629
double lwgeom_mindistance2d(const LWGEOM *lw1, const LWGEOM *lw2)
Function initializing min distance calculation.
Definition measures.c:212
int lw_dist2d_ptarray_ptarray(POINTARRAY *l1, POINTARRAY *l2, DISTPTS *dl)
test each segment of l1 against each segment of l2.
Definition measures.c:1197
int lw_dist2d_ptarray_ptarrayarc(const POINTARRAY *pa, const POINTARRAY *pb, DISTPTS *dl)
Test each segment of pa against each arc of pb for distance.
Definition measures.c:1246
int lw_dist2d_selected_seg_seg(const POINT2D *A, const POINT2D *B, const POINT2D *C, const POINT2D *D, DISTPTS *dl)
This is the same function as lw_dist2d_seg_seg but without any calculations to determine intersection...
Definition measures.c:2176
int lw_dist2d_line_poly(LWLINE *line, LWPOLY *poly, DISTPTS *dl)
line to polygon calculation Brute force.
Definition measures.c:748
int lw_dist2d_distribute_bruteforce(const LWGEOM *lwg1, const LWGEOM *lwg2, DISTPTS *dl)
Definition measures.c:393
static const POINT2D * lw_curvering_getfirstpoint2d_cp(LWGEOM *geom)
Definition measures.c:883
int lw_dist2d_ptarrayarc_ptarrayarc(const POINTARRAY *pa, const POINTARRAY *pb, DISTPTS *dl)
Test each arc of pa against each arc of pb for distance.
Definition measures.c:1300
double lwgeom_maxdistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
Function handling max distance calculations and dfullywithin calculations.
Definition measures.c:192
int project_pt_pt(const POINT4D *A, const POINT4D *B, double distance, POINT4D *R)
Azimuth is angle in radians from vertical axis.
Definition measures.c:2445
LWGEOM * lw_dist2d_distancepoint(const LWGEOM *lw1, const LWGEOM *lw2, int32_t srid, int mode)
Function initializing closestpoint calculations.
Definition measures.c:143
LWGEOM * lwgeom_furthest_line(const LWGEOM *lw1, const LWGEOM *lw2)
Definition measures.c:49
double distance2d_sqr_pt_seg(const POINT2D *C, const POINT2D *A, const POINT2D *B)
Definition measures.c:2354
int lw_dist2d_distribute_fast(LWGEOM *lwg1, LWGEOM *lwg2, DISTPTS *dl)
Geometries are distributed for the new faster distance-calculations.
Definition measures.c:555
int lw_dist2d_pt_pt(const POINT2D *thep1, const POINT2D *thep2, DISTPTS *dl)
Compares incoming points and stores the points closest to each other or most far away from each other...
Definition measures.c:2312
int lw_dist2d_circstring_circstring(LWCIRCSTRING *line1, LWCIRCSTRING *line2, DISTPTS *dl)
Definition measures.c:1048
int lw_dist2d_pt_ptarrayarc(const POINT2D *p, const POINTARRAY *pa, DISTPTS *dl)
Search all the arcs of pointarray to see which one is closest to p1 Returns minimum distance between ...
Definition measures.c:1148
int lw_dist2d_point_poly(LWPOINT *point, LWPOLY *poly, DISTPTS *dl)
Definition measures.c:656
LWGEOM * lwgeom_furthest_point(const LWGEOM *lw1, const LWGEOM *lw2)
Definition measures.c:61
static int lw_dist2d_circle_intersections(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, const POINT2D *center_A, double radius_A, const POINT2D *P, POINT2D *I, uint32_t *ni)
Calculates the intersection points of a circle and line.
Definition measures.c:1566
int lw_dist2d_tri_poly(LWTRIANGLE *tri, LWPOLY *poly, DISTPTS *dl)
Definition measures.c:835
int lw_dist2d_point_point(LWPOINT *point1, LWPOINT *point2, DISTPTS *dl)
point to point calculation
Definition measures.c:612
#define DIST_MIN
Definition measures.h:44
#define DIST_MAX
Definition measures.h:43
POINT2D p1
Definition measures.h:52
POINT2D p2
Definition measures.h:53
int twisted
Definition measures.h:55
double tolerance
Definition measures.h:56
int mode
Definition measures.h:54
double distance
Definition measures.h:51
Structure used in distance-calculations.
Definition measures.h:50
double ymax
Definition liblwgeom.h:357
double xmax
Definition liblwgeom.h:355
double ymin
Definition liblwgeom.h:356
double xmin
Definition liblwgeom.h:354
double themeasure
Definition measures.h:61
POINTARRAY * points
Definition liblwgeom.h:507
uint32_t ngeoms
Definition liblwgeom.h:580
LWGEOM ** geoms
Definition liblwgeom.h:575
LWGEOM ** geoms
Definition liblwgeom.h:589
LWGEOM ** rings
Definition liblwgeom.h:603
uint32_t nrings
Definition liblwgeom.h:608
uint8_t type
Definition liblwgeom.h:462
GBOX * bbox
Definition liblwgeom.h:458
int32_t srid
Definition liblwgeom.h:460
POINTARRAY * points
Definition liblwgeom.h:483
POINTARRAY * point
Definition liblwgeom.h:471
POINTARRAY ** rings
Definition liblwgeom.h:519
uint32_t nrings
Definition liblwgeom.h:524
POINTARRAY * points
Definition liblwgeom.h:495
double y
Definition liblwgeom.h:390
double x
Definition liblwgeom.h:390
double m
Definition liblwgeom.h:414
double x
Definition liblwgeom.h:414
double z
Definition liblwgeom.h:414
double y
Definition liblwgeom.h:414
uint32_t npoints
Definition liblwgeom.h:427