PostGIS 3.0.6dev-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 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 %d 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 %d", point_size, where);
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), 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 = %d", 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
509ptarray_addPoint(const POINTARRAY *pa, uint8_t *p, size_t pdims, uint32_t where)
510{
511 POINTARRAY *ret;
512 POINT4D pbuf;
513 size_t ptsize = ptarray_point_size(pa);
514
515 LWDEBUGF(3, "pa %x p %x size %d where %d",
516 pa, p, pdims, where);
517
518 if ( pdims < 2 || pdims > 4 )
519 {
520 lwerror("ptarray_addPoint: point dimension out of range (%d)",
521 pdims);
522 return NULL;
523 }
524
525 if ( where > pa->npoints )
526 {
527 lwerror("ptarray_addPoint: offset out of range (%d)",
528 where);
529 return NULL;
530 }
531
532 LWDEBUG(3, "called with a %dD point");
533
534 pbuf.x = pbuf.y = pbuf.z = pbuf.m = 0.0;
535 memcpy((uint8_t *)&pbuf, p, pdims*sizeof(double));
536
537 LWDEBUG(3, "initialized point buffer");
538
540 FLAGS_GET_M(pa->flags), pa->npoints+1);
541
542
543 if ( where )
544 {
545 memcpy(getPoint_internal(ret, 0), getPoint_internal(pa, 0), ptsize*where);
546 }
547
548 memcpy(getPoint_internal(ret, where), (uint8_t *)&pbuf, ptsize);
549
550 if ( where+1 != ret->npoints )
551 {
552 memcpy(getPoint_internal(ret, where+1),
553 getPoint_internal(pa, where),
554 ptsize*(pa->npoints-where));
555 }
556
557 return ret;
558}
559
561ptarray_removePoint(POINTARRAY *pa, uint32_t which)
562{
563 POINTARRAY *ret;
564 size_t ptsize = ptarray_point_size(pa);
565
566 LWDEBUGF(3, "pa %x which %d", pa, which);
567
568#if PARANOIA_LEVEL > 0
569 if ( which > pa->npoints-1 )
570 {
571 lwerror("%s [%d] offset (%d) out of range (%d..%d)", __FILE__, __LINE__,
572 which, 0, pa->npoints-1);
573 return NULL;
574 }
575
576 if ( pa->npoints < 3 )
577 {
578 lwerror("%s [%d] can't remove a point from a 2-vertex POINTARRAY", __FILE__, __LINE__);
579 return NULL;
580 }
581#endif
582
584 FLAGS_GET_M(pa->flags), pa->npoints-1);
585
586 /* copy initial part */
587 if ( which )
588 {
589 memcpy(getPoint_internal(ret, 0), getPoint_internal(pa, 0), ptsize*which);
590 }
591
592 /* copy final part */
593 if ( which < pa->npoints-1 )
594 {
595 memcpy(getPoint_internal(ret, which), getPoint_internal(pa, which+1),
596 ptsize*(pa->npoints-which-1));
597 }
598
599 return ret;
600}
601
604{
605 POINTARRAY *pa;
606 size_t ptsize = ptarray_point_size(pa1);
607
608 if (FLAGS_GET_ZM(pa1->flags) != FLAGS_GET_ZM(pa2->flags))
609 lwerror("ptarray_cat: Mixed dimension");
610
612 FLAGS_GET_M(pa1->flags),
613 pa1->npoints + pa2->npoints);
614
615 memcpy( getPoint_internal(pa, 0),
616 getPoint_internal(pa1, 0),
617 ptsize*(pa1->npoints));
618
619 memcpy( getPoint_internal(pa, pa1->npoints),
620 getPoint_internal(pa2, 0),
621 ptsize*(pa2->npoints));
622
623 ptarray_free(pa1);
624 ptarray_free(pa2);
625
626 return pa;
627}
628
629
635{
636 POINTARRAY *out = lwalloc(sizeof(POINTARRAY));
637
638 LWDEBUG(3, "ptarray_clone_deep called.");
639
640 out->flags = in->flags;
641 out->npoints = in->npoints;
642 out->maxpoints = in->npoints;
643
644 FLAGS_SET_READONLY(out->flags, 0);
645
646 if (!in->npoints)
647 {
648 // Avoid calling lwalloc of 0 bytes
649 out->serialized_pointlist = NULL;
650 }
651 else
652 {
653 size_t size = in->npoints * ptarray_point_size(in);
654 out->serialized_pointlist = lwalloc(size);
655 memcpy(out->serialized_pointlist, in->serialized_pointlist, size);
656 }
657
658 return out;
659}
660
666{
667 POINTARRAY *out = lwalloc(sizeof(POINTARRAY));
668
669 LWDEBUG(3, "ptarray_clone called.");
670
671 out->flags = in->flags;
672 out->npoints = in->npoints;
673 out->maxpoints = in->maxpoints;
674
675 FLAGS_SET_READONLY(out->flags, 1);
676
678
679 return out;
680}
681
686int
688{
689 if (!in)
690 {
691 lwerror("ptarray_is_closed: called with null point array");
692 return 0;
693 }
694 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
695
696 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), ptarray_point_size(in));
697}
698
699
700int
702{
703 if (!in)
704 {
705 lwerror("ptarray_is_closed_2d: called with null point array");
706 return 0;
707 }
708 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
709
710 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), sizeof(POINT2D) );
711}
712
713int
715{
716 if (!in)
717 {
718 lwerror("ptarray_is_closed_3d: called with null point array");
719 return 0;
720 }
721 if (in->npoints <= 1 ) return in->npoints; /* single-point are closed, empty not closed */
722
723 return 0 == memcmp(getPoint_internal(in, 0), getPoint_internal(in, in->npoints-1), sizeof(POINT3D) );
724}
725
726int
728{
729 if ( FLAGS_GET_Z(in->flags) )
730 return ptarray_is_closed_3d(in);
731 else
732 return ptarray_is_closed_2d(in);
733}
734
739int
741{
742 return ptarray_contains_point_partial(pa, pt, LW_TRUE, NULL);
743}
744
745int
746ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
747{
748 int wn = 0;
749 uint32_t i;
750 double side;
751 const POINT2D *seg1;
752 const POINT2D *seg2;
753 double ymin, ymax;
754
755 seg1 = getPoint2d_cp(pa, 0);
756 seg2 = getPoint2d_cp(pa, pa->npoints-1);
757 if ( check_closed && ! p2d_same(seg1, seg2) )
758 lwerror("ptarray_contains_point called on unclosed ring");
759
760 for ( i=1; i < pa->npoints; i++ )
761 {
762 seg2 = getPoint2d_cp(pa, i);
763
764 /* Zero length segments are ignored. */
765 if ( seg1->x == seg2->x && seg1->y == seg2->y )
766 {
767 seg1 = seg2;
768 continue;
769 }
770
771 ymin = FP_MIN(seg1->y, seg2->y);
772 ymax = FP_MAX(seg1->y, seg2->y);
773
774 /* Only test segments in our vertical range */
775 if ( pt->y > ymax || pt->y < ymin )
776 {
777 seg1 = seg2;
778 continue;
779 }
780
781 side = lw_segment_side(seg1, seg2, pt);
782
783 /*
784 * A point on the boundary of a ring is not contained.
785 * WAS: if (fabs(side) < 1e-12), see #852
786 */
787 if ( (side == 0) && lw_pt_in_seg(pt, seg1, seg2) )
788 {
789 return LW_BOUNDARY;
790 }
791
792 /*
793 * If the point is to the left of the line, and it's rising,
794 * then the line is to the right of the point and
795 * circling counter-clockwise, so increment.
796 */
797 if ( (side < 0) && (seg1->y <= pt->y) && (pt->y < seg2->y) )
798 {
799 wn++;
800 }
801
802 /*
803 * If the point is to the right of the line, and it's falling,
804 * then the line is to the right of the point and circling
805 * clockwise, so decrement.
806 */
807 else if ( (side > 0) && (seg2->y <= pt->y) && (pt->y < seg1->y) )
808 {
809 wn--;
810 }
811
812 seg1 = seg2;
813 }
814
815 /* Sent out the winding number for calls that are building on this as a primitive */
816 if ( winding_number )
817 *winding_number = wn;
818
819 /* Outside */
820 if (wn == 0)
821 {
822 return LW_OUTSIDE;
823 }
824
825 /* Inside */
826 return LW_INSIDE;
827}
828
838int
840{
841 return ptarrayarc_contains_point_partial(pa, pt, LW_TRUE /* Check closed*/, NULL);
842}
843
844int
845ptarrayarc_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
846{
847 int wn = 0;
848 uint32_t i;
849 int side;
850 const POINT2D *seg1;
851 const POINT2D *seg2;
852 const POINT2D *seg3;
853 GBOX gbox;
854
855 /* Check for not an arc ring (always have odd # of points) */
856 if ( (pa->npoints % 2) == 0 )
857 {
858 lwerror("ptarrayarc_contains_point called with even number of points");
859 return LW_OUTSIDE;
860 }
861
862 /* Check for not an arc ring (always have >= 3 points) */
863 if ( pa->npoints < 3 )
864 {
865 lwerror("ptarrayarc_contains_point called too-short pointarray");
866 return LW_OUTSIDE;
867 }
868
869 /* Check for unclosed case */
870 seg1 = getPoint2d_cp(pa, 0);
871 seg3 = getPoint2d_cp(pa, pa->npoints-1);
872 if ( check_closed && ! p2d_same(seg1, seg3) )
873 {
874 lwerror("ptarrayarc_contains_point called on unclosed ring");
875 return LW_OUTSIDE;
876 }
877 /* OK, it's closed. Is it just one circle? */
878 else if ( p2d_same(seg1, seg3) && pa->npoints == 3 )
879 {
880 double radius, d;
881 POINT2D c;
882 seg2 = getPoint2d_cp(pa, 1);
883
884 /* Wait, it's just a point, so it can't contain anything */
885 if ( lw_arc_is_pt(seg1, seg2, seg3) )
886 return LW_OUTSIDE;
887
888 /* See if the point is within the circle radius */
889 radius = lw_arc_center(seg1, seg2, seg3, &c);
890 d = distance2d_pt_pt(pt, &c);
891 if ( FP_EQUALS(d, radius) )
892 return LW_BOUNDARY; /* Boundary of circle */
893 else if ( d < radius )
894 return LW_INSIDE; /* Inside circle */
895 else
896 return LW_OUTSIDE; /* Outside circle */
897 }
898 else if ( p2d_same(seg1, pt) || p2d_same(seg3, pt) )
899 {
900 return LW_BOUNDARY; /* Boundary case */
901 }
902
903 /* Start on the ring */
904 seg1 = getPoint2d_cp(pa, 0);
905 for ( i=1; i < pa->npoints; i += 2 )
906 {
907 seg2 = getPoint2d_cp(pa, i);
908 seg3 = getPoint2d_cp(pa, i+1);
909
910 /* Catch an easy boundary case */
911 if( p2d_same(seg3, pt) )
912 return LW_BOUNDARY;
913
914 /* Skip arcs that have no size */
915 if ( lw_arc_is_pt(seg1, seg2, seg3) )
916 {
917 seg1 = seg3;
918 continue;
919 }
920
921 /* Only test segments in our vertical range */
922 lw_arc_calculate_gbox_cartesian_2d(seg1, seg2, seg3, &gbox);
923 if ( pt->y > gbox.ymax || pt->y < gbox.ymin )
924 {
925 seg1 = seg3;
926 continue;
927 }
928
929 /* Outside of horizontal range, and not between end points we also skip */
930 if ( (pt->x > gbox.xmax || pt->x < gbox.xmin) &&
931 (pt->y > FP_MAX(seg1->y, seg3->y) || pt->y < FP_MIN(seg1->y, seg3->y)) )
932 {
933 seg1 = seg3;
934 continue;
935 }
936
937 side = lw_arc_side(seg1, seg2, seg3, pt);
938
939 /* On the boundary */
940 if ( (side == 0) && lw_pt_in_arc(pt, seg1, seg2, seg3) )
941 {
942 return LW_BOUNDARY;
943 }
944
945 /* Going "up"! Point to left of arc. */
946 if ( side < 0 && (seg1->y <= pt->y) && (pt->y < seg3->y) )
947 {
948 wn++;
949 }
950
951 /* Going "down"! */
952 if ( side > 0 && (seg3->y <= pt->y) && (pt->y < seg1->y) )
953 {
954 wn--;
955 }
956
957 /* Inside the arc! */
958 if ( pt->x <= gbox.xmax && pt->x >= gbox.xmin )
959 {
960 POINT2D C;
961 double radius = lw_arc_center(seg1, seg2, seg3, &C);
962 double d = distance2d_pt_pt(pt, &C);
963
964 /* On the boundary! */
965 if ( d == radius )
966 return LW_BOUNDARY;
967
968 /* Within the arc! */
969 if ( d < radius )
970 {
971 /* Left side, increment winding number */
972 if ( side < 0 )
973 wn++;
974 /* Right side, decrement winding number */
975 if ( side > 0 )
976 wn--;
977 }
978 }
979
980 seg1 = seg3;
981 }
982
983 /* Sent out the winding number for calls that are building on this as a primitive */
984 if ( winding_number )
985 *winding_number = wn;
986
987 /* Outside */
988 if (wn == 0)
989 {
990 return LW_OUTSIDE;
991 }
992
993 /* Inside */
994 return LW_INSIDE;
995}
996
1002double
1004{
1005 const POINT2D *P1;
1006 const POINT2D *P2;
1007 const POINT2D *P3;
1008 double sum = 0.0;
1009 double x0, x, y1, y2;
1010 uint32_t i;
1011
1012 if (! pa || pa->npoints < 3 )
1013 return 0.0;
1014
1015 P1 = getPoint2d_cp(pa, 0);
1016 P2 = getPoint2d_cp(pa, 1);
1017 x0 = P1->x;
1018 for ( i = 2; i < pa->npoints; i++ )
1019 {
1020 P3 = getPoint2d_cp(pa, i);
1021 x = P2->x - x0;
1022 y1 = P3->y;
1023 y2 = P1->y;
1024 sum += x * (y2-y1);
1025
1026 /* Move forwards! */
1027 P1 = P2;
1028 P2 = P3;
1029 }
1030 return sum / 2.0;
1031}
1032
1033int
1035{
1036 double area = 0;
1037 area = ptarray_signed_area(pa);
1038 if ( area > 0 ) return LW_FALSE;
1039 else return LW_TRUE;
1040}
1041
1043ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm)
1044{
1045 /* TODO handle zero-length point arrays */
1046 uint32_t i;
1047 int in_hasz = FLAGS_GET_Z(pa->flags);
1048 int in_hasm = FLAGS_GET_M(pa->flags);
1049 POINT4D pt;
1050 POINTARRAY *pa_out = ptarray_construct_empty(hasz, hasm, pa->npoints);
1051
1052 for( i = 0; i < pa->npoints; i++ )
1053 {
1054 getPoint4d_p(pa, i, &pt);
1055 if( hasz && ! in_hasz )
1056 pt.z = 0.0;
1057 if( hasm && ! in_hasm )
1058 pt.m = 0.0;
1059 ptarray_append_point(pa_out, &pt, LW_TRUE);
1060 }
1061
1062 return pa_out;
1063}
1064
1065POINTARRAY *
1066ptarray_substring(POINTARRAY *ipa, double from, double to, double tolerance)
1067{
1068 POINTARRAY *dpa;
1069 POINT4D pt;
1070 POINT4D p1, p2;
1071 POINT4D *p1ptr=&p1; /* don't break strict-aliasing rule */
1072 POINT4D *p2ptr=&p2;
1073 int nsegs, i;
1074 double length, slength, tlength;
1075 int state = 0; /* 0=before, 1=inside */
1076
1077 /*
1078 * Create a dynamic pointarray with an initial capacity
1079 * equal to full copy of input points
1080 */
1082
1083 /* Compute total line length */
1084 length = ptarray_length_2d(ipa);
1085
1086
1087 LWDEBUGF(3, "Total length: %g", length);
1088
1089
1090 /* Get 'from' and 'to' lengths */
1091 from = length*from;
1092 to = length*to;
1093
1094
1095 LWDEBUGF(3, "From/To: %g/%g", from, to);
1096
1097
1098 tlength = 0;
1099 getPoint4d_p(ipa, 0, &p1);
1100 nsegs = ipa->npoints - 1;
1101 for ( i = 0; i < nsegs; i++ )
1102 {
1103 double dseg;
1104
1105 getPoint4d_p(ipa, i+1, &p2);
1106
1107
1108 LWDEBUGF(3 ,"Segment %d: (%g,%g,%g,%g)-(%g,%g,%g,%g)",
1109 i, p1.x, p1.y, p1.z, p1.m, p2.x, p2.y, p2.z, p2.m);
1110
1111
1112 /* Find the length of this segment */
1113 slength = distance2d_pt_pt((POINT2D *)p1ptr, (POINT2D *)p2ptr);
1114
1115 /*
1116 * We are before requested start.
1117 */
1118 if ( state == 0 ) /* before */
1119 {
1120
1121 LWDEBUG(3, " Before start");
1122
1123 if ( fabs ( from - ( tlength + slength ) ) <= tolerance )
1124 {
1125
1126 LWDEBUG(3, " Second point is our start");
1127
1128 /*
1129 * Second point is our start
1130 */
1131 ptarray_append_point(dpa, &p2, LW_FALSE);
1132 state=1; /* we're inside now */
1133 goto END;
1134 }
1135
1136 else if ( fabs(from - tlength) <= tolerance )
1137 {
1138
1139 LWDEBUG(3, " First point is our start");
1140
1141 /*
1142 * First point is our start
1143 */
1144 ptarray_append_point(dpa, &p1, LW_FALSE);
1145
1146 /*
1147 * We're inside now, but will check
1148 * 'to' point as well
1149 */
1150 state=1;
1151 }
1152
1153 /*
1154 * Didn't reach the 'from' point,
1155 * nothing to do
1156 */
1157 else if ( from > tlength + slength ) goto END;
1158
1159 else /* tlength < from < tlength+slength */
1160 {
1161
1162 LWDEBUG(3, " Seg contains first point");
1163
1164 /*
1165 * Our start is between first and
1166 * second point
1167 */
1168 dseg = (from - tlength) / slength;
1169
1170 interpolate_point4d(&p1, &p2, &pt, dseg);
1171
1172 ptarray_append_point(dpa, &pt, LW_FALSE);
1173
1174 /*
1175 * We're inside now, but will check
1176 * 'to' point as well
1177 */
1178 state=1;
1179 }
1180 }
1181
1182 if ( state == 1 ) /* inside */
1183 {
1184
1185 LWDEBUG(3, " Inside");
1186
1187 /*
1188 * 'to' point is our second point.
1189 */
1190 if ( fabs(to - ( tlength + slength ) ) <= tolerance )
1191 {
1192
1193 LWDEBUG(3, " Second point is our end");
1194
1195 ptarray_append_point(dpa, &p2, LW_FALSE);
1196 break; /* substring complete */
1197 }
1198
1199 /*
1200 * 'to' point is our first point.
1201 * (should only happen if 'to' is 0)
1202 */
1203 else if ( fabs(to - tlength) <= tolerance )
1204 {
1205
1206 LWDEBUG(3, " First point is our end");
1207
1208 ptarray_append_point(dpa, &p1, LW_FALSE);
1209
1210 break; /* substring complete */
1211 }
1212
1213 /*
1214 * Didn't reach the 'end' point,
1215 * just copy second point
1216 */
1217 else if ( to > tlength + slength )
1218 {
1219 ptarray_append_point(dpa, &p2, LW_FALSE);
1220 goto END;
1221 }
1222
1223 /*
1224 * 'to' point falls on this segment
1225 * Interpolate and break.
1226 */
1227 else if ( to < tlength + slength )
1228 {
1229
1230 LWDEBUG(3, " Seg contains our end");
1231
1232 dseg = (to - tlength) / slength;
1233 interpolate_point4d(&p1, &p2, &pt, dseg);
1234
1235 ptarray_append_point(dpa, &pt, LW_FALSE);
1236
1237 break;
1238 }
1239
1240 else
1241 {
1242 LWDEBUG(3, "Unhandled case");
1243 }
1244 }
1245
1246
1247END:
1248
1249 tlength += slength;
1250 memcpy(&p1, &p2, sizeof(POINT4D));
1251 }
1252
1253 LWDEBUGF(3, "Out of loop, ptarray has %d points", dpa->npoints);
1254
1255 return dpa;
1256}
1257
1258/*
1259 * Write into the *ret argument coordinates of the closes point on
1260 * the given segment to the reference input point.
1261 */
1262void
1263closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
1264{
1265 double r;
1266
1267 if ( FP_EQUALS(A->x, B->x) && FP_EQUALS(A->y, B->y) )
1268 {
1269 *ret = *A;
1270 return;
1271 }
1272
1273 /*
1274 * We use comp.graphics.algorithms Frequently Asked Questions method
1275 *
1276 * (1) AC dot AB
1277 * r = ----------
1278 * ||AB||^2
1279 * r has the following meaning:
1280 * r=0 P = A
1281 * r=1 P = B
1282 * r<0 P is on the backward extension of AB
1283 * r>1 P is on the forward extension of AB
1284 * 0<r<1 P is interior to AB
1285 *
1286 */
1287 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) );
1288
1289 if (r<0)
1290 {
1291 *ret = *A;
1292 return;
1293 }
1294 if (r>1)
1295 {
1296 *ret = *B;
1297 return;
1298 }
1299
1300 ret->x = A->x + ( (B->x - A->x) * r );
1301 ret->y = A->y + ( (B->y - A->y) * r );
1302 ret->z = A->z + ( (B->z - A->z) * r );
1303 ret->m = A->m + ( (B->m - A->m) * r );
1304}
1305
1306/*
1307 * Given a point, returns the location of closest point on pointarray
1308 * and, optionally, it's actual distance from the point array.
1309 */
1310double
1311ptarray_locate_point(const POINTARRAY *pa, const POINT4D *p4d, double *mindistout, POINT4D *proj4d)
1312{
1313 double mindist=DBL_MAX;
1314 double tlen, plen;
1315 uint32_t t, seg=0;
1316 POINT4D start4d, end4d, projtmp;
1317 POINT2D proj, p;
1318 const POINT2D *start = NULL, *end = NULL;
1319
1320 /* Initialize our 2D copy of the input parameter */
1321 p.x = p4d->x;
1322 p.y = p4d->y;
1323
1324 if ( ! proj4d ) proj4d = &projtmp;
1325
1326 /* Check for special cases (length 0 and 1) */
1327 if ( pa->npoints <= 1 )
1328 {
1329 if ( pa->npoints == 1 )
1330 {
1331 getPoint4d_p(pa, 0, proj4d);
1332 if ( mindistout )
1333 *mindistout = distance2d_pt_pt(&p, getPoint2d_cp(pa, 0));
1334 }
1335 return 0.0;
1336 }
1337
1338 start = getPoint2d_cp(pa, 0);
1339 /* Loop through pointarray looking for nearest segment */
1340 for (t=1; t<pa->npoints; t++)
1341 {
1342 double dist_sqr;
1343 end = getPoint2d_cp(pa, t);
1344 dist_sqr = distance2d_sqr_pt_seg(&p, start, end);
1345
1346 if (dist_sqr < mindist)
1347 {
1348 mindist = dist_sqr;
1349 seg=t-1;
1350 if ( mindist == 0 )
1351 {
1352 LWDEBUG(3, "Breaking on mindist=0");
1353 break;
1354 }
1355 }
1356
1357 start = end;
1358 }
1359 mindist = sqrt(mindist);
1360
1361 if ( mindistout ) *mindistout = mindist;
1362
1363 LWDEBUGF(3, "Closest segment: %d", seg);
1364 LWDEBUGF(3, "mindist: %g", mindist);
1365
1366 /*
1367 * We need to project the
1368 * point on the closest segment.
1369 */
1370 getPoint4d_p(pa, seg, &start4d);
1371 getPoint4d_p(pa, seg+1, &end4d);
1372 closest_point_on_segment(p4d, &start4d, &end4d, proj4d);
1373
1374 /* Copy 4D values into 2D holder */
1375 proj.x = proj4d->x;
1376 proj.y = proj4d->y;
1377
1378 LWDEBUGF(3, "Closest segment:%d, npoints:%d", seg, pa->npoints);
1379
1380 /* For robustness, force 1 when closest point == endpoint */
1381 if ( (seg >= (pa->npoints-2)) && p2d_same(&proj, end) )
1382 {
1383 return 1.0;
1384 }
1385
1386 LWDEBUGF(3, "Closest point on segment: %g,%g", proj.x, proj.y);
1387
1388 tlen = ptarray_length_2d(pa);
1389
1390 LWDEBUGF(3, "tlen %g", tlen);
1391
1392 /* Location of any point on a zero-length line is 0 */
1393 /* See http://trac.osgeo.org/postgis/ticket/1772#comment:2 */
1394 if ( tlen == 0 ) return 0;
1395
1396 plen=0;
1397 start = getPoint2d_cp(pa, 0);
1398 for (t=0; t<seg; t++, start=end)
1399 {
1400 end = getPoint2d_cp(pa, t+1);
1401 plen += distance2d_pt_pt(start, end);
1402
1403 LWDEBUGF(4, "Segment %d made plen %g", t, plen);
1404 }
1405
1406 plen+=distance2d_pt_pt(&proj, start);
1407
1408 LWDEBUGF(3, "plen %g, tlen %g", plen, tlen);
1409
1410 return plen/tlen;
1411}
1412
1422void
1424{
1425 uint32_t i;
1426 double x;
1427
1428 for (i=0; i<pa->npoints; i++)
1429 {
1430 memcpy(&x, getPoint_internal(pa, i), sizeof(double));
1431 if ( x < 0 ) x+= 360;
1432 else if ( x > 180 ) x -= 360;
1433 memcpy(getPoint_internal(pa, i), &x, sizeof(double));
1434 }
1435}
1436
1437
1438/*
1439 * Returns a POINTARRAY with consecutive equal points
1440 * removed. Equality test on all dimensions of input.
1441 *
1442 * Always returns a newly allocated object.
1443 */
1444static POINTARRAY *
1445ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
1446{
1447 POINTARRAY *out = ptarray_clone_deep(in);
1448 ptarray_remove_repeated_points_in_place(out, tolerance, minpoints);
1449 return out;
1450}
1451
1452POINTARRAY *
1453ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
1454{
1455 return ptarray_remove_repeated_points_minpoints(in, tolerance, 2);
1456}
1457
1458
1459void
1460ptarray_remove_repeated_points_in_place(POINTARRAY *pa, double tolerance, uint32_t min_points)
1461{
1462 uint32_t i;
1463 double tolsq = tolerance * tolerance;
1464 const POINT2D *last = NULL;
1465 const POINT2D *pt;
1466 uint32_t n_points = pa->npoints;
1467 uint32_t n_points_out = 1;
1468 size_t pt_size = ptarray_point_size(pa);
1469
1470 double dsq = FLT_MAX;
1471
1472 /* No-op on short inputs */
1473 if ( n_points <= min_points ) return;
1474
1475 last = getPoint2d_cp(pa, 0);
1476 void *p_to = ((char *)last) + pt_size;
1477 for (i = 1; i < n_points; i++)
1478 {
1479 int last_point = (i == n_points - 1);
1480
1481 /* Look straight into the abyss */
1482 pt = getPoint2d_cp(pa, i);
1483
1484 /* Don't drop points if we are running short of points */
1485 if (n_points + n_points_out > min_points + i)
1486 {
1487 if (tolerance > 0.0)
1488 {
1489 /* Only drop points that are within our tolerance */
1490 dsq = distance2d_sqr_pt_pt(last, pt);
1491 /* Allow any point but the last one to be dropped */
1492 if (!last_point && dsq <= tolsq)
1493 {
1494 continue;
1495 }
1496 }
1497 else
1498 {
1499 /* At tolerance zero, only skip exact dupes */
1500 if (memcmp((char*)pt, (char*)last, pt_size) == 0)
1501 continue;
1502 }
1503
1504 /* Got to last point, and it's not very different from */
1505 /* the point that preceded it. We want to keep the last */
1506 /* point, not the second-to-last one, so we pull our write */
1507 /* index back one value */
1508 if (last_point && n_points_out > 1 && tolerance > 0.0 && dsq <= tolsq)
1509 {
1510 n_points_out--;
1511 p_to -= pt_size;
1512 }
1513 }
1514
1515 /* Compact all remaining values to front of array */
1516 memcpy(p_to, pt, pt_size);
1517 n_points_out++;
1518 p_to += pt_size;
1519 last = pt;
1520 }
1521 /* Adjust array length */
1522 pa->npoints = n_points_out;
1523 return;
1524}
1525
1526/* Out of the points in pa [itfist .. itlast], finds the one that's farthest away from
1527 * the segment determined by pts[itfist] and pts[itlast].
1528 * Returns itfirst if no point was found futher away than max_distance_sqr
1529 */
1530static uint32_t
1531ptarray_dp_findsplit_in_place(const POINTARRAY *pts, uint32_t it_first, uint32_t it_last, double max_distance_sqr)
1532{
1533 uint32_t split = it_first;
1534 if ((it_first - it_last) < 2)
1535 return it_first;
1536
1537 const POINT2D *A = getPoint2d_cp(pts, it_first);
1538 const POINT2D *B = getPoint2d_cp(pts, it_last);
1539
1540 if (distance2d_sqr_pt_pt(A, B) < DBL_EPSILON)
1541 {
1542 /* If p1 == p2, we can just calculate the distance from each point to A */
1543 for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1544 {
1545 const POINT2D *pk = getPoint2d_cp(pts, itk);
1546 double distance_sqr = distance2d_sqr_pt_pt(pk, A);
1547 if (distance_sqr > max_distance_sqr)
1548 {
1549 split = itk;
1550 max_distance_sqr = distance_sqr;
1551 }
1552 }
1553 return split;
1554 }
1555
1556 /* This is based on distance2d_sqr_pt_seg, but heavily inlined here to avoid recalculations */
1557 double ba_x = (B->x - A->x);
1558 double ba_y = (B->y - A->y);
1559 double ab_length_sqr = (ba_x * ba_x + ba_y * ba_y);
1560 /* To avoid the division by ab_length_sqr in the 3rd path, we normalize here
1561 * and multiply in the first two paths [(dot_ac_ab < 0) and (> ab_length_sqr)] */
1562 max_distance_sqr *= ab_length_sqr;
1563 for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1564 {
1565 const POINT2D *C = getPoint2d_cp(pts, itk);
1566 double distance_sqr;
1567 double ca_x = (C->x - A->x);
1568 double ca_y = (C->y - A->y);
1569 double dot_ac_ab = (ca_x * ba_x + ca_y * ba_y);
1570
1571 if (dot_ac_ab <= 0.0)
1572 {
1573 distance_sqr = distance2d_sqr_pt_pt(C, A) * ab_length_sqr;
1574 }
1575 else if (dot_ac_ab >= ab_length_sqr)
1576 {
1577 distance_sqr = distance2d_sqr_pt_pt(C, B) * ab_length_sqr;
1578 }
1579 else
1580 {
1581 double s_numerator = ca_x * ba_y - ca_y * ba_x;
1582 distance_sqr = s_numerator * s_numerator; /* Missing division by ab_length_sqr on purpose */
1583 }
1584
1585 if (distance_sqr > max_distance_sqr)
1586 {
1587 split = itk;
1588 max_distance_sqr = distance_sqr;
1589 }
1590 }
1591 return split;
1592}
1593
1594void
1595ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
1596{
1597 /* Do not try to simplify really short things */
1598 if (pa->npoints < 3 || pa->npoints <= minpts)
1599 return;
1600
1601 /* We use this array to keep track of the points we are keeping, so
1602 * we store just TRUE / FALSE in their position */
1603 uint8_t *kept_points = lwalloc(sizeof(uint8_t) * pa->npoints);
1604 memset(kept_points, LW_FALSE, sizeof(uint8_t) * pa->npoints);
1605 kept_points[0] = LW_TRUE;
1606 kept_points[pa->npoints - 1] = LW_TRUE;
1607 uint32_t keptn = 2;
1608
1609 /* We use this array as a stack to store the iterators that we are going to need
1610 * in the following steps.
1611 * This is ~10% faster than iterating over @kept_points looking for them
1612 */
1613 uint32_t *iterator_stack = lwalloc(sizeof(uint32_t) * pa->npoints);
1614 iterator_stack[0] = 0;
1615 uint32_t iterator_stack_size = 1;
1616
1617 uint32_t it_first = 0;
1618 uint32_t it_last = pa->npoints - 1;
1619
1620 const double tolerance_sqr = tolerance * tolerance;
1621 /* For the first @minpts points we ignore the tolerance */
1622 double it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
1623
1624 while (iterator_stack_size)
1625 {
1626 uint32_t split = ptarray_dp_findsplit_in_place(pa, it_first, it_last, it_tol);
1627 if (split == it_first)
1628 {
1629 it_first = it_last;
1630 it_last = iterator_stack[--iterator_stack_size];
1631 }
1632 else
1633 {
1634 kept_points[split] = LW_TRUE;
1635 keptn++;
1636
1637 iterator_stack[iterator_stack_size++] = it_last;
1638 it_last = split;
1639 it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
1640 }
1641 }
1642
1643 const size_t pt_size = ptarray_point_size(pa);
1644 /* The first point is already in place, so we don't need to copy it */
1645 size_t kept_it = 1;
1646 if (keptn == 2)
1647 {
1648 /* If there are 2 points remaining, it has to be first and last as
1649 * we added those at the start */
1650 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1651 pa->serialized_pointlist + pt_size * (pa->npoints - 1),
1652 pt_size);
1653 }
1654 else
1655 {
1656 for (uint32_t i = 1; i < pa->npoints; i++)
1657 {
1658 if (kept_points[i])
1659 {
1660 memcpy(pa->serialized_pointlist + pt_size * kept_it,
1661 pa->serialized_pointlist + pt_size * i,
1662 pt_size);
1663 kept_it++;
1664 }
1665 }
1666 }
1667 pa->npoints = keptn;
1668
1669 lwfree(kept_points);
1670 lwfree(iterator_stack);
1671}
1672
1673/************************************************************************/
1674
1680double
1682{
1683 double dist = 0.0;
1684 uint32_t i;
1685 const POINT2D *a1;
1686 const POINT2D *a2;
1687 const POINT2D *a3;
1688
1689 if ( pts->npoints % 2 != 1 )
1690 lwerror("arc point array with even number of points");
1691
1692 a1 = getPoint2d_cp(pts, 0);
1693
1694 for ( i=2; i < pts->npoints; i += 2 )
1695 {
1696 a2 = getPoint2d_cp(pts, i-1);
1697 a3 = getPoint2d_cp(pts, i);
1698 dist += lw_arc_length(a1, a2, a3);
1699 a1 = a3;
1700 }
1701 return dist;
1702}
1703
1707double
1709{
1710 double dist = 0.0;
1711 uint32_t i;
1712 const POINT2D *frm;
1713 const POINT2D *to;
1714
1715 if ( pts->npoints < 2 ) return 0.0;
1716
1717 frm = getPoint2d_cp(pts, 0);
1718
1719 for ( i=1; i < pts->npoints; i++ )
1720 {
1721 to = getPoint2d_cp(pts, i);
1722
1723 dist += sqrt( ((frm->x - to->x)*(frm->x - to->x)) +
1724 ((frm->y - to->y)*(frm->y - to->y)) );
1725
1726 frm = to;
1727 }
1728 return dist;
1729}
1730
1735double
1737{
1738 double dist = 0.0;
1739 uint32_t i;
1740 POINT3DZ frm;
1741 POINT3DZ to;
1742
1743 if ( pts->npoints < 2 ) return 0.0;
1744
1745 /* compute 2d length if 3d is not available */
1746 if ( ! FLAGS_GET_Z(pts->flags) ) return ptarray_length_2d(pts);
1747
1748 getPoint3dz_p(pts, 0, &frm);
1749 for ( i=1; i < pts->npoints; i++ )
1750 {
1751 getPoint3dz_p(pts, i, &to);
1752 dist += sqrt( ((frm.x - to.x)*(frm.x - to.x)) +
1753 ((frm.y - to.y)*(frm.y - to.y)) +
1754 ((frm.z - to.z)*(frm.z - to.z)) );
1755 frm = to;
1756 }
1757 return dist;
1758}
1759
1760
1761
1765void
1767{
1768 uint32_t i;
1769 double x,y,z;
1770 POINT4D p4d;
1771
1772 LWDEBUG(2, "lwgeom_affine_ptarray start");
1773
1774 if ( FLAGS_GET_Z(pa->flags) )
1775 {
1776 LWDEBUG(3, " has z");
1777
1778 for (i=0; i<pa->npoints; i++)
1779 {
1780 getPoint4d_p(pa, i, &p4d);
1781 x = p4d.x;
1782 y = p4d.y;
1783 z = p4d.z;
1784 p4d.x = a->afac * x + a->bfac * y + a->cfac * z + a->xoff;
1785 p4d.y = a->dfac * x + a->efac * y + a->ffac * z + a->yoff;
1786 p4d.z = a->gfac * x + a->hfac * y + a->ifac * z + a->zoff;
1787 ptarray_set_point4d(pa, i, &p4d);
1788
1789 LWDEBUGF(3, " POINT %g %g %g => %g %g %g", x, y, z, p4d.x, p4d.y, p4d.z);
1790 }
1791 }
1792 else
1793 {
1794 LWDEBUG(3, " doesn't have z");
1795
1796 for (i=0; i<pa->npoints; i++)
1797 {
1798 getPoint4d_p(pa, i, &p4d);
1799 x = p4d.x;
1800 y = p4d.y;
1801 p4d.x = a->afac * x + a->bfac * y + a->xoff;
1802 p4d.y = a->dfac * x + a->efac * y + a->yoff;
1803 ptarray_set_point4d(pa, i, &p4d);
1804
1805 LWDEBUGF(3, " POINT %g %g => %g %g", x, y, p4d.x, p4d.y);
1806 }
1807 }
1808
1809 LWDEBUG(3, "lwgeom_affine_ptarray end");
1810
1811}
1812
1817#if 0
1818static int gluInvertMatrix(const double *m, double *invOut)
1819{
1820 double inv[16], det;
1821 int i;
1822
1823 inv[0] = m[5] * m[10] * m[15] -
1824 m[5] * m[11] * m[14] -
1825 m[9] * m[6] * m[15] +
1826 m[9] * m[7] * m[14] +
1827 m[13] * m[6] * m[11] -
1828 m[13] * m[7] * m[10];
1829
1830 inv[4] = -m[4] * m[10] * m[15] +
1831 m[4] * m[11] * m[14] +
1832 m[8] * m[6] * m[15] -
1833 m[8] * m[7] * m[14] -
1834 m[12] * m[6] * m[11] +
1835 m[12] * m[7] * m[10];
1836
1837 inv[8] = m[4] * m[9] * m[15] -
1838 m[4] * m[11] * m[13] -
1839 m[8] * m[5] * m[15] +
1840 m[8] * m[7] * m[13] +
1841 m[12] * m[5] * m[11] -
1842 m[12] * m[7] * m[9];
1843
1844 inv[12] = -m[4] * m[9] * m[14] +
1845 m[4] * m[10] * m[13] +
1846 m[8] * m[5] * m[14] -
1847 m[8] * m[6] * m[13] -
1848 m[12] * m[5] * m[10] +
1849 m[12] * m[6] * m[9];
1850
1851 inv[1] = -m[1] * m[10] * m[15] +
1852 m[1] * m[11] * m[14] +
1853 m[9] * m[2] * m[15] -
1854 m[9] * m[3] * m[14] -
1855 m[13] * m[2] * m[11] +
1856 m[13] * m[3] * m[10];
1857
1858 inv[5] = m[0] * m[10] * m[15] -
1859 m[0] * m[11] * m[14] -
1860 m[8] * m[2] * m[15] +
1861 m[8] * m[3] * m[14] +
1862 m[12] * m[2] * m[11] -
1863 m[12] * m[3] * m[10];
1864
1865 inv[9] = -m[0] * m[9] * m[15] +
1866 m[0] * m[11] * m[13] +
1867 m[8] * m[1] * m[15] -
1868 m[8] * m[3] * m[13] -
1869 m[12] * m[1] * m[11] +
1870 m[12] * m[3] * m[9];
1871
1872 inv[13] = m[0] * m[9] * m[14] -
1873 m[0] * m[10] * m[13] -
1874 m[8] * m[1] * m[14] +
1875 m[8] * m[2] * m[13] +
1876 m[12] * m[1] * m[10] -
1877 m[12] * m[2] * m[9];
1878
1879 inv[2] = m[1] * m[6] * m[15] -
1880 m[1] * m[7] * m[14] -
1881 m[5] * m[2] * m[15] +
1882 m[5] * m[3] * m[14] +
1883 m[13] * m[2] * m[7] -
1884 m[13] * m[3] * m[6];
1885
1886 inv[6] = -m[0] * m[6] * m[15] +
1887 m[0] * m[7] * m[14] +
1888 m[4] * m[2] * m[15] -
1889 m[4] * m[3] * m[14] -
1890 m[12] * m[2] * m[7] +
1891 m[12] * m[3] * m[6];
1892
1893 inv[10] = m[0] * m[5] * m[15] -
1894 m[0] * m[7] * m[13] -
1895 m[4] * m[1] * m[15] +
1896 m[4] * m[3] * m[13] +
1897 m[12] * m[1] * m[7] -
1898 m[12] * m[3] * m[5];
1899
1900 inv[14] = -m[0] * m[5] * m[14] +
1901 m[0] * m[6] * m[13] +
1902 m[4] * m[1] * m[14] -
1903 m[4] * m[2] * m[13] -
1904 m[12] * m[1] * m[6] +
1905 m[12] * m[2] * m[5];
1906
1907 inv[3] = -m[1] * m[6] * m[11] +
1908 m[1] * m[7] * m[10] +
1909 m[5] * m[2] * m[11] -
1910 m[5] * m[3] * m[10] -
1911 m[9] * m[2] * m[7] +
1912 m[9] * m[3] * m[6];
1913
1914 inv[7] = m[0] * m[6] * m[11] -
1915 m[0] * m[7] * m[10] -
1916 m[4] * m[2] * m[11] +
1917 m[4] * m[3] * m[10] +
1918 m[8] * m[2] * m[7] -
1919 m[8] * m[3] * m[6];
1920
1921 inv[11] = -m[0] * m[5] * m[11] +
1922 m[0] * m[7] * m[9] +
1923 m[4] * m[1] * m[11] -
1924 m[4] * m[3] * m[9] -
1925 m[8] * m[1] * m[7] +
1926 m[8] * m[3] * m[5];
1927
1928 inv[15] = m[0] * m[5] * m[10] -
1929 m[0] * m[6] * m[9] -
1930 m[4] * m[1] * m[10] +
1931 m[4] * m[2] * m[9] +
1932 m[8] * m[1] * m[6] -
1933 m[8] * m[2] * m[5];
1934
1935 det = m[0] * inv[0] + m[1] * inv[4] + m[2] * inv[8] + m[3] * inv[12];
1936
1937 if (det == 0)
1938 return LW_FALSE;
1939
1940 det = 1.0 / det;
1941
1942 for (i = 0; i < 16; i++)
1943 invOut[i] = inv[i] * det;
1944
1945 return LW_TRUE;
1946}
1947#endif
1948
1952void
1954{
1955 uint32_t i;
1956 POINT4D p4d;
1957 LWDEBUG(3, "ptarray_scale start");
1958 for (i=0; i<pa->npoints; i++)
1959 {
1960 getPoint4d_p(pa, i, &p4d);
1961 p4d.x *= fact->x;
1962 p4d.y *= fact->y;
1963 p4d.z *= fact->z;
1964 p4d.m *= fact->m;
1965 ptarray_set_point4d(pa, i, &p4d);
1966 }
1967 LWDEBUG(3, "ptarray_scale end");
1968}
1969
1970int
1972{
1973 return getPoint4d_p(pa, 0, pt);
1974}
1975
1976
1977/*
1978 * Stick an array of points to the given gridspec.
1979 * Return "gridded" points in *outpts and their number in *outptsn.
1980 *
1981 * Two consecutive points falling on the same grid cell are collapsed
1982 * into one single point.
1983 *
1984 */
1985void
1987{
1988 uint32_t i, j = 0;
1989 POINT4D *p, *p_out = NULL;
1990 int ndims = FLAGS_NDIMS(pa->flags);
1991 int has_z = FLAGS_GET_Z(pa->flags);
1992 int has_m = FLAGS_GET_M(pa->flags);
1993
1994 LWDEBUGF(2, "%s called on %p", __func__, pa);
1995
1996 for (i = 0; i < pa->npoints; i++)
1997 {
1998 /* Look straight into the abyss */
1999 p = (POINT4D*)(getPoint_internal(pa, i));
2000
2001 if (grid->xsize > 0)
2002 {
2003 p->x = rint((p->x - grid->ipx)/grid->xsize) * grid->xsize + grid->ipx;
2004 }
2005
2006 if (grid->ysize > 0)
2007 {
2008 p->y = rint((p->y - grid->ipy)/grid->ysize) * grid->ysize + grid->ipy;
2009 }
2010
2011 /* Read and round this point */
2012 /* Z is always in third position */
2013 if (has_z)
2014 {
2015 if (grid->zsize > 0)
2016 p->z = rint((p->z - grid->ipz)/grid->zsize) * grid->zsize + grid->ipz;
2017 }
2018 /* M might be in 3rd or 4th position */
2019 if (has_m)
2020 {
2021 /* In POINT M, M is in 3rd position */
2022 if (grid->msize > 0 && !has_z)
2023 p->z = rint((p->z - grid->ipm)/grid->msize) * grid->msize + grid->ipm;
2024 /* In POINT ZM, M is in 4th position */
2025 if (grid->msize > 0 && has_z)
2026 p->m = rint((p->m - grid->ipm)/grid->msize) * grid->msize + grid->ipm;
2027 }
2028
2029 /* Skip duplicates */
2030 if ( p_out && FP_EQUALS(p_out->x, p->x) && FP_EQUALS(p_out->y, p->y)
2031 && (ndims > 2 ? FP_EQUALS(p_out->z, p->z) : 1)
2032 && (ndims > 3 ? FP_EQUALS(p_out->m, p->m) : 1) )
2033 {
2034 continue;
2035 }
2036
2037 /* Write rounded values into the next available point */
2038 p_out = (POINT4D*)(getPoint_internal(pa, j++));
2039 p_out->x = p->x;
2040 p_out->y = p->y;
2041 if (ndims > 2)
2042 p_out->z = p->z;
2043 if (ndims > 3)
2044 p_out->m = p->m;
2045 }
2046
2047 /* Update output ptarray length */
2048 pa->npoints = j;
2049 return;
2050}
2051
2052
2053int
2055{
2056 const POINT2D *pt;
2057 int n = 0;
2058 uint32_t i;
2059 for ( i = 0; i < pa->npoints; i++ )
2060 {
2061 pt = getPoint2d_cp(pa, i);
2062 if ( gbox_contains_point2d(gbox, pt) )
2063 n++;
2064 }
2065 return n;
2066}
2067
2068
char * r
Definition cu_in_wkt.c:24
int gbox_contains_point2d(const GBOX *g, const POINT2D *p)
Definition gbox.c:350
int lw_arc_calculate_gbox_cartesian_2d(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, GBOX *gbox)
Definition gbox.c:453
#define LW_FALSE
Definition liblwgeom.h:108
void * lwrealloc(void *mem, size_t size)
Definition lwutil.c:235
double distance2d_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition measures.c:2397
#define LW_FAILURE
Definition liblwgeom.h:110
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:656
#define LW_SUCCESS
Definition liblwgeom.h:111
int getPoint2d_p(const POINTARRAY *pa, uint32_t n, POINT2D *point)
Definition lwgeom_api.c:349
#define FLAGS_GET_READONLY(flags)
Definition liblwgeom.h:183
#define FLAGS_GET_Z(flags)
Definition liblwgeom.h:179
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:242
#define FLAGS_NDIMS(flags)
Definition liblwgeom.h:193
#define FLAGS_GET_M(flags)
Definition liblwgeom.h:180
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:2407
#define FLAGS_GET_ZM(flags)
Definition liblwgeom.h:194
#define FLAGS_SET_READONLY(flags, value)
Definition liblwgeom.h:190
lwflags_t lwflags(int hasz, int hasm, int geodetic)
Construct a new flags bitmask.
Definition lwutil.c:471
#define LW_TRUE
Return types for functions with status returns.
Definition liblwgeom.h:107
void ptarray_set_point4d(POINTARRAY *pa, uint32_t n, const POINT4D *p4d)
Definition lwgeom_api.c:376
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)
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.
int lw_arc_side(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, const POINT2D *Q)
#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.
Definition lwalgorithm.c:96
#define LW_OUTSIDE
int lw_segment_side(const POINT2D *p1, const POINT2D *p2, const POINT2D *q)
lw_segment_side()
Definition lwalgorithm.c:65
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:86
int p2d_same(const POINT2D *p1, const POINT2D *p2)
Definition lwalgorithm.c:50
static double det(double m00, double m01, double m10, double m11)
#define LWDEBUG(level, msg)
Definition lwgeom_log.h:83
#define LWDEBUGF(level, msg,...)
Definition lwgeom_log.h:88
void lwerror(const char *fmt,...)
Write a notice out to the error handler.
Definition lwutil.c:190
void lwnotice(const char *fmt,...)
Write a notice out to the notice handler.
Definition lwutil.c:177
#define INT32_MAX
static uint8_t * getPoint_internal(const POINTARRAY *pa, uint32_t n)
Definition lwinline.h:67
static double distance2d_sqr_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition lwinline.h:35
static size_t ptarray_point_size(const POINTARRAY *pa)
Definition lwinline.h:48
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:91
int ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
Definition ptarray.c:746
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:1531
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:1423
POINTARRAY * ptarray_clone_deep(const POINTARRAY *in)
Deep clone a pointarray (also clones serialized pointlist)
Definition ptarray.c:634
POINTARRAY * ptarray_clone(const POINTARRAY *in)
Clone a POINTARRAY object.
Definition ptarray.c:665
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 ptarrayarc_contains_point(const POINTARRAY *pa, const POINT2D *pt)
For POINTARRAYs representing CIRCULARSTRINGS.
Definition ptarray.c:839
int ptarray_is_closed_2d(const POINTARRAY *in)
Definition ptarray.c:701
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:1066
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:727
double ptarray_length(const POINTARRAY *pts)
Find the 3d/2d length of the given POINTARRAY (depending on its dimensionality)
Definition ptarray.c:1736
double ptarray_signed_area(const POINTARRAY *pa)
Returns the area in cartesian units.
Definition ptarray.c:1003
static POINTARRAY * ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
Definition ptarray.c:1445
void closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
Definition ptarray.c:1263
int ptarray_startpoint(const POINTARRAY *pa, POINT4D *pt)
Definition ptarray.c:1971
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:687
void ptarray_grid_in_place(POINTARRAY *pa, const gridspec *grid)
Snap to grid.
Definition ptarray.c:1986
double ptarray_length_2d(const POINTARRAY *pts)
Find the 2d length of the given POINTARRAY (even if it's 3d)
Definition ptarray.c:1708
char ptarray_same(const POINTARRAY *pa1, const POINTARRAY *pa2)
Definition ptarray.c:484
int ptarray_is_closed_3d(const POINTARRAY *in)
Definition ptarray.c:714
void ptarray_remove_repeated_points_in_place(POINTARRAY *pa, double tolerance, uint32_t min_points)
Definition ptarray.c:1460
void ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
Definition ptarray.c:1595
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
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:1681
POINTARRAY * ptarray_merge(POINTARRAY *pa1, POINTARRAY *pa2)
Merge two given POINTARRAY and returns a pointer on the new aggregate one.
Definition ptarray.c:603
void ptarray_affine(POINTARRAY *pa, const AFFINE *a)
Affine transform a pointarray.
Definition ptarray.c:1766
int ptarray_contains_point(const POINTARRAY *pa, const POINT2D *pt)
Return 1 if the point is inside the POINTARRAY, -1 if it is outside, and 0 if it is on the boundary.
Definition ptarray.c:740
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:561
int ptarrayarc_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
Definition ptarray.c:845
void ptarray_free(POINTARRAY *pa)
Definition ptarray.c:327
int ptarray_isccw(const POINTARRAY *pa)
Definition ptarray.c:1034
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:509
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:1953
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:2054
POINTARRAY * ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
Definition ptarray.c:1453
POINTARRAY * ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm)
Definition ptarray.c:1043
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:1311
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:318
double zoff
Definition liblwgeom.h:318
double bfac
Definition liblwgeom.h:318
double ifac
Definition liblwgeom.h:318
double xoff
Definition liblwgeom.h:318
double dfac
Definition liblwgeom.h:318
double afac
Definition liblwgeom.h:318
double ffac
Definition liblwgeom.h:318
double cfac
Definition liblwgeom.h:318
double hfac
Definition liblwgeom.h:318
double efac
Definition liblwgeom.h:318
double yoff
Definition liblwgeom.h:318
double ymax
Definition liblwgeom.h:343
double xmax
Definition liblwgeom.h:341
double ymin
Definition liblwgeom.h:342
double xmin
Definition liblwgeom.h:340
double y
Definition liblwgeom.h:376
double x
Definition liblwgeom.h:376
double z
Definition liblwgeom.h:382
double x
Definition liblwgeom.h:382
double y
Definition liblwgeom.h:382
double m
Definition liblwgeom.h:400
double x
Definition liblwgeom.h:400
double z
Definition liblwgeom.h:400
double y
Definition liblwgeom.h:400
lwflags_t flags
Definition liblwgeom.h:417
uint32_t maxpoints
Definition liblwgeom.h:414
uint32_t npoints
Definition liblwgeom.h:413
uint8_t * serialized_pointlist
Definition liblwgeom.h:420
double ipm
Definition liblwgeom.h:1345
double zsize
Definition liblwgeom.h:1348
double ysize
Definition liblwgeom.h:1347
double xsize
Definition liblwgeom.h:1346
double ipx
Definition liblwgeom.h:1342
double msize
Definition liblwgeom.h:1349
double ipy
Definition liblwgeom.h:1343
double ipz
Definition liblwgeom.h:1344
Snap-to-grid.
Definition liblwgeom.h:1341