PostGIS 3.6.2dev-r@@SVN_REVISION@@
Loading...
Searching...
No Matches
ptarray.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 (C) 2012-2021 Sandro Santilli <strk@kbt.io>
22 * Copyright (C) 2001-2006 Refractions Research Inc.
23 *
24 **********************************************************************/
25
26
27#include <stdio.h>
28#include <string.h>
29#include <stdlib.h>
30
31#include "../postgis_config.h"
32/*#define POSTGIS_DEBUG_LEVEL 4*/
33#include "liblwgeom_internal.h"
34#include "lwgeom_log.h"
35
36int
38{
39 if ( ! pa ) return LW_FALSE;
40 return FLAGS_GET_Z(pa->flags);
41}
42
43int
45{
46 if ( ! pa ) return LW_FALSE;
47 return FLAGS_GET_M(pa->flags);
48}
49
51ptarray_construct(char hasz, char hasm, uint32_t npoints)
52{
53 POINTARRAY *pa = ptarray_construct_empty(hasz, hasm, npoints);
54 pa->npoints = npoints;
55 return pa;
56}
57
59ptarray_construct_empty(char hasz, char hasm, uint32_t maxpoints)
60{
61 POINTARRAY *pa = lwalloc(sizeof(POINTARRAY));
62 pa->serialized_pointlist = NULL;
63
64 /* Set our dimensionality info on the bitmap */
65 pa->flags = lwflags(hasz, hasm, 0);
66
67 /* We will be allocating a bit of room */
68 pa->npoints = 0;
69 pa->maxpoints = maxpoints;
70
71 /* Allocate the coordinate array */
72 if ( maxpoints > 0 )
73 pa->serialized_pointlist = lwalloc(maxpoints * ptarray_point_size(pa));
74 else
75 pa->serialized_pointlist = NULL;
76
77 return pa;
78}
79
80/*
81* Add a point into a pointarray. Only adds as many dimensions as the
82* pointarray supports.
83*/
84int
85ptarray_insert_point(POINTARRAY *pa, const POINT4D *p, uint32_t where)
86{
87 if (!pa || !p)
88 return LW_FAILURE;
89 size_t point_size = ptarray_point_size(pa);
90 LWDEBUGF(5,"pa = %p; p = %p; where = %d", pa, p, where);
91 LWDEBUGF(5,"pa->npoints = %d; pa->maxpoints = %d", pa->npoints, pa->maxpoints);
92
93 if ( FLAGS_GET_READONLY(pa->flags) )
94 {
95 lwerror("ptarray_insert_point: called on read-only point array");
96 return LW_FAILURE;
97 }
98
99 /* Error on invalid offset value */
100 if ( where > pa->npoints )
101 {
102 lwerror("ptarray_insert_point: offset out of range (%d)", where);
103 return LW_FAILURE;
104 }
105
106 /* If we have no storage, let's allocate some */
107 if( pa->maxpoints == 0 || ! pa->serialized_pointlist )
108 {
109 pa->maxpoints = 32;
110 pa->npoints = 0;
112 }
113
114 /* Error out if we have a bad situation */
115 if ( pa->npoints > pa->maxpoints )
116 {
117 lwerror("npoints (%d) is greater than maxpoints (%d)", pa->npoints, pa->maxpoints);
118 return LW_FAILURE;
119 }
120
121 /* Check if we have enough storage, add more if necessary */
122 if( pa->npoints == pa->maxpoints )
123 {
124 pa->maxpoints *= 2;
126 }
127
128 /* Make space to insert the new point */
129 if( where < pa->npoints )
130 {
131 size_t copy_size = point_size * (pa->npoints - where);
132 memmove(getPoint_internal(pa, where+1), getPoint_internal(pa, where), copy_size);
133 LWDEBUGF(5,"copying %zu bytes to start vertex %d from start vertex %d", copy_size, where+1, where);
134 }
135
136 /* We have one more point */
137 ++pa->npoints;
138
139 /* Copy the new point into the gap */
140 ptarray_set_point4d(pa, where, p);
141 LWDEBUGF(5,"copying new point to start vertex %zu", point_size);
142
143 return LW_SUCCESS;
144}
145
146int
147ptarray_append_point(POINTARRAY *pa, const POINT4D *pt, int repeated_points)
148{
149 /* Check for pathology */
150 if( ! pa || ! pt )
151 {
152 lwerror("ptarray_append_point: null input");
153 return LW_FAILURE;
154 }
155
156 /* Check for duplicate end point */
157 if ( repeated_points == LW_FALSE && pa->npoints > 0 )
158 {
159 POINT4D tmp;
160 getPoint4d_p(pa, pa->npoints-1, &tmp);
161 LWDEBUGF(4,"checking for duplicate end point (pt = POINT(%g %g) pa->npoints-q = POINT(%g %g))",pt->x,pt->y,tmp.x,tmp.y);
162
163 /* Return LW_SUCCESS and do nothing else if previous point in list is equal to this one */
164 if ( (pt->x == tmp.x) && (pt->y == tmp.y) &&
165 (FLAGS_GET_Z(pa->flags) ? pt->z == tmp.z : 1) &&
166 (FLAGS_GET_M(pa->flags) ? pt->m == tmp.m : 1) )
167 {
168 return LW_SUCCESS;
169 }
170 }
171
172 /* Append is just a special case of insert */
173 return ptarray_insert_point(pa, pt, pa->npoints);
174}
175
176int
177ptarray_append_ptarray(POINTARRAY *pa1, POINTARRAY *pa2, double gap_tolerance)
178{
179 unsigned int poff = 0;
180 unsigned int npoints;
181 unsigned int ncap;
182 unsigned int ptsize;
183
184 /* Check for pathology */
185 if( ! pa1 || ! pa2 )
186 {
187 lwerror("ptarray_append_ptarray: null input");
188 return LW_FAILURE;
189 }
190
191 npoints = pa2->npoints;
192
193 if ( ! npoints ) return LW_SUCCESS; /* nothing more to do */
194
195 if( FLAGS_GET_READONLY(pa1->flags) )
196 {
197 lwerror("ptarray_append_ptarray: target pointarray is read-only");
198 return LW_FAILURE;
199 }
200
201 if( FLAGS_GET_ZM(pa1->flags) != FLAGS_GET_ZM(pa2->flags) )
202 {
203 lwerror("ptarray_append_ptarray: appending mixed dimensionality is not allowed");
204 return LW_FAILURE;
205 }
206
207 ptsize = ptarray_point_size(pa1);
208
209 /* Check for duplicate end point */
210 if ( pa1->npoints )
211 {
212 POINT2D tmp1, tmp2;
213 getPoint2d_p(pa1, pa1->npoints-1, &tmp1);
214 getPoint2d_p(pa2, 0, &tmp2);
215
216 /* If the end point and start point are the same, then don't copy start point */
217 if (p2d_same(&tmp1, &tmp2)) {
218 poff = 1;
219 --npoints;
220 }
221 else if ( gap_tolerance == 0 || ( gap_tolerance > 0 &&
222 distance2d_pt_pt(&tmp1, &tmp2) > gap_tolerance ) )
223 {
224 lwerror("Second line start point too far from first line end point");
225 return LW_FAILURE;
226 }
227 }
228
229 /* Check if we need extra space */
230 ncap = pa1->npoints + npoints;
231 if ( pa1->maxpoints < ncap )
232 {
233 if ( pa1->maxpoints )
234 {
235 pa1->maxpoints = ncap > pa1->maxpoints*2 ?
236 ncap : pa1->maxpoints*2;
237 pa1->serialized_pointlist = lwrealloc(pa1->serialized_pointlist, (size_t)ptsize * pa1->maxpoints);
238 }
239 else
240 {
241 pa1->maxpoints = ncap;
242 pa1->serialized_pointlist = lwalloc((size_t)ptsize * pa1->maxpoints);
243 }
244 }
245
246 memcpy(getPoint_internal(pa1, pa1->npoints),
247 getPoint_internal(pa2, poff), (size_t)ptsize * npoints);
248
249 pa1->npoints = ncap;
250
251 return LW_SUCCESS;
252}
253
254/*
255* Add a point into a pointarray. Only adds as many dimensions as the
256* pointarray supports.
257*/
258int
260{
261 /* Check for pathology */
262 if( ! pa )
263 {
264 lwerror("ptarray_remove_point: null input");
265 return LW_FAILURE;
266 }
267
268 /* Error on invalid offset value */
269 if ( where >= pa->npoints )
270 {
271 lwerror("ptarray_remove_point: offset out of range (%d)", where);
272 return LW_FAILURE;
273 }
274
275 /* If the point is any but the last, we need to copy the data back one point */
276 if (where < pa->npoints - 1)
277 memmove(getPoint_internal(pa, where),
278 getPoint_internal(pa, where + 1),
279 ptarray_point_size(pa) * (pa->npoints - where - 1));
280
281 /* We have one less point */
282 pa->npoints--;
283
284 return LW_SUCCESS;
285}
286
291POINTARRAY* ptarray_construct_reference_data(char hasz, char hasm, uint32_t npoints, uint8_t *ptlist)
292{
293 POINTARRAY *pa = lwalloc(sizeof(POINTARRAY));
294 LWDEBUGF(5, "hasz = %d, hasm = %d, npoints = %d, ptlist = %p", hasz, hasm, npoints, ptlist);
295 pa->flags = lwflags(hasz, hasm, 0);
296 FLAGS_SET_READONLY(pa->flags, 1); /* We don't own this memory, so we can't alter or free it. */
297 pa->npoints = npoints;
298 pa->maxpoints = npoints;
299 pa->serialized_pointlist = ptlist;
300 return pa;
301}
302
303
305ptarray_construct_copy_data(char hasz, char hasm, uint32_t npoints, const uint8_t *ptlist)
306{
307 POINTARRAY *pa = lwalloc(sizeof(POINTARRAY));
308
309 pa->flags = lwflags(hasz, hasm, 0);
310 pa->npoints = npoints;
311 pa->maxpoints = npoints;
312
313 if ( npoints > 0 )
314 {
316 memcpy(pa->serialized_pointlist, ptlist, ptarray_point_size(pa) * npoints);
317 }
318 else
319 {
320 pa->serialized_pointlist = NULL;
321 }
322
323 return pa;
324}
325
326void
328{
329 if (pa)
330 {
333 lwfree(pa);
334 }
335}
336
337
338void
340{
341 if (!pa->npoints)
342 return;
343 uint32_t i;
344 uint32_t last = pa->npoints - 1;
345 uint32_t mid = pa->npoints / 2;
346
347 double *d = (double*)(pa->serialized_pointlist);
348 int j;
349 int ndims = FLAGS_NDIMS(pa->flags);
350 for (i = 0; i < mid; i++)
351 {
352 for (j = 0; j < ndims; j++)
353 {
354 double buf;
355 buf = d[i*ndims+j];
356 d[i*ndims+j] = d[(last-i)*ndims+j];
357 d[(last-i)*ndims+j] = buf;
358 }
359 }
360 return;
361}
362
363
369{
370 uint32_t i;
371 double d;
372 POINT4D p;
373
374 for (i=0 ; i < pa->npoints ; i++)
375 {
376 getPoint4d_p(pa, i, &p);
377 d = p.y;
378 p.y = p.x;
379 p.x = d;
380 ptarray_set_point4d(pa, i, &p);
381 }
382
383 return pa;
384}
385
386void
388{
389 uint32_t i;
390 double d, *dp1, *dp2;
391 POINT4D p;
392
393 dp1 = ((double*)&p)+(unsigned)o1;
394 dp2 = ((double*)&p)+(unsigned)o2;
395 for (i=0 ; i < pa->npoints ; i++)
396 {
397 getPoint4d_p(pa, i, &p);
398 d = *dp2;
399 *dp2 = *dp1;
400 *dp1 = d;
401 ptarray_set_point4d(pa, i, &p);
402 }
403}
404
413ptarray_segmentize2d(const POINTARRAY *ipa, double dist)
414{
415 double segdist;
416 POINT4D p1, p2;
417 POINT4D pbuf;
418 POINTARRAY *opa;
419 uint32_t i, j, nseg;
420 int hasz = FLAGS_GET_Z(ipa->flags);
421 int hasm = FLAGS_GET_M(ipa->flags);
422
423 pbuf.x = pbuf.y = pbuf.z = pbuf.m = 0;
424
425 /* Initial storage */
426 opa = ptarray_construct_empty(hasz, hasm, ipa->npoints);
427
428 /* Add first point */
429 getPoint4d_p(ipa, 0, &p1);
431
432 /* Loop on all other input points */
433 for (i = 1; i < ipa->npoints; i++)
434 {
435 /*
436 * We use these pointers to avoid
437 * "strict-aliasing rules break" warning raised
438 * by gcc (3.3 and up).
439 *
440 * It looks that casting a variable address (also
441 * referred to as "type-punned pointer")
442 * breaks those "strict" rules.
443 */
444 POINT4D *p1ptr=&p1, *p2ptr=&p2;
445 double segments;
446
447 getPoint4d_p(ipa, i, &p2);
448
449 segdist = distance2d_pt_pt((POINT2D *)p1ptr, (POINT2D *)p2ptr);
450 /* Split input segment into shorter even chunks */
451 segments = ceil(segdist / dist);
452
453 /* Uses INT32_MAX instead of UINT32_MAX to be safe that it fits */
454 if (segments >= INT32_MAX)
455 {
456 lwnotice("%s:%d - %s: Too many segments required (%e)",
457 __FILE__, __LINE__,__func__, segments);
458 ptarray_free(opa);
459 return NULL;
460 }
461 nseg = segments;
462
463 for (j = 1; j < nseg; j++)
464 {
465 pbuf.x = p1.x + (p2.x - p1.x) * j / nseg;
466 pbuf.y = p1.y + (p2.y - p1.y) * j / nseg;
467 if (hasz)
468 pbuf.z = p1.z + (p2.z - p1.z) * j / nseg;
469 if (hasm)
470 pbuf.m = p1.m + (p2.m - p1.m) * j / nseg;
471 ptarray_append_point(opa, &pbuf, LW_FALSE);
472 LW_ON_INTERRUPT(ptarray_free(opa); return NULL);
473 }
474
475 ptarray_append_point(opa, &p2, (ipa->npoints == 2) ? LW_TRUE : LW_FALSE);
476 p1 = p2;
477 LW_ON_INTERRUPT(ptarray_free(opa); return NULL);
478 }
479
480 return opa;
481}
482
483char
484ptarray_same(const POINTARRAY *pa1, const POINTARRAY *pa2)
485{
486 uint32_t i;
487 size_t ptsize;
488
489 if ( FLAGS_GET_ZM(pa1->flags) != FLAGS_GET_ZM(pa2->flags) ) return LW_FALSE;
490 LWDEBUG(5,"dimensions are the same");
491
492 if ( pa1->npoints != pa2->npoints ) return LW_FALSE;
493 LWDEBUG(5,"npoints are the same");
494
495 ptsize = ptarray_point_size(pa1);
496 LWDEBUGF(5, "ptsize = %zu", ptsize);
497
498 for (i=0; i<pa1->npoints; i++)
499 {
500 if ( memcmp(getPoint_internal(pa1, i), getPoint_internal(pa2, i), ptsize) )
501 return LW_FALSE;
502 LWDEBUGF(5,"point #%d is the same",i);
503 }
504
505 return LW_TRUE;
506}
507
508char
509ptarray_same2d(const POINTARRAY *pa1, const POINTARRAY *pa2)
510{
511 uint32_t i;
512
513 if ( FLAGS_GET_ZM(pa1->flags) != FLAGS_GET_ZM(pa2->flags) ) return LW_FALSE;
514 LWDEBUG(5,"dimensions are the same");
515
516 if ( pa1->npoints != pa2->npoints ) return LW_FALSE;
517 LWDEBUG(5,"npoints are the same");
518
519 for (i=0; i<pa1->npoints; i++)
520 {
521 if ( memcmp(getPoint_internal(pa1, i), getPoint_internal(pa2, i), sizeof(POINT2D)) )
522 return LW_FALSE;
523 LWDEBUGF(5,"point #%d is the same",i);
524 }
525
526 return LW_TRUE;
527}
528
530ptarray_addPoint(const POINTARRAY *pa, uint8_t *p, size_t pdims, uint32_t where)
531{
532 POINTARRAY *ret;
533 POINT4D pbuf;
534 size_t ptsize = ptarray_point_size(pa);
535
536 LWDEBUGF(3, "pa %p p %p size %zu where %u",
537 pa, p, pdims, where);
538
539 if ( pdims < 2 || pdims > 4 )
540 {
541 lwerror("ptarray_addPoint: point dimension out of range (%zu)",
542 pdims);
543 return NULL;
544 }
545
546 if ( where > pa->npoints )
547 {
548 lwerror("ptarray_addPoint: offset out of range (%d)",
549 where);
550 return NULL;
551 }
552
553 LWDEBUG(3, "called with a point");
554
555 pbuf.x = pbuf.y = pbuf.z = pbuf.m = 0.0;
556 memcpy((uint8_t *)&pbuf, p, pdims*sizeof(double));
557
558 LWDEBUG(3, "initialized point buffer");
559
561 FLAGS_GET_M(pa->flags), pa->npoints+1);
562
563
564 if ( where )
565 {
566 memcpy(getPoint_internal(ret, 0), getPoint_internal(pa, 0), ptsize*where);
567 }
568
569 memcpy(getPoint_internal(ret, where), (uint8_t *)&pbuf, ptsize);
570
571 if ( where+1 != ret->npoints )
572 {
573 memcpy(getPoint_internal(ret, where+1),
574 getPoint_internal(pa, where),
575 ptsize*(pa->npoints-where));
576 }
577
578 return ret;
579}
580
582ptarray_removePoint(POINTARRAY *pa, uint32_t which)
583{
584 POINTARRAY *ret;
585 size_t ptsize = ptarray_point_size(pa);
586
587 LWDEBUGF(3, "pa %p which %u", pa, which);
588
589 assert(which <= pa->npoints-1);
590 assert(pa->npoints >= 3);
591
593 FLAGS_GET_M(pa->flags), pa->npoints-1);
594
595 /* copy initial part */
596 if ( which )
597 {
598 memcpy(getPoint_internal(ret, 0), getPoint_internal(pa, 0), ptsize*which);
599 }
600
601 /* copy final part */
602 if ( which < pa->npoints-1 )
603 {
604 memcpy(getPoint_internal(ret, which), getPoint_internal(pa, which+1),
605 ptsize*(pa->npoints-which-1));
606 }
607
608 return ret;
609}
610
613{
614 POINTARRAY *pa;
615 size_t ptsize = ptarray_point_size(pa1);
616
617 if (FLAGS_GET_ZM(pa1->flags) != FLAGS_GET_ZM(pa2->flags))
618 lwerror("ptarray_cat: Mixed dimension");
619
621 FLAGS_GET_M(pa1->flags),
622 pa1->npoints + pa2->npoints);
623
624 memcpy( getPoint_internal(pa, 0),
625 getPoint_internal(pa1, 0),
626 ptsize*(pa1->npoints));
627
628 memcpy( getPoint_internal(pa, pa1->npoints),
629 getPoint_internal(pa2, 0),
630 ptsize*(pa2->npoints));
631
632 ptarray_free(pa1);
633 ptarray_free(pa2);
634
635 return pa;
636}
637
638
644{
645 POINTARRAY *out = lwalloc(sizeof(POINTARRAY));
646
647 LWDEBUG(3, "ptarray_clone_deep called.");
648
649 out->flags = in->flags;
650 out->npoints = in->npoints;
651 out->maxpoints = in->npoints;
652
653 FLAGS_SET_READONLY(out->flags, 0);
654
655 if (!in->npoints)
656 {
657 // Avoid calling lwalloc of 0 bytes
658 out->serialized_pointlist = NULL;
659 }
660 else
661 {
662 size_t size = in->npoints * ptarray_point_size(in);
663 out->serialized_pointlist = lwalloc(size);
664 memcpy(out->serialized_pointlist, in->serialized_pointlist, size);
665 }
666
667 return out;
668}
669
675{
676 POINTARRAY *out = lwalloc(sizeof(POINTARRAY));
677
678 LWDEBUG(3, "ptarray_clone called.");
679
680 out->flags = in->flags;
681 out->npoints = in->npoints;
682 out->maxpoints = in->maxpoints;
683
684 FLAGS_SET_READONLY(out->flags, 1);
685
687
688 return out;
689}
690
695int
697{
698 if (!in)
699 {
700 lwerror("ptarray_is_closed: called with null point array");
701 return 0;
702 }
703 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
704
705 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), ptarray_point_size(in));
706}
707
708
709int
711{
712 if (!in)
713 {
714 lwerror("ptarray_is_closed_2d: called with null point array");
715 return 0;
716 }
717 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
718
719 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), sizeof(POINT2D) );
720}
721
722int
724{
725 if (!in)
726 {
727 lwerror("ptarray_is_closed_3d: called with null point array");
728 return 0;
729 }
730 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
731
732 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), sizeof(POINT3D) );
733}
734
735int
737{
738 if ( FLAGS_GET_Z(in->flags) )
739 return ptarray_is_closed_3d(in);
740 else
741 return ptarray_is_closed_2d(in);
742}
743
754int
756{
757 const POINT2D *seg1, *seg2;
758 int wn = 0;
759
760 seg1 = getPoint2d_cp(pa, 0);
761 seg2 = getPoint2d_cp(pa, pa->npoints-1);
762 if (!p2d_same(seg1, seg2))
763 lwerror("ptarray_contains_point called on unclosed ring");
764
765 for (uint32_t i = 1; i < pa->npoints; i++)
766 {
767 double side, ymin, ymax;
768
769 seg2 = getPoint2d_cp(pa, i);
770
771 /* Zero length segments are ignored. */
772 if (p2d_same(seg1, seg2))
773 {
774 seg1 = seg2;
775 continue;
776 }
777
778 ymin = FP_MIN(seg1->y, seg2->y);
779 ymax = FP_MAX(seg1->y, seg2->y);
780
781 /* Only test segments in our vertical range */
782 if (pt->y > ymax || pt->y < ymin)
783 {
784 seg1 = seg2;
785 continue;
786 }
787
788 side = lw_segment_side(seg1, seg2, pt);
789
790 /*
791 * A point on the boundary of a ring is not contained.
792 * WAS: if (fabs(side) < 1e-12), see #852
793 */
794 if ((side == 0) && lw_pt_in_seg(pt, seg1, seg2))
795 {
796 return LW_BOUNDARY;
797 }
798
799 /*
800 * If the point is to the left of the line, and it's rising,
801 * then the line is to the right of the point and
802 * circling counter-clockwise, so increment.
803 */
804 if ((side < 0) && (seg1->y <= pt->y) && (pt->y < seg2->y))
805 {
806 wn++;
807 }
808
809 /*
810 * If the point is to the right of the line, and it's falling,
811 * then the line is to the right of the point and circling
812 * clockwise, so decrement.
813 */
814 else if ( (side > 0) && (seg2->y <= pt->y) && (pt->y < seg1->y) )
815 {
816 wn--;
817 }
818
819 seg1 = seg2;
820 }
821
822 /* wn == 0 => Outside, wn != 0 => Inside */
823 return wn == 0 ? LW_OUTSIDE : LW_INSIDE;
824}
825
826int
827ptarray_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
828{
829 // A valid linestring must have at least 2 point
830 if (pa->npoints < 2)
831 lwerror("%s called on invalid linestring", __func__);
832
833 int intersections = 0;
834 double px = p->x;
835 double py = p->y;
836
837 // Iterate through each edge of the polygon
838 for (uint32_t i = 0; i < pa->npoints-1; ++i)
839 {
840 const POINT2D* p1 = getPoint2d_cp(pa, i);
841 const POINT2D* p2 = getPoint2d_cp(pa, i+1);
842
843 /* Skip zero-length edges */
844 if (p2d_same(p1, p2))
845 continue;
846
847 /* --- Step 1: Check if the point is ON the boundary edge --- */
848 if (lw_pt_on_segment(p1, p2, p))
849 {
850 *on_boundary = LW_TRUE;
851 return 0;
852 }
853
854 /* --- Step 2: Perform the Ray Casting intersection test --- */
855
856 /*
857 * Check if the horizontal ray from p intersects the edge (p1, p2).
858 * This is the core condition for handling vertices correctly:
859 * - One vertex must be strictly above the ray (py < vertex.y)
860 * - The other must be on or below the ray (py >= vertex.y)
861 */
862 if (((p1->y <= py) && (py < p2->y)) || ((p2->y <= py) && (py < p1->y)))
863 {
864 /*
865 * Calculate the x-coordinate where the edge intersects the ray's horizontal line.
866 * Formula: x = x1 + (y - y1) * (x2 - x1) / (y2 - y1)
867 */
868 double x_intersection = p1->x + (py - p1->y) * (p2->x - p1->x) / (p2->y - p1->y);
869
870 /*
871 * If the intersection point is to the right of our test point,
872 * it's a valid "crossing".
873 */
874 if (x_intersection > px)
875 {
876 intersections++;
877 }
878 }
879 }
880
881 return intersections;
882}
883
884
902static int
903circle_raycast_intersections(const POINT2D *center, double radius, double ray, POINT2D *i0, POINT2D *i1)
904{
905 // Calculate the vertical distance from the circle's center to the horizontal line.
906 double dy = ray - center->y;
907
908 // If the absolute vertical distance is greater than the radius, there are no intersections.
909 if (fabs(dy) > radius)
910 return 0;
911
912 // Use the Pythagorean theorem to find the horizontal distance (dx) from the
913 // center's x-coordinate to the intersection points.
914 // dx^2 + dy^2 = radius^2 => dx^2 = radius^2 - dy^2
915 double dx_squared = radius * radius - dy * dy;
916
917 // Case 1: One intersection (tangent)
918 // This occurs when the line just touches the top or bottom of the circle.
919 // dx_squared will be zero. We check against a small epsilon for floating-point safety.
920 if (FP_EQUALS(dx_squared, 0.0))
921 {
922 i0->x = center->x;
923 i0->y = ray;
924 return 1;
925 }
926
927 // Case 2: Two intersections
928 // The line cuts through the circle.
929 double dx = sqrt(dx_squared);
930
931 // The first intersection point has the smaller x-value.
932 i0->x = center->x - dx;
933 i0->y = ray;
934
935 // The second intersection point has the larger x-value.
936 i1->x = center->x + dx;
937 i1->y = ray;
938
939 return 2;
940}
941
942
943int
944ptarrayarc_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
945{
946 int intersections = 0;
947 double px = p->x;
948 double py = p->y;
949
950 assert(on_boundary);
951
952 // A valid circular arc must have at least 3 vertices (circle).
953 if (pa->npoints < 3)
954 lwerror("%s called on invalid circularstring", __func__);
955
956 if (pa->npoints % 2 == 0)
957 lwerror("%s called with even number of points", __func__);
958
959 // Iterate through each arc of the circularstring
960 for (uint32_t i = 1; i < pa->npoints-1; i +=2)
961 {
962 const POINT2D* p0 = getPoint2d_cp(pa, i-1);
963 const POINT2D* p1 = getPoint2d_cp(pa, i);
964 const POINT2D* p2 = getPoint2d_cp(pa, i+1);
965 POINT2D center = {0,0};
966 double radius, d;
967 GBOX gbox;
968
969 // Skip zero-length arc
970 if (lw_arc_is_pt(p0, p1, p2))
971 continue;
972
973 // --- Step 1: Check if the point is ON the boundary edge ---
974 if (p2d_same(p0, p) || p2d_same(p1, p) || p2d_same(p2, p))
975 {
976 *on_boundary = LW_TRUE;
977 return 0;
978 }
979
980 // Calculate some important pieces
981 radius = lw_arc_center(p0, p1, p2, &center);
982
983 d = distance2d_pt_pt(p, &center);
984 if (FP_EQUALS(d, radius) && lw_pt_in_arc(p, p0, p1, p2))
985 {
986 *on_boundary = LW_TRUE;
987 return 0;
988 }
989
990 // --- Step 2: Perform the Ray Casting intersection test ---
991
992 // Only process arcs that our ray crosses
993 lw_arc_calculate_gbox_cartesian_2d(p0, p1, p2, &gbox);
994 if ((gbox.ymin <= py) && (py < gbox.ymax))
995 {
996 // Point of intersection on the circle that defines the arc
997 POINT2D i0, i1;
998
999 // How many points of intersection are there
1000 int iCount = circle_raycast_intersections(&center, radius, py, &i0, &i1);
1001
1002 // Nothing to see here
1003 if (iCount == 0)
1004 continue;
1005
1006 // Cannot think of a case where a grazing is not a
1007 // no-op
1008 if (iCount == 1)
1009 continue;
1010
1011 // So we must have 2 intersections
1012 // Only increment the counter for intersections to the right
1013 // of the test point
1014 if (i0.x > px && lw_pt_in_arc(&i0, p0, p1, p2))
1015 intersections++;
1016
1017 if (i1.x > px && lw_pt_in_arc(&i1, p0, p1, p2))
1018 intersections++;
1019 }
1020 }
1021
1022 return intersections;
1023}
1024
1025/*
1026 * The following is based on the "Fast Winding Number Inclusion of a Point
1027 * in a Polygon" algorithm by Dan Sunday.
1028 * http://softsurfer.com/Archive/algorithm_0103/algorithm_0103.htm#Winding%20Number
1029 */
1030int
1031ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
1032{
1033 int wn = 0;
1034 uint32_t i;
1035 double side;
1036 const POINT2D *seg1, *seg2;
1037 double ymin, ymax;
1038
1039 seg1 = getPoint2d_cp(pa, 0);
1040 seg2 = getPoint2d_cp(pa, pa->npoints-1);
1041 if ( check_closed && ! p2d_same(seg1, seg2) )
1042 lwerror("ptarray_contains_point called on unclosed ring");
1043
1044 for ( i=1; i < pa->npoints; i++ )
1045 {
1046 seg2 = getPoint2d_cp(pa, i);
1047
1048 /* Zero length segments are ignored. */
1049 if ( seg1->x == seg2->x && seg1->y == seg2->y )
1050 {
1051 seg1 = seg2;
1052 continue;
1053 }
1054
1055 ymin = FP_MIN(seg1->y, seg2->y);
1056 ymax = FP_MAX(seg1->y, seg2->y);
1057
1058 /* Only test segments in our vertical range */
1059 if ( pt->y > ymax || pt->y < ymin )
1060 {
1061 seg1 = seg2;
1062 continue;
1063 }
1064
1065 side = lw_segment_side(seg1, seg2, pt);
1066
1067 /*
1068 * A point on the boundary of a ring is not contained.
1069 * WAS: if (fabs(side) < 1e-12), see #852
1070 */
1071 if ( (side == 0) && lw_pt_in_seg(pt, seg1, seg2) )
1072 {
1073 return LW_BOUNDARY;
1074 }
1075
1076 /*
1077 * If the point is to the left of the line, and it's rising,
1078 * then the line is to the right of the point and
1079 * circling counter-clockwise, so increment.
1080 */
1081 if ( (side < 0) && (seg1->y <= pt->y) && (pt->y < seg2->y) )
1082 {
1083 wn++;
1084 }
1085
1086 /*
1087 * If the point is to the right of the line, and it's falling,
1088 * then the line is to the right of the point and circling
1089 * clockwise, so decrement.
1090 */
1091 else if ( (side > 0) && (seg2->y <= pt->y) && (pt->y < seg1->y) )
1092 {
1093 wn--;
1094 }
1095
1096 seg1 = seg2;
1097 }
1098
1099 /* Sent out the winding number for calls that are building on this as a primitive */
1100 if ( winding_number )
1101 *winding_number = wn;
1102
1103 /* Outside */
1104 if (wn == 0)
1105 {
1106 return LW_OUTSIDE;
1107 }
1108
1109 /* Inside */
1110 return LW_INSIDE;
1111}
1112
1122int
1124{
1125 int on_boundary = LW_FALSE;
1126 int intersections;
1127 if (!ptarray_is_closed_2d(pa))
1128 lwerror("%s called on unclosed ring", __func__);
1129
1130 intersections = ptarrayarc_raycast_intersections(pa, pt, &on_boundary);
1131 if (on_boundary)
1132 return LW_BOUNDARY;
1133 else
1134 return (intersections % 2) ? LW_INSIDE : LW_OUTSIDE;
1135}
1136
1142double
1144{
1145 const POINT2D *P1;
1146 const POINT2D *P2;
1147 const POINT2D *P3;
1148 double sum = 0.0;
1149 double x0, x, y1, y2;
1150 uint32_t i;
1151
1152 if (! pa || pa->npoints < 3 )
1153 return 0.0;
1154
1155 P1 = getPoint2d_cp(pa, 0);
1156 P2 = getPoint2d_cp(pa, 1);
1157 x0 = P1->x;
1158 for ( i = 2; i < pa->npoints; i++ )
1159 {
1160 P3 = getPoint2d_cp(pa, i);
1161 x = P2->x - x0;
1162 y1 = P3->y;
1163 y2 = P1->y;
1164 sum += x * (y2-y1);
1165
1166 /* Move forwards! */
1167 P1 = P2;
1168 P2 = P3;
1169 }
1170 return sum / 2.0;
1171}
1172
1173int
1175{
1176 double area = 0;
1177 area = ptarray_signed_area(pa);
1178 if ( area > 0 ) return LW_FALSE;
1179 else return LW_TRUE;
1180}
1181
1183ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm, double zval, double mval)
1184{
1185 /* TODO handle zero-length point arrays */
1186 uint32_t i;
1187 int in_hasz = FLAGS_GET_Z(pa->flags);
1188 int in_hasm = FLAGS_GET_M(pa->flags);
1189 POINT4D pt;
1190 POINTARRAY *pa_out = ptarray_construct_empty(hasz, hasm, pa->npoints);
1191
1192 for( i = 0; i < pa->npoints; i++ )
1193 {
1194 getPoint4d_p(pa, i, &pt);
1195 if( hasz && ! in_hasz )
1196 pt.z = zval;
1197 if( hasm && ! in_hasm )
1198 pt.m = mval;
1199 ptarray_append_point(pa_out, &pt, LW_TRUE);
1200 }
1201
1202 return pa_out;
1203}
1204
1205POINTARRAY *
1206ptarray_substring(POINTARRAY *ipa, double from, double to, double tolerance)
1207{
1208 POINTARRAY *dpa;
1209 POINT4D pt;
1210 POINT4D p1, p2;
1211 POINT4D *p1ptr=&p1; /* don't break strict-aliasing rule */
1212 POINT4D *p2ptr=&p2;
1213 int nsegs, i;
1214 double length, slength, tlength;
1215 int state = 0; /* 0=before, 1=inside */
1216
1217 /*
1218 * Create a dynamic pointarray with an initial capacity
1219 * equal to full copy of input points
1220 */
1222
1223 /* Compute total line length */
1224 length = ptarray_length_2d(ipa);
1225
1226
1227 LWDEBUGF(3, "Total length: %g", length);
1228
1229
1230 /* Get 'from' and 'to' lengths */
1231 from = length*from;
1232 to = length*to;
1233
1234
1235 LWDEBUGF(3, "From/To: %g/%g", from, to);
1236
1237
1238 tlength = 0;
1239 getPoint4d_p(ipa, 0, &p1);
1240 nsegs = ipa->npoints - 1;
1241 for ( i = 0; i < nsegs; i++ )
1242 {
1243 double dseg;
1244
1245 getPoint4d_p(ipa, i+1, &p2);
1246
1247
1248 LWDEBUGF(3 ,"Segment %d: (%g,%g,%g,%g)-(%g,%g,%g,%g)",
1249 i, p1.x, p1.y, p1.z, p1.m, p2.x, p2.y, p2.z, p2.m);
1250
1251
1252 /* Find the length of this segment */
1253 slength = distance2d_pt_pt((POINT2D *)p1ptr, (POINT2D *)p2ptr);
1254
1255 /*
1256 * We are before requested start.
1257 */
1258 if ( state == 0 ) /* before */
1259 {
1260
1261 LWDEBUG(3, " Before start");
1262
1263 if ( fabs ( from - ( tlength + slength ) ) <= tolerance )
1264 {
1265
1266 LWDEBUG(3, " Second point is our start");
1267
1268 /*
1269 * Second point is our start
1270 */
1271 ptarray_append_point(dpa, &p2, LW_FALSE);
1272 state=1; /* we're inside now */
1273 goto END;
1274 }
1275
1276 else if ( fabs(from - tlength) <= tolerance )
1277 {
1278
1279 LWDEBUG(3, " First point is our start");
1280
1281 /*
1282 * First point is our start
1283 */
1284 ptarray_append_point(dpa, &p1, LW_FALSE);
1285
1286 /*
1287 * We're inside now, but will check
1288 * 'to' point as well
1289 */
1290 state=1;
1291 }
1292
1293 /*
1294 * Didn't reach the 'from' point,
1295 * nothing to do
1296 */
1297 else if ( from > tlength + slength ) goto END;
1298
1299 else /* tlength < from < tlength+slength */
1300 {
1301
1302 LWDEBUG(3, " Seg contains first point");
1303
1304 /*
1305 * Our start is between first and
1306 * second point
1307 */
1308 dseg = (from - tlength) / slength;
1309
1310 interpolate_point4d(&p1, &p2, &pt, dseg);
1311
1312 ptarray_append_point(dpa, &pt, LW_FALSE);
1313
1314 /*
1315 * We're inside now, but will check
1316 * 'to' point as well
1317 */
1318 state=1;
1319 }
1320 }
1321
1322 if ( state == 1 ) /* inside */
1323 {
1324
1325 LWDEBUG(3, " Inside");
1326
1327 /*
1328 * 'to' point is our second point.
1329 */
1330 if ( fabs(to - ( tlength + slength ) ) <= tolerance )
1331 {
1332
1333 LWDEBUG(3, " Second point is our end");
1334
1335 ptarray_append_point(dpa, &p2, LW_FALSE);
1336 break; /* substring complete */
1337 }
1338
1339 /*
1340 * 'to' point is our first point.
1341 * (should only happen if 'to' is 0)
1342 */
1343 else if ( fabs(to - tlength) <= tolerance )
1344 {
1345
1346 LWDEBUG(3, " First point is our end");
1347
1348 ptarray_append_point(dpa, &p1, LW_FALSE);
1349
1350 break; /* substring complete */
1351 }
1352
1353 /*
1354 * Didn't reach the 'end' point,
1355 * just copy second point
1356 */
1357 else if ( to > tlength + slength )
1358 {
1359 ptarray_append_point(dpa, &p2, LW_FALSE);
1360 goto END;
1361 }
1362
1363 /*
1364 * 'to' point falls on this segment
1365 * Interpolate and break.
1366 */
1367 else if ( to < tlength + slength )
1368 {
1369
1370 LWDEBUG(3, " Seg contains our end");
1371
1372 dseg = (to - tlength) / slength;
1373 interpolate_point4d(&p1, &p2, &pt, dseg);
1374
1375 ptarray_append_point(dpa, &pt, LW_FALSE);
1376
1377 break;
1378 }
1379
1380 else
1381 {
1382 LWDEBUG(3, "Unhandled case");
1383 }
1384 }
1385
1386
1387END:
1388
1389 tlength += slength;
1390 memcpy(&p1, &p2, sizeof(POINT4D));
1391 }
1392
1393 LWDEBUGF(3, "Out of loop, ptarray has %d points", dpa->npoints);
1394
1395 return dpa;
1396}
1397
1398/*
1399 * Write into the *ret argument coordinates of the closes point on
1400 * the given segment to the reference input point.
1401 */
1402void
1403closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
1404{
1405 double r;
1406
1407 if ( FP_EQUALS(A->x, B->x) && FP_EQUALS(A->y, B->y) )
1408 {
1409 *ret = *A;
1410 return;
1411 }
1412
1413 /*
1414 * We use comp.graphics.algorithms Frequently Asked Questions method
1415 *
1416 * (1) AC dot AB
1417 * r = ----------
1418 * ||AB||^2
1419 * r has the following meaning:
1420 * r=0 P = A
1421 * r=1 P = B
1422 * r<0 P is on the backward extension of AB
1423 * r>1 P is on the forward extension of AB
1424 * 0<r<1 P is interior to AB
1425 *
1426 */
1427 r = ( (p->x-A->x) * (B->x-A->x) + (p->y-A->y) * (B->y-A->y) )/( (B->x-A->x)*(B->x-A->x) +(B->y-A->y)*(B->y-A->y) );
1428
1429 if (r<=0)
1430 {
1431 *ret = *A;
1432 return;
1433 }
1434 if (r>=1)
1435 {
1436 *ret = *B;
1437 return;
1438 }
1439
1440 ret->x = A->x + ( (B->x - A->x) * r );
1441 ret->y = A->y + ( (B->y - A->y) * r );
1442 ret->z = A->z + ( (B->z - A->z) * r );
1443 ret->m = A->m + ( (B->m - A->m) * r );
1444}
1445
1446int
1447ptarray_closest_segment_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
1448{
1449 const POINT2D *start = getPoint2d_cp(pa, 0), *end = NULL;
1450 uint32_t t, seg=0;
1451 double mindist=DBL_MAX;
1452
1453 /* Loop through pointarray looking for nearest segment */
1454 for (t=1; t<pa->npoints; t++)
1455 {
1456 double dist_sqr;
1457 end = getPoint2d_cp(pa, t);
1458 dist_sqr = distance2d_sqr_pt_seg(qp, start, end);
1459
1460 if (dist_sqr < mindist)
1461 {
1462 mindist = dist_sqr;
1463 seg=t-1;
1464 if ( mindist == 0 )
1465 {
1466 LWDEBUG(3, "Breaking on mindist=0");
1467 break;
1468 }
1469 }
1470
1471 start = end;
1472 }
1473
1474 if ( dist ) *dist = sqrt(mindist);
1475 return seg;
1476}
1477
1478
1479int
1480ptarray_closest_vertex_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
1481{
1482 uint32_t t, pn=0;
1483 const POINT2D *p;
1484 double mindist = DBL_MAX;
1485
1486 /* Loop through pointarray looking for nearest segment */
1487 for (t=0; t<pa->npoints; t++)
1488 {
1489 double dist_sqr;
1490 p = getPoint2d_cp(pa, t);
1491 dist_sqr = distance2d_sqr_pt_pt(p, qp);
1492
1493 if (dist_sqr < mindist)
1494 {
1495 mindist = dist_sqr;
1496 pn = t;
1497 if ( mindist == 0 )
1498 {
1499 LWDEBUG(3, "Breaking on mindist=0");
1500 break;
1501 }
1502 }
1503 }
1504 if ( dist ) *dist = sqrt(mindist);
1505 return pn;
1506}
1507
1508/*
1509 * Given a point, returns the location of closest point on pointarray
1510 * and, optionally, it's actual distance from the point array.
1511 */
1512double
1513ptarray_locate_point(const POINTARRAY *pa, const POINT4D *p4d, double *mindistout, POINT4D *proj4d)
1514{
1515 double mindist=DBL_MAX;
1516 double tlen, plen;
1517 uint32_t t, seg=0;
1518 POINT4D start4d, end4d, projtmp;
1519 POINT2D proj, p;
1520 const POINT2D *start = NULL, *end = NULL;
1521
1522 /* Initialize our 2D copy of the input parameter */
1523 p.x = p4d->x;
1524 p.y = p4d->y;
1525
1526 if ( ! proj4d ) proj4d = &projtmp;
1527
1528 /* Check for special cases (length 0 and 1) */
1529 if ( pa->npoints <= 1 )
1530 {
1531 if ( pa->npoints == 1 )
1532 {
1533 getPoint4d_p(pa, 0, proj4d);
1534 if ( mindistout )
1535 *mindistout = distance2d_pt_pt(&p, getPoint2d_cp(pa, 0));
1536 }
1537 return 0.0;
1538 }
1539
1540 start = getPoint2d_cp(pa, 0);
1541 /* Loop through pointarray looking for nearest segment */
1542 for (t=1; t<pa->npoints; t++)
1543 {
1544 double dist_sqr;
1545 end = getPoint2d_cp(pa, t);
1546 dist_sqr = distance2d_sqr_pt_seg(&p, start, end);
1547
1548 if (dist_sqr < mindist)
1549 {
1550 mindist = dist_sqr;
1551 seg=t-1;
1552 if ( mindist == 0 )
1553 {
1554 LWDEBUG(3, "Breaking on mindist=0");
1555 break;
1556 }
1557 }
1558
1559 start = end;
1560 }
1561 mindist = sqrt(mindist);
1562
1563 if ( mindistout ) *mindistout = mindist;
1564
1565 LWDEBUGF(3, "Closest segment: %d", seg);
1566 LWDEBUGF(3, "mindist: %g", mindist);
1567
1568 /*
1569 * We need to project the
1570 * point on the closest segment.
1571 */
1572 getPoint4d_p(pa, seg, &start4d);
1573 getPoint4d_p(pa, seg+1, &end4d);
1574 closest_point_on_segment(p4d, &start4d, &end4d, proj4d);
1575
1576 /* Copy 4D values into 2D holder */
1577 proj.x = proj4d->x;
1578 proj.y = proj4d->y;
1579
1580 LWDEBUGF(3, "Closest segment:%d, npoints:%d", seg, pa->npoints);
1581
1582 /* For robustness, force 1 when closest point == endpoint */
1583 if ( (seg >= (pa->npoints-2)) && p2d_same(&proj, end) )
1584 {
1585 return 1.0;
1586 }
1587
1588 LWDEBUGF(3, "Closest point on segment: %g,%g", proj.x, proj.y);
1589
1590 tlen = ptarray_length_2d(pa);
1591
1592 LWDEBUGF(3, "tlen %g", tlen);
1593
1594 /* Location of any point on a zero-length line is 0 */
1595 /* See http://trac.osgeo.org/postgis/ticket/1772#comment:2 */
1596 if ( tlen == 0 ) return 0;
1597
1598 plen=0;
1599 start = getPoint2d_cp(pa, 0);
1600 for (t=0; t<seg; t++, start=end)
1601 {
1602 end = getPoint2d_cp(pa, t+1);
1603 plen += distance2d_pt_pt(start, end);
1604
1605 LWDEBUGF(4, "Segment %d made plen %g", t, plen);
1606 }
1607
1608 plen+=distance2d_pt_pt(&proj, start);
1609
1610 LWDEBUGF(3, "plen %g, tlen %g", plen, tlen);
1611
1612 return plen/tlen;
1613}
1614
1624void
1626{
1627 uint32_t i;
1628 double x;
1629
1630 for (i=0; i<pa->npoints; i++)
1631 {
1632 memcpy(&x, getPoint_internal(pa, i), sizeof(double));
1633 if ( x < 0 ) x+= 360;
1634 else if ( x > 180 ) x -= 360;
1635 memcpy(getPoint_internal(pa, i), &x, sizeof(double));
1636 }
1637}
1638
1639
1640/*
1641 * Returns a POINTARRAY with consecutive equal points
1642 * removed. Equality test on all dimensions of input.
1643 *
1644 * Always returns a newly allocated object.
1645 */
1646static POINTARRAY *
1647ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
1648{
1649 POINTARRAY *out = ptarray_clone_deep(in);
1650 ptarray_remove_repeated_points_in_place(out, tolerance, minpoints);
1651 return out;
1652}
1653
1654POINTARRAY *
1655ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
1656{
1657 return ptarray_remove_repeated_points_minpoints(in, tolerance, 2);
1658}
1659
1660
1661void
1662ptarray_remove_repeated_points_in_place(POINTARRAY *pa, double tolerance, uint32_t min_points)
1663{
1664 uint32_t i;
1665 double tolsq = tolerance * tolerance;
1666 const POINT2D *last = NULL;
1667 const POINT2D *pt;
1668 uint32_t n_points = pa->npoints;
1669 uint32_t n_points_out = 1;
1670 size_t pt_size = ptarray_point_size(pa);
1671
1672 double dsq = FLT_MAX;
1673
1674 /* No-op on short inputs */
1675 if ( n_points <= min_points ) return;
1676
1677 last = getPoint2d_cp(pa, 0);
1678 void *p_to = ((char *)last) + pt_size;
1679 for (i = 1; i < n_points; i++)
1680 {
1681 int last_point = (i == n_points - 1);
1682
1683 /* Look straight into the abyss */
1684 pt = getPoint2d_cp(pa, i);
1685
1686 /* Don't drop points if we are running short of points */
1687 if (n_points + n_points_out > min_points + i)
1688 {
1689 if (tolerance > 0.0)
1690 {
1691 /* Only drop points that are within our tolerance */
1692 dsq = distance2d_sqr_pt_pt(last, pt);
1693 /* Allow any point but the last one to be dropped */
1694 if (!last_point && dsq <= tolsq)
1695 {
1696 continue;
1697 }
1698 }
1699 else
1700 {
1701 /* At tolerance zero, only skip exact dupes */
1702 if (memcmp((char*)pt, (char*)last, pt_size) == 0)
1703 continue;
1704 }
1705
1706 /* Got to last point, and it's not very different from */
1707 /* the point that preceded it. We want to keep the last */
1708 /* point, not the second-to-last one, so we pull our write */
1709 /* index back one value */
1710 if (last_point && n_points_out > 1 && tolerance > 0.0 && dsq <= tolsq)
1711 {
1712 n_points_out--;
1713 p_to = (char*)p_to - pt_size;
1714 }
1715 }
1716
1717 /* Compact all remaining values to front of array */
1718 memcpy(p_to, pt, pt_size);
1719 n_points_out++;
1720 p_to = (char*)p_to + pt_size;
1721 last = pt;
1722 }
1723 /* Adjust array length */
1724 pa->npoints = n_points_out;
1725 return;
1726}
1727
1728/* Out of the points in pa [itfist .. itlast], finds the one that's farthest away from
1729 * the segment determined by pts[itfist] and pts[itlast].
1730 * Returns itfirst if no point was found further away than max_distance_sqr
1731 */
1732static uint32_t
1733ptarray_dp_findsplit_in_place(const POINTARRAY *pts, uint32_t it_first, uint32_t it_last, double max_distance_sqr)
1734{
1735 uint32_t split = it_first;
1736 if ((it_first - it_last) < 2)
1737 return it_first;
1738
1739 const POINT2D *A = getPoint2d_cp(pts, it_first);
1740 const POINT2D *B = getPoint2d_cp(pts, it_last);
1741
1742 if (distance2d_sqr_pt_pt(A, B) < DBL_EPSILON)
1743 {
1744 /* If p1 == p2, we can just calculate the distance from each point to A */
1745 for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1746 {
1747 const POINT2D *pk = getPoint2d_cp(pts, itk);
1748 double distance_sqr = distance2d_sqr_pt_pt(pk, A);
1749 if (distance_sqr > max_distance_sqr)
1750 {
1751 split = itk;
1752 max_distance_sqr = distance_sqr;
1753 }
1754 }
1755 return split;
1756 }
1757
1758 /* This is based on distance2d_sqr_pt_seg, but heavily inlined here to avoid recalculations */
1759 double ba_x = (B->x - A->x);
1760 double ba_y = (B->y - A->y);
1761 double ab_length_sqr = (ba_x * ba_x + ba_y * ba_y);
1762 /* To avoid the division by ab_length_sqr in the 3rd path, we normalize here
1763 * and multiply in the first two paths [(dot_ac_ab < 0) and (> ab_length_sqr)] */
1764 max_distance_sqr *= ab_length_sqr;
1765 for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1766 {
1767 const POINT2D *C = getPoint2d_cp(pts, itk);
1768 double distance_sqr;
1769 double ca_x = (C->x - A->x);
1770 double ca_y = (C->y - A->y);
1771 double dot_ac_ab = (ca_x * ba_x + ca_y * ba_y);
1772
1773 if (dot_ac_ab <= 0.0)
1774 {
1775 distance_sqr = distance2d_sqr_pt_pt(C, A) * ab_length_sqr;
1776 }
1777 else if (dot_ac_ab >= ab_length_sqr)
1778 {
1779 distance_sqr = distance2d_sqr_pt_pt(C, B) * ab_length_sqr;
1780 }
1781 else
1782 {
1783 double s_numerator = ca_x * ba_y - ca_y * ba_x;
1784 distance_sqr = s_numerator * s_numerator; /* Missing division by ab_length_sqr on purpose */
1785 }
1786
1787 if (distance_sqr > max_distance_sqr)
1788 {
1789 split = itk;
1790 max_distance_sqr = distance_sqr;
1791 }
1792 }
1793 return split;
1794}
1795
1796/* O(N) simplification for tolerance = 0 */
1797static void
1799{
1800 uint32_t kept_it = 0;
1801 uint32_t last_it = pa->npoints - 1;
1802 const POINT2D *kept_pt = getPoint2d_cp(pa, 0);
1803 const size_t pt_size = ptarray_point_size(pa);
1804
1805 for (uint32_t i = 1; i < last_it; i++)
1806 {
1807 const POINT2D *curr_pt = getPoint2d_cp(pa, i);
1808 const POINT2D *next_pt = getPoint2d_cp(pa, i + 1);
1809
1810 double ba_x = next_pt->x - kept_pt->x;
1811 double ba_y = next_pt->y - kept_pt->y;
1812 double ab_length_sqr = ba_x * ba_x + ba_y * ba_y;
1813
1814 double ca_x = curr_pt->x - kept_pt->x;
1815 double ca_y = curr_pt->y - kept_pt->y;
1816 double dot_ac_ab = ca_x * ba_x + ca_y * ba_y;
1817 double s_numerator = ca_x * ba_y - ca_y * ba_x;
1818
1819 if (p2d_same(kept_pt, next_pt) ||
1820 dot_ac_ab < 0.0 ||
1821 dot_ac_ab > ab_length_sqr ||
1822 s_numerator != 0)
1823
1824 {
1825 kept_it++;
1826 kept_pt = curr_pt;
1827 if (kept_it != i)
1828 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1829 pa->serialized_pointlist + pt_size * i,
1830 pt_size);
1831 }
1832 }
1833
1834 /* Append last point */
1835 kept_it++;
1836 if (kept_it != last_it)
1837 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1838 pa->serialized_pointlist + pt_size * last_it,
1839 pt_size);
1840 pa->npoints = kept_it + 1;
1841}
1842
1843void
1844ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
1845{
1846 /* Do not try to simplify really short things */
1847 if (pa->npoints < 3 || pa->npoints <= minpts)
1848 return;
1849
1850 if (tolerance == 0 && minpts <= 2)
1851 {
1853 return;
1854 }
1855
1856 /* We use this array to keep track of the points we are keeping, so
1857 * we store just TRUE / FALSE in their position */
1858 uint8_t *kept_points = lwalloc(sizeof(uint8_t) * pa->npoints);
1859 memset(kept_points, LW_FALSE, sizeof(uint8_t) * pa->npoints);
1860 kept_points[0] = LW_TRUE;
1861 kept_points[pa->npoints - 1] = LW_TRUE;
1862 uint32_t keptn = 2;
1863
1864 /* We use this array as a stack to store the iterators that we are going to need
1865 * in the following steps.
1866 * This is ~10% faster than iterating over @kept_points looking for them
1867 */
1868 uint32_t *iterator_stack = lwalloc(sizeof(uint32_t) * pa->npoints);
1869 iterator_stack[0] = 0;
1870 uint32_t iterator_stack_size = 1;
1871
1872 uint32_t it_first = 0;
1873 uint32_t it_last = pa->npoints - 1;
1874
1875 const double tolerance_sqr = tolerance * tolerance;
1876 /* For the first @minpts points we ignore the tolerance */
1877 double it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
1878
1879 while (iterator_stack_size)
1880 {
1881 uint32_t split = ptarray_dp_findsplit_in_place(pa, it_first, it_last, it_tol);
1882 if (split == it_first)
1883 {
1884 it_first = it_last;
1885 it_last = iterator_stack[--iterator_stack_size];
1886 }
1887 else
1888 {
1889 kept_points[split] = LW_TRUE;
1890 keptn++;
1891
1892 iterator_stack[iterator_stack_size++] = it_last;
1893 it_last = split;
1894 it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
1895 }
1896 }
1897
1898 const size_t pt_size = ptarray_point_size(pa);
1899 /* The first point is already in place, so we don't need to copy it */
1900 size_t kept_it = 1;
1901 if (keptn == 2)
1902 {
1903 /* If there are 2 points remaining, it has to be first and last as
1904 * we added those at the start */
1905 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1906 pa->serialized_pointlist + pt_size * (pa->npoints - 1),
1907 pt_size);
1908 }
1909 else if (pa->npoints != keptn) /* We don't need to move any points if we are keeping them all */
1910 {
1911 for (uint32_t i = 1; i < pa->npoints; i++)
1912 {
1913 if (kept_points[i])
1914 {
1915 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1916 pa->serialized_pointlist + pt_size * i,
1917 pt_size);
1918 kept_it++;
1919 }
1920 }
1921 }
1922 pa->npoints = keptn;
1923
1924 lwfree(kept_points);
1925 lwfree(iterator_stack);
1926}
1927
1928/************************************************************************/
1929
1935double
1937{
1938 double dist = 0.0;
1939 uint32_t i;
1940 const POINT2D *a1;
1941 const POINT2D *a2;
1942 const POINT2D *a3;
1943
1944 if ( pts->npoints % 2 != 1 )
1945 lwerror("arc point array with even number of points");
1946
1947 a1 = getPoint2d_cp(pts, 0);
1948
1949 for ( i=2; i < pts->npoints; i += 2 )
1950 {
1951 a2 = getPoint2d_cp(pts, i-1);
1952 a3 = getPoint2d_cp(pts, i);
1953 dist += lw_arc_length(a1, a2, a3);
1954 a1 = a3;
1955 }
1956 return dist;
1957}
1958
1962double
1964{
1965 double dist = 0.0;
1966 uint32_t i;
1967 const POINT2D *frm;
1968 const POINT2D *to;
1969
1970 if ( pts->npoints < 2 ) return 0.0;
1971
1972 frm = getPoint2d_cp(pts, 0);
1973
1974 for ( i=1; i < pts->npoints; i++ )
1975 {
1976 to = getPoint2d_cp(pts, i);
1977
1978 dist += sqrt( ((frm->x - to->x)*(frm->x - to->x)) +
1979 ((frm->y - to->y)*(frm->y - to->y)) );
1980
1981 frm = to;
1982 }
1983 return dist;
1984}
1985
1990double
1992{
1993 double dist = 0.0;
1994 uint32_t i;
1995 POINT3DZ frm;
1996 POINT3DZ to;
1997
1998 if ( pts->npoints < 2 ) return 0.0;
1999
2000 /* compute 2d length if 3d is not available */
2001 if ( ! FLAGS_GET_Z(pts->flags) ) return ptarray_length_2d(pts);
2002
2003 getPoint3dz_p(pts, 0, &frm);
2004 for ( i=1; i < pts->npoints; i++ )
2005 {
2006 getPoint3dz_p(pts, i, &to);
2007 dist += sqrt( ((frm.x - to.x)*(frm.x - to.x)) +
2008 ((frm.y - to.y)*(frm.y - to.y)) +
2009 ((frm.z - to.z)*(frm.z - to.z)) );
2010 frm = to;
2011 }
2012 return dist;
2013}
2014
2015
2016
2020void
2022{
2023 if (FLAGS_GET_Z(pa->flags))
2024 {
2025 for (uint32_t i = 0; i < pa->npoints; i++)
2026 {
2027 POINT4D *p4d = (POINT4D *)(getPoint_internal(pa, i));
2028 double x = p4d->x;
2029 double y = p4d->y;
2030 double z = p4d->z;
2031 p4d->x = a->afac * x + a->bfac * y + a->cfac * z + a->xoff;
2032 p4d->y = a->dfac * x + a->efac * y + a->ffac * z + a->yoff;
2033 p4d->z = a->gfac * x + a->hfac * y + a->ifac * z + a->zoff;
2034 }
2035 }
2036 else
2037 {
2038 for (uint32_t i = 0; i < pa->npoints; i++)
2039 {
2040 POINT2D *pt = (POINT2D *)(getPoint_internal(pa, i));
2041 double x = pt->x;
2042 double y = pt->y;
2043 pt->x = a->afac * x + a->bfac * y + a->xoff;
2044 pt->y = a->dfac * x + a->efac * y + a->yoff;
2045 }
2046 }
2047}
2048
2053#if 0
2054static int gluInvertMatrix(const double *m, double *invOut)
2055{
2056 double inv[16], det;
2057 int i;
2058
2059 inv[0] = m[5] * m[10] * m[15] -
2060 m[5] * m[11] * m[14] -
2061 m[9] * m[6] * m[15] +
2062 m[9] * m[7] * m[14] +
2063 m[13] * m[6] * m[11] -
2064 m[13] * m[7] * m[10];
2065
2066 inv[4] = -m[4] * m[10] * m[15] +
2067 m[4] * m[11] * m[14] +
2068 m[8] * m[6] * m[15] -
2069 m[8] * m[7] * m[14] -
2070 m[12] * m[6] * m[11] +
2071 m[12] * m[7] * m[10];
2072
2073 inv[8] = m[4] * m[9] * m[15] -
2074 m[4] * m[11] * m[13] -
2075 m[8] * m[5] * m[15] +
2076 m[8] * m[7] * m[13] +
2077 m[12] * m[5] * m[11] -
2078 m[12] * m[7] * m[9];
2079
2080 inv[12] = -m[4] * m[9] * m[14] +
2081 m[4] * m[10] * m[13] +
2082 m[8] * m[5] * m[14] -
2083 m[8] * m[6] * m[13] -
2084 m[12] * m[5] * m[10] +
2085 m[12] * m[6] * m[9];
2086
2087 inv[1] = -m[1] * m[10] * m[15] +
2088 m[1] * m[11] * m[14] +
2089 m[9] * m[2] * m[15] -
2090 m[9] * m[3] * m[14] -
2091 m[13] * m[2] * m[11] +
2092 m[13] * m[3] * m[10];
2093
2094 inv[5] = m[0] * m[10] * m[15] -
2095 m[0] * m[11] * m[14] -
2096 m[8] * m[2] * m[15] +
2097 m[8] * m[3] * m[14] +
2098 m[12] * m[2] * m[11] -
2099 m[12] * m[3] * m[10];
2100
2101 inv[9] = -m[0] * m[9] * m[15] +
2102 m[0] * m[11] * m[13] +
2103 m[8] * m[1] * m[15] -
2104 m[8] * m[3] * m[13] -
2105 m[12] * m[1] * m[11] +
2106 m[12] * m[3] * m[9];
2107
2108 inv[13] = m[0] * m[9] * m[14] -
2109 m[0] * m[10] * m[13] -
2110 m[8] * m[1] * m[14] +
2111 m[8] * m[2] * m[13] +
2112 m[12] * m[1] * m[10] -
2113 m[12] * m[2] * m[9];
2114
2115 inv[2] = m[1] * m[6] * m[15] -
2116 m[1] * m[7] * m[14] -
2117 m[5] * m[2] * m[15] +
2118 m[5] * m[3] * m[14] +
2119 m[13] * m[2] * m[7] -
2120 m[13] * m[3] * m[6];
2121
2122 inv[6] = -m[0] * m[6] * m[15] +
2123 m[0] * m[7] * m[14] +
2124 m[4] * m[2] * m[15] -
2125 m[4] * m[3] * m[14] -
2126 m[12] * m[2] * m[7] +
2127 m[12] * m[3] * m[6];
2128
2129 inv[10] = m[0] * m[5] * m[15] -
2130 m[0] * m[7] * m[13] -
2131 m[4] * m[1] * m[15] +
2132 m[4] * m[3] * m[13] +
2133 m[12] * m[1] * m[7] -
2134 m[12] * m[3] * m[5];
2135
2136 inv[14] = -m[0] * m[5] * m[14] +
2137 m[0] * m[6] * m[13] +
2138 m[4] * m[1] * m[14] -
2139 m[4] * m[2] * m[13] -
2140 m[12] * m[1] * m[6] +
2141 m[12] * m[2] * m[5];
2142
2143 inv[3] = -m[1] * m[6] * m[11] +
2144 m[1] * m[7] * m[10] +
2145 m[5] * m[2] * m[11] -
2146 m[5] * m[3] * m[10] -
2147 m[9] * m[2] * m[7] +
2148 m[9] * m[3] * m[6];
2149
2150 inv[7] = m[0] * m[6] * m[11] -
2151 m[0] * m[7] * m[10] -
2152 m[4] * m[2] * m[11] +
2153 m[4] * m[3] * m[10] +
2154 m[8] * m[2] * m[7] -
2155 m[8] * m[3] * m[6];
2156
2157 inv[11] = -m[0] * m[5] * m[11] +
2158 m[0] * m[7] * m[9] +
2159 m[4] * m[1] * m[11] -
2160 m[4] * m[3] * m[9] -
2161 m[8] * m[1] * m[7] +
2162 m[8] * m[3] * m[5];
2163
2164 inv[15] = m[0] * m[5] * m[10] -
2165 m[0] * m[6] * m[9] -
2166 m[4] * m[1] * m[10] +
2167 m[4] * m[2] * m[9] +
2168 m[8] * m[1] * m[6] -
2169 m[8] * m[2] * m[5];
2170
2171 det = m[0] * inv[0] + m[1] * inv[4] + m[2] * inv[8] + m[3] * inv[12];
2172
2173 if (det == 0)
2174 return LW_FALSE;
2175
2176 det = 1.0 / det;
2177
2178 for (i = 0; i < 16; i++)
2179 invOut[i] = inv[i] * det;
2180
2181 return LW_TRUE;
2182}
2183#endif
2184
2188void
2190{
2191 uint32_t i;
2192 POINT4D p4d;
2193 LWDEBUG(3, "ptarray_scale start");
2194 for (i=0; i<pa->npoints; i++)
2195 {
2196 getPoint4d_p(pa, i, &p4d);
2197 p4d.x *= fact->x;
2198 p4d.y *= fact->y;
2199 p4d.z *= fact->z;
2200 p4d.m *= fact->m;
2201 ptarray_set_point4d(pa, i, &p4d);
2202 }
2203 LWDEBUG(3, "ptarray_scale end");
2204}
2205
2206int
2208{
2209 return getPoint4d_p(pa, 0, pt);
2210}
2211
2212
2213/*
2214 * Stick an array of points to the given gridspec.
2215 * Return "gridded" points in *outpts and their number in *outptsn.
2216 *
2217 * Two consecutive points falling on the same grid cell are collapsed
2218 * into one single point.
2219 *
2220 */
2221void
2223{
2224 uint32_t j = 0;
2225 POINT4D *p, *p_out = NULL;
2226 double x, y, z = 0, m = 0;
2227 uint32_t ndims = FLAGS_NDIMS(pa->flags);
2228 uint32_t has_z = FLAGS_GET_Z(pa->flags);
2229 uint32_t has_m = FLAGS_GET_M(pa->flags);
2230
2231 for (uint32_t i = 0; i < pa->npoints; i++)
2232 {
2233 /* Look straight into the abyss */
2234 p = (POINT4D *)(getPoint_internal(pa, i));
2235 x = p->x;
2236 y = p->y;
2237 if (ndims > 2)
2238 z = p->z;
2239 if (ndims > 3)
2240 m = p->m;
2241
2242 /*
2243 * See https://github.com/libgeos/geos/pull/956
2244 * We use scale for rounding when gridsize is < 1 and
2245 * gridsize for rounding when scale < 1.
2246 */
2247 if (grid->xsize > 0) {
2248 if (grid->xsize < 1)
2249 x = rint((x - grid->ipx) * grid->xscale) / grid->xscale + grid->ipx;
2250 else
2251 x = rint((x - grid->ipx) / grid->xsize) * grid->xsize + grid->ipx;
2252 }
2253
2254 if (grid->ysize > 0) {
2255 if (grid->ysize < 1)
2256 y = rint((y - grid->ipy) * grid->yscale) / grid->yscale + grid->ipy;
2257 else
2258 y = rint((y - grid->ipy) / grid->ysize) * grid->ysize + grid->ipy;
2259 }
2260
2261 /* Read and round this point */
2262 /* Z is always in third position */
2263 if (has_z && grid->zsize > 0)
2264 z = rint((z - grid->ipz) / grid->zsize) * grid->zsize + grid->ipz;
2265
2266 /* M might be in 3rd or 4th position */
2267 if (has_m && grid->msize > 0)
2268 {
2269 /* In POINT ZM, M is in 4th position, in POINT M, M is in 3rd position which is Z in POINT4D */
2270 if (has_z)
2271 m = rint((m - grid->ipm) / grid->msize) * grid->msize + grid->ipm;
2272 else
2273 z = rint((z - grid->ipm) / grid->msize) * grid->msize + grid->ipm;
2274 }
2275
2276 /* Skip duplicates */
2277 if (p_out &&
2278 p_out->x == x &&
2279 p_out->y == y &&
2280 (ndims > 2 ? p_out->z == z : 1) &&
2281 (ndims > 3 ? p_out->m == m : 1))
2282 {
2283 continue;
2284 }
2285
2286 /* Write rounded values into the next available point */
2287 p_out = (POINT4D *)(getPoint_internal(pa, j++));
2288 p_out->x = x;
2289 p_out->y = y;
2290 if (ndims > 2)
2291 p_out->z = z;
2292 if (ndims > 3)
2293 p_out->m = m;
2294 }
2295
2296 /* Update output ptarray length */
2297 pa->npoints = j;
2298 return;
2299}
2300
2301int
2303{
2304 const POINT2D *pt;
2305 int n = 0;
2306 uint32_t i;
2307 for ( i = 0; i < pa->npoints; i++ )
2308 {
2309 pt = getPoint2d_cp(pa, i);
2310 if ( gbox_contains_point2d(gbox, pt) )
2311 n++;
2312 }
2313 return n;
2314}
2315
2316
2317/*
2318 * Reorder the vertices of a closed pointarray so that the
2319 * given point is the first/last one.
2320 *
2321 * Error out if pointarray is not closed or it does not
2322 * contain the given point.
2323 */
2324int
2326{
2327 POINTARRAY *tmp;
2328 int found;
2329 uint32_t it;
2330 int ptsize;
2331
2332 if ( ! ptarray_is_closed_2d(pa) )
2333 {
2334 lwerror("ptarray_scroll_in_place: input POINTARRAY is not closed");
2335 return LW_FAILURE;
2336 }
2337
2338 ptsize = ptarray_point_size(pa);
2339
2340 /* Find the point in the array */
2341 found = 0;
2342 for ( it = 0; it < pa->npoints; ++it )
2343 {
2344 if ( ! memcmp(getPoint_internal(pa, it), pt, ptsize) )
2345 {
2346 found = 1;
2347 break;
2348 }
2349 }
2350
2351 if ( ! found )
2352 {
2353 lwerror("ptarray_scroll_in_place: input POINTARRAY does not contain the given point");
2354 return LW_FAILURE;
2355 }
2356
2357 if ( 0 == it )
2358 {
2359 /* Point is already the start/end point, just clone the input */
2360 return LW_SUCCESS;
2361 }
2362
2363 /* TODO: reduce allocations */
2365
2366 memset(getPoint_internal(tmp, 0), 0, (size_t)ptsize * pa->npoints);
2367 /* Copy the block from found point to last point into the output array */
2368 memcpy(
2369 getPoint_internal(tmp, 0),
2370 getPoint_internal(pa, it),
2371 (size_t)ptsize * ( pa->npoints - it )
2372 );
2373
2374 /* Copy the block from second point to the found point into the last portion of the
2375 * return */
2376 memcpy(
2377 getPoint_internal(tmp, pa->npoints - it),
2378 getPoint_internal(pa, 1),
2379 (size_t)ptsize * ( it )
2380 );
2381
2382 /* Copy the resulting pointarray back to source one */
2383 memcpy(
2384 getPoint_internal(pa, 0),
2385 getPoint_internal(tmp, 0),
2386 (size_t)ptsize * ( pa->npoints )
2387 );
2388
2389 ptarray_free(tmp);
2390
2391 return LW_SUCCESS;
2392}
char * r
Definition cu_in_wkt.c:24
int gbox_contains_point2d(const GBOX *g, const POINT2D *p)
Definition gbox.c:362
int lw_arc_calculate_gbox_cartesian_2d(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, GBOX *gbox)
Definition gbox.c:465
#define LW_FALSE
Definition liblwgeom.h:94
void * lwrealloc(void *mem, size_t size)
Definition lwutil.c:242
double distance2d_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition measures.c:2344
#define LW_FAILURE
Definition liblwgeom.h:96
void interpolate_point4d(const POINT4D *A, const POINT4D *B, POINT4D *I, double F)
Find interpolation point I between point A and point B so that the len(AI) == len(AB)*F and I falls o...
Definition lwgeom_api.c:646
#define LW_SUCCESS
Definition liblwgeom.h:97
int getPoint2d_p(const POINTARRAY *pa, uint32_t n, POINT2D *point)
Definition lwgeom_api.c:342
#define FLAGS_GET_READONLY(flags)
Definition liblwgeom.h:169
#define FLAGS_GET_Z(flags)
Definition liblwgeom.h:165
void * lwalloc(size_t size)
Definition lwutil.c:227
int getPoint3dz_p(const POINTARRAY *pa, uint32_t n, POINT3DZ *point)
Definition lwgeom_api.c:215
void lwfree(void *mem)
Definition lwutil.c:248
#define FLAGS_NDIMS(flags)
Definition liblwgeom.h:179
#define FLAGS_GET_M(flags)
Definition liblwgeom.h:166
int getPoint4d_p(const POINTARRAY *pa, uint32_t n, POINT4D *point)
Definition lwgeom_api.c:125
double distance2d_sqr_pt_seg(const POINT2D *p, const POINT2D *A, const POINT2D *B)
Definition measures.c:2354
#define FLAGS_GET_ZM(flags)
Definition liblwgeom.h:180
#define FLAGS_SET_READONLY(flags, value)
Definition liblwgeom.h:176
lwflags_t lwflags(int hasz, int hasm, int geodetic)
Construct a new flags bitmask.
Definition lwutil.c:477
#define LW_TRUE
Return types for functions with status returns.
Definition liblwgeom.h:93
void ptarray_set_point4d(POINTARRAY *pa, uint32_t n, const POINT4D *p4d)
Definition lwgeom_api.c:369
enum LWORD_T LWORD
Ordinate names.
#define LW_INSIDE
Constants for point-in-polygon return values.
#define LW_BOUNDARY
double lw_arc_length(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3)
Returns the length of a circular arc segment.
#define LW_ON_INTERRUPT(x)
#define FP_MAX(A, B)
#define FP_MIN(A, B)
int lw_pt_on_segment(const POINT2D *p1, const POINT2D *p2, const POINT2D *p)
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 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 lw_segment_side(const POINT2D *p1, const POINT2D *p2, const POINT2D *q)
lw_segment_side()
Definition lwalgorithm.c:70
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
int p2d_same(const POINT2D *p1, const POINT2D *p2)
Definition lwalgorithm.c:57
static double det(double m00, double m01, double m10, double m11)
#define LWDEBUG(level, msg)
Definition lwgeom_log.h:101
#define LWDEBUGF(level, msg,...)
Definition lwgeom_log.h:106
void lwnotice(const char *fmt,...) __attribute__((format(printf
Write a notice out to the notice handler.
void void lwerror(const char *fmt,...) __attribute__((format(printf
Write a notice out to the error handler.
#define INT32_MAX
static uint8_t * getPoint_internal(const POINTARRAY *pa, uint32_t n)
Definition lwinline.h:75
static double distance2d_sqr_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition lwinline.h:33
static size_t ptarray_point_size(const POINTARRAY *pa)
Definition lwinline.h:56
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
int ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
Definition ptarray.c:1031
static uint32_t ptarray_dp_findsplit_in_place(const POINTARRAY *pts, uint32_t it_first, uint32_t it_last, double max_distance_sqr)
Definition ptarray.c:1733
int ptarray_remove_point(POINTARRAY *pa, uint32_t where)
Remove a point from an existing POINTARRAY.
Definition ptarray.c:259
void ptarray_longitude_shift(POINTARRAY *pa)
Longitude shift for a pointarray.
Definition ptarray.c:1625
int ptarray_closest_segment_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
Definition ptarray.c:1447
POINTARRAY * ptarray_clone_deep(const POINTARRAY *in)
Deep clone a pointarray (also clones serialized pointlist)
Definition ptarray.c:643
POINTARRAY * ptarray_clone(const POINTARRAY *in)
Clone a POINTARRAY object.
Definition ptarray.c:674
int ptarray_append_ptarray(POINTARRAY *pa1, POINTARRAY *pa2, double gap_tolerance)
Append a POINTARRAY, pa2 to the end of an existing POINTARRAY, pa1.
Definition ptarray.c:177
void ptarray_reverse_in_place(POINTARRAY *pa)
Definition ptarray.c:339
int ptarray_closest_vertex_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
Definition ptarray.c:1480
int ptarrayarc_contains_point(const POINTARRAY *pa, const POINT2D *pt)
For POINTARRAYs representing CIRCULARSTRINGS.
Definition ptarray.c:1123
int ptarray_is_closed_2d(const POINTARRAY *in)
Definition ptarray.c:710
POINTARRAY * ptarray_substring(POINTARRAY *ipa, double from, double to, double tolerance)
@d1 start location (distance from start / total distance) @d2 end location (distance from start / tot...
Definition ptarray.c:1206
POINTARRAY * ptarray_construct_reference_data(char hasz, char hasm, uint32_t npoints, uint8_t *ptlist)
Build a new POINTARRAY, but on top of someone else's ordinate array.
Definition ptarray.c:291
int ptarray_is_closed_z(const POINTARRAY *in)
Definition ptarray.c:736
int ptarray_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
Definition ptarray.c:827
double ptarray_length(const POINTARRAY *pts)
Find the 3d/2d length of the given POINTARRAY (depending on its dimensionality)
Definition ptarray.c:1991
double ptarray_signed_area(const POINTARRAY *pa)
Returns the area in cartesian units.
Definition ptarray.c:1143
static POINTARRAY * ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
Definition ptarray.c:1647
void closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
Definition ptarray.c:1403
int ptarray_startpoint(const POINTARRAY *pa, POINT4D *pt)
Definition ptarray.c:2207
POINTARRAY * ptarray_construct_empty(char hasz, char hasm, uint32_t maxpoints)
Create a new POINTARRAY with no points.
Definition ptarray.c:59
int ptarray_is_closed(const POINTARRAY *in)
Check for ring closure using whatever dimensionality is declared on the pointarray.
Definition ptarray.c:696
double ptarray_length_2d(const POINTARRAY *pts)
Find the 2d length of the given POINTARRAY (even if it's 3d)
Definition ptarray.c:1963
char ptarray_same(const POINTARRAY *pa1, const POINTARRAY *pa2)
Definition ptarray.c:484
int ptarray_is_closed_3d(const POINTARRAY *in)
Definition ptarray.c:723
void ptarray_remove_repeated_points_in_place(POINTARRAY *pa, double tolerance, uint32_t min_points)
Definition ptarray.c:1662
int ptarrayarc_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
Definition ptarray.c:944
int ptarray_scroll_in_place(POINTARRAY *pa, const POINT4D *pt)
Definition ptarray.c:2325
void ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
Definition ptarray.c:1844
POINTARRAY * ptarray_flip_coordinates(POINTARRAY *pa)
Reverse X and Y axis on a given POINTARRAY.
Definition ptarray.c:368
int ptarray_insert_point(POINTARRAY *pa, const POINT4D *p, uint32_t where)
Insert a point into an existing POINTARRAY.
Definition ptarray.c:85
void ptarray_grid_in_place(POINTARRAY *pa, gridspec *grid)
Snap to grid.
Definition ptarray.c:2222
double ptarray_arc_length_2d(const POINTARRAY *pts)
Find the 2d length of the given POINTARRAY, using circular arc interpolation between each coordinate ...
Definition ptarray.c:1936
static int circle_raycast_intersections(const POINT2D *center, double radius, double ray, POINT2D *i0, POINT2D *i1)
Calculates the intersection points of a circle and a horizontal line.
Definition ptarray.c:903
POINTARRAY * ptarray_merge(POINTARRAY *pa1, POINTARRAY *pa2)
Merge two given POINTARRAY and returns a pointer on the new aggregate one.
Definition ptarray.c:612
void ptarray_affine(POINTARRAY *pa, const AFFINE *a)
Affine transform a pointarray.
Definition ptarray.c:2021
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
POINTARRAY * ptarray_segmentize2d(const POINTARRAY *ipa, double dist)
Returns a modified POINTARRAY so that no segment is longer than the given distance (computed using 2d...
Definition ptarray.c:413
POINTARRAY * ptarray_removePoint(POINTARRAY *pa, uint32_t which)
Remove a point from a pointarray.
Definition ptarray.c:582
void ptarray_free(POINTARRAY *pa)
Definition ptarray.c:327
static void ptarray_simplify_in_place_tolerance0(POINTARRAY *pa)
Definition ptarray.c:1798
int ptarray_isccw(const POINTARRAY *pa)
Definition ptarray.c:1174
int ptarray_has_z(const POINTARRAY *pa)
Definition ptarray.c:37
int ptarray_append_point(POINTARRAY *pa, const POINT4D *pt, int repeated_points)
Append a point to the end of an existing POINTARRAY If allow_duplicate is LW_FALSE,...
Definition ptarray.c:147
POINTARRAY * ptarray_addPoint(const POINTARRAY *pa, uint8_t *p, size_t pdims, uint32_t where)
Add a point in a pointarray.
Definition ptarray.c:530
void ptarray_scale(POINTARRAY *pa, const POINT4D *fact)
WARNING, make sure you send in only 16-member double arrays or obviously things will go pear-shaped f...
Definition ptarray.c:2189
void ptarray_swap_ordinates(POINTARRAY *pa, LWORD o1, LWORD o2)
Swap ordinate values o1 and o2 on a given POINTARRAY.
Definition ptarray.c:387
POINTARRAY * ptarray_construct_copy_data(char hasz, char hasm, uint32_t npoints, const uint8_t *ptlist)
Construct a new POINTARRAY, copying in the data from ptlist.
Definition ptarray.c:305
int ptarray_npoints_in_rect(const POINTARRAY *pa, const GBOX *gbox)
Definition ptarray.c:2302
char ptarray_same2d(const POINTARRAY *pa1, const POINTARRAY *pa2)
Definition ptarray.c:509
POINTARRAY * ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm, double zval, double mval)
Definition ptarray.c:1183
POINTARRAY * ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
Definition ptarray.c:1655
int ptarray_has_m(const POINTARRAY *pa)
Definition ptarray.c:44
double ptarray_locate_point(const POINTARRAY *pa, const POINT4D *p4d, double *mindistout, POINT4D *proj4d)
Definition ptarray.c:1513
POINTARRAY * ptarray_construct(char hasz, char hasm, uint32_t npoints)
Construct an empty pointarray, allocating storage and setting the npoints, but not filling in any inf...
Definition ptarray.c:51
double gfac
Definition liblwgeom.h:332
double zoff
Definition liblwgeom.h:332
double bfac
Definition liblwgeom.h:332
double ifac
Definition liblwgeom.h:332
double xoff
Definition liblwgeom.h:332
double dfac
Definition liblwgeom.h:332
double afac
Definition liblwgeom.h:332
double ffac
Definition liblwgeom.h:332
double cfac
Definition liblwgeom.h:332
double hfac
Definition liblwgeom.h:332
double efac
Definition liblwgeom.h:332
double yoff
Definition liblwgeom.h:332
double ymax
Definition liblwgeom.h:357
double ymin
Definition liblwgeom.h:356
double y
Definition liblwgeom.h:390
double x
Definition liblwgeom.h:390
double z
Definition liblwgeom.h:396
double x
Definition liblwgeom.h:396
double y
Definition liblwgeom.h:396
double m
Definition liblwgeom.h:414
double x
Definition liblwgeom.h:414
double z
Definition liblwgeom.h:414
double y
Definition liblwgeom.h:414
lwflags_t flags
Definition liblwgeom.h:431
uint32_t maxpoints
Definition liblwgeom.h:428
uint32_t npoints
Definition liblwgeom.h:427
uint8_t * serialized_pointlist
Definition liblwgeom.h:434
double xscale
Definition liblwgeom.h:1407
double ipm
Definition liblwgeom.h:1402
double zsize
Definition liblwgeom.h:1405
double ysize
Definition liblwgeom.h:1404
double xsize
Definition liblwgeom.h:1403
double yscale
Definition liblwgeom.h:1408
double ipx
Definition liblwgeom.h:1399
double msize
Definition liblwgeom.h:1406
double ipy
Definition liblwgeom.h:1400
double ipz
Definition liblwgeom.h:1401
Snap-to-grid.
Definition liblwgeom.h:1398