PostGIS  3.2.2dev-r@@SVN_REVISION@@
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 
36 int
38 {
39  if ( ! pa ) return LW_FALSE;
40  return FLAGS_GET_Z(pa->flags);
41 }
42 
43 int
45 {
46  if ( ! pa ) return LW_FALSE;
47  return FLAGS_GET_M(pa->flags);
48 }
49 
51 ptarray_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 
59 ptarray_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 */
84 int
85 ptarray_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 
146 int
147 ptarray_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 
176 int
177 ptarray_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 */
258 int
259 ptarray_remove_point(POINTARRAY *pa, uint32_t where)
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 
291 POINTARRAY* 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 
304 POINTARRAY*
305 ptarray_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  {
315  pa->serialized_pointlist = lwalloc(ptarray_point_size(pa) * npoints);
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 
326 void
328 {
329  if (pa)
330  {
331  if (pa->serialized_pointlist && (!FLAGS_GET_READONLY(pa->flags)))
333  lwfree(pa);
334  }
335 }
336 
337 
338 void
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 
367 POINTARRAY*
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 
386 void
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 
412 POINTARRAY *
413 ptarray_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);
430  ptarray_append_point(opa, &p1, LW_FALSE);
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 
483 char
484 ptarray_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 
508 POINTARRAY *
509 ptarray_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 
560 POINTARRAY *
561 ptarray_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 
602 POINTARRAY *
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 
611  pa = ptarray_construct( FLAGS_GET_Z(pa1->flags),
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 
633 POINTARRAY *
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 
664 POINTARRAY *
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 
686 int
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 
700 int
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 
713 int
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 
726 int
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 
745 int
747 {
748  const POINT2D *seg1, *seg2;
749  int wn = 0;
750 
751  seg1 = getPoint2d_cp(pa, 0);
752  seg2 = getPoint2d_cp(pa, pa->npoints-1);
753  if (!p2d_same(seg1, seg2))
754  lwerror("ptarray_contains_point called on unclosed ring");
755 
756  for (uint32_t i = 1; i < pa->npoints; i++)
757  {
758  double side, ymin, ymax;
759 
760  seg2 = getPoint2d_cp(pa, i);
761 
762  /* Zero length segments are ignored. */
763  if (p2d_same(seg1, seg2))
764  {
765  seg1 = seg2;
766  continue;
767  }
768 
769  ymin = FP_MIN(seg1->y, seg2->y);
770  ymax = FP_MAX(seg1->y, seg2->y);
771 
772  /* Only test segments in our vertical range */
773  if (pt->y > ymax || pt->y < ymin)
774  {
775  seg1 = seg2;
776  continue;
777  }
778 
779  side = lw_segment_side(seg1, seg2, pt);
780 
781  /*
782  * A point on the boundary of a ring is not contained.
783  * WAS: if (fabs(side) < 1e-12), see #852
784  */
785  if ((side == 0) && lw_pt_in_seg(pt, seg1, seg2))
786  {
787  return LW_BOUNDARY;
788  }
789 
790  /*
791  * If the point is to the left of the line, and it's rising,
792  * then the line is to the right of the point and
793  * circling counter-clockwise, so increment.
794  */
795  if ((side < 0) && (seg1->y <= pt->y) && (pt->y < seg2->y))
796  {
797  wn++;
798  }
799 
800  /*
801  * If the point is to the right of the line, and it's falling,
802  * then the line is to the right of the point and circling
803  * clockwise, so decrement.
804  */
805  else if ( (side > 0) && (seg2->y <= pt->y) && (pt->y < seg1->y) )
806  {
807  wn--;
808  }
809 
810  seg1 = seg2;
811  }
812 
813  /* wn == 0 => Outside, wn != 0 => Inside */
814  return wn == 0 ? LW_OUTSIDE : LW_INSIDE;
815 }
816 
817 int
818 ptarray_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
819 {
820  // A valid linestring must have at least 2 point
821  if (pa->npoints < 2)
822  lwerror("%s called on invalid linestring", __func__);
823 
824  int intersections = 0;
825  double px = p->x;
826  double py = p->y;
827 
828  // Iterate through each edge of the polygon
829  for (uint32_t i = 0; i < pa->npoints-1; ++i)
830  {
831  const POINT2D* p1 = getPoint2d_cp(pa, i);
832  const POINT2D* p2 = getPoint2d_cp(pa, i+1);
833 
834  /* Skip zero-length edges */
835  if (p2d_same(p1, p2))
836  continue;
837 
838  /* --- Step 1: Check if the point is ON the boundary edge --- */
839  if (lw_pt_on_segment(p1, p2, p))
840  {
841  *on_boundary = LW_TRUE;
842  return 0;
843  }
844 
845  /* --- Step 2: Perform the Ray Casting intersection test --- */
846 
847  /*
848  * Check if the horizontal ray from p intersects the edge (p1, p2).
849  * This is the core condition for handling vertices correctly:
850  * - One vertex must be strictly above the ray (py < vertex.y)
851  * - The other must be on or below the ray (py >= vertex.y)
852  */
853  if (((p1->y <= py) && (py < p2->y)) || ((p2->y <= py) && (py < p1->y)))
854  {
855  /*
856  * Calculate the x-coordinate where the edge intersects the ray's horizontal line.
857  * Formula: x = x1 + (y - y1) * (x2 - x1) / (y2 - y1)
858  */
859  double x_intersection = p1->x + (py - p1->y) * (p2->x - p1->x) / (p2->y - p1->y);
860 
861  /*
862  * If the intersection point is to the right of our test point,
863  * it's a valid "crossing".
864  */
865  if (x_intersection > px)
866  {
867  intersections++;
868  }
869  }
870  }
871 
872  return intersections;
873 }
874 
875 
893 static int
894 circle_raycast_intersections(const POINT2D *center, double radius, double ray, POINT2D *i0, POINT2D *i1)
895 {
896  // Calculate the vertical distance from the circle's center to the horizontal line.
897  double dy = ray - center->y;
898 
899  // If the absolute vertical distance is greater than the radius, there are no intersections.
900  if (fabs(dy) > radius)
901  return 0;
902 
903  // Use the Pythagorean theorem to find the horizontal distance (dx) from the
904  // center's x-coordinate to the intersection points.
905  // dx^2 + dy^2 = radius^2 => dx^2 = radius^2 - dy^2
906  double dx_squared = radius * radius - dy * dy;
907 
908  // Case 1: One intersection (tangent)
909  // This occurs when the line just touches the top or bottom of the circle.
910  // dx_squared will be zero. We check against a small epsilon for floating-point safety.
911  if (FP_EQUALS(dx_squared, 0.0))
912  {
913  i0->x = center->x;
914  i0->y = ray;
915  return 1;
916  }
917 
918  // Case 2: Two intersections
919  // The line cuts through the circle.
920  double dx = sqrt(dx_squared);
921 
922  // The first intersection point has the smaller x-value.
923  i0->x = center->x - dx;
924  i0->y = ray;
925 
926  // The second intersection point has the larger x-value.
927  i1->x = center->x + dx;
928  i1->y = ray;
929 
930  return 2;
931 }
932 
933 
934 int
935 ptarrayarc_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
936 {
937  int intersections = 0;
938  double px = p->x;
939  double py = p->y;
940 
941  assert(on_boundary);
942 
943  // A valid circular arc must have at least 3 vertices (circle).
944  if (pa->npoints < 3)
945  lwerror("%s called on invalid circularstring", __func__);
946 
947  if (pa->npoints % 2 == 0)
948  lwerror("%s called with even number of points", __func__);
949 
950  // Iterate through each arc of the circularstring
951  for (uint32_t i = 1; i < pa->npoints-1; i +=2)
952  {
953  const POINT2D* p0 = getPoint2d_cp(pa, i-1);
954  const POINT2D* p1 = getPoint2d_cp(pa, i);
955  const POINT2D* p2 = getPoint2d_cp(pa, i+1);
956  POINT2D center = {0,0};
957  double radius, d;
958  GBOX gbox;
959 
960  // Skip zero-length arc
961  if (lw_arc_is_pt(p0, p1, p2))
962  continue;
963 
964  // --- Step 1: Check if the point is ON the boundary edge ---
965  if (p2d_same(p0, p) || p2d_same(p1, p) || p2d_same(p2, p))
966  {
967  *on_boundary = LW_TRUE;
968  return 0;
969  }
970 
971  // Calculate some important pieces
972  radius = lw_arc_center(p0, p1, p2, &center);
973 
974  d = distance2d_pt_pt(p, &center);
975  if (FP_EQUALS(d, radius) && lw_pt_in_arc(p, p0, p1, p2))
976  {
977  *on_boundary = LW_TRUE;
978  return 0;
979  }
980 
981  // --- Step 2: Perform the Ray Casting intersection test ---
982 
983  // Only process arcs that our ray crosses
984  lw_arc_calculate_gbox_cartesian_2d(p0, p1, p2, &gbox);
985  if ((gbox.ymin <= py) && (py < gbox.ymax))
986  {
987  // Point of intersection on the circle that defines the arc
988  POINT2D i0, i1;
989 
990  // How many points of intersection are there
991  int iCount = circle_raycast_intersections(&center, radius, py, &i0, &i1);
992 
993  // Nothing to see here
994  if (iCount == 0)
995  continue;
996 
997  // Cannot think of a case where a grazing is not a
998  // no-op
999  if (iCount == 1)
1000  continue;
1001 
1002  // So we must have 2 intersections
1003  // Only increment the counter for intersections to the right
1004  // of the test point
1005  if (i0.x > px && lw_pt_in_arc(&i0, p0, p1, p2))
1006  intersections++;
1007 
1008  if (i1.x > px && lw_pt_in_arc(&i1, p0, p1, p2))
1009  intersections++;
1010  }
1011  }
1012 
1013  return intersections;
1014 }
1015 
1016 int
1017 ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
1018 {
1019  int wn = 0;
1020  uint32_t i;
1021  double side;
1022  const POINT2D *seg1;
1023  const POINT2D *seg2;
1024  double ymin, ymax;
1025 
1026  seg1 = getPoint2d_cp(pa, 0);
1027  seg2 = getPoint2d_cp(pa, pa->npoints-1);
1028  if ( check_closed && ! p2d_same(seg1, seg2) )
1029  lwerror("ptarray_contains_point called on unclosed ring");
1030 
1031  for ( i=1; i < pa->npoints; i++ )
1032  {
1033  seg2 = getPoint2d_cp(pa, i);
1034 
1035  /* Zero length segments are ignored. */
1036  if ( seg1->x == seg2->x && seg1->y == seg2->y )
1037  {
1038  seg1 = seg2;
1039  continue;
1040  }
1041 
1042  ymin = FP_MIN(seg1->y, seg2->y);
1043  ymax = FP_MAX(seg1->y, seg2->y);
1044 
1045  /* Only test segments in our vertical range */
1046  if ( pt->y > ymax || pt->y < ymin )
1047  {
1048  seg1 = seg2;
1049  continue;
1050  }
1051 
1052  side = lw_segment_side(seg1, seg2, pt);
1053 
1054  /*
1055  * A point on the boundary of a ring is not contained.
1056  * WAS: if (fabs(side) < 1e-12), see #852
1057  */
1058  if ( (side == 0) && lw_pt_in_seg(pt, seg1, seg2) )
1059  {
1060  return LW_BOUNDARY;
1061  }
1062 
1063  /*
1064  * If the point is to the left of the line, and it's rising,
1065  * then the line is to the right of the point and
1066  * circling counter-clockwise, so increment.
1067  */
1068  if ( (side < 0) && (seg1->y <= pt->y) && (pt->y < seg2->y) )
1069  {
1070  wn++;
1071  }
1072 
1073  /*
1074  * If the point is to the right of the line, and it's falling,
1075  * then the line is to the right of the point and circling
1076  * clockwise, so decrement.
1077  */
1078  else if ( (side > 0) && (seg2->y <= pt->y) && (pt->y < seg1->y) )
1079  {
1080  wn--;
1081  }
1082 
1083  seg1 = seg2;
1084  }
1085 
1086  /* Sent out the winding number for calls that are building on this as a primitive */
1087  if ( winding_number )
1088  *winding_number = wn;
1089 
1090  /* Outside */
1091  if (wn == 0)
1092  {
1093  return LW_OUTSIDE;
1094  }
1095 
1096  /* Inside */
1097  return LW_INSIDE;
1098 }
1099 
1109 int
1111 {
1112  int on_boundary = LW_FALSE;
1113  int intersections;
1114  if (!ptarray_is_closed_2d(pa))
1115  lwerror("%s called on unclosed ring", __func__);
1116 
1117  intersections = ptarrayarc_raycast_intersections(pa, pt, &on_boundary);
1118  if (on_boundary)
1119  return LW_BOUNDARY;
1120  else
1121  return (intersections % 2) ? LW_INSIDE : LW_OUTSIDE;
1122 }
1123 
1124 int
1125 ptarrayarc_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
1126 {
1127  int wn = 0;
1128  uint32_t i;
1129  int side;
1130  const POINT2D *seg1;
1131  const POINT2D *seg2;
1132  const POINT2D *seg3;
1133  GBOX gbox;
1134 
1135  /* Check for not an arc ring (always have odd # of points) */
1136  if ( (pa->npoints % 2) == 0 )
1137  {
1138  lwerror("ptarrayarc_contains_point called with even number of points");
1139  return LW_OUTSIDE;
1140  }
1141 
1142  /* Check for not an arc ring (always have >= 3 points) */
1143  if ( pa->npoints < 3 )
1144  {
1145  lwerror("ptarrayarc_contains_point called too-short pointarray");
1146  return LW_OUTSIDE;
1147  }
1148 
1149  /* Check for unclosed case */
1150  seg1 = getPoint2d_cp(pa, 0);
1151  seg3 = getPoint2d_cp(pa, pa->npoints-1);
1152  if ( check_closed && ! p2d_same(seg1, seg3) )
1153  {
1154  lwerror("ptarrayarc_contains_point called on unclosed ring");
1155  return LW_OUTSIDE;
1156  }
1157  /* OK, it's closed. Is it just one circle? */
1158  else if ( p2d_same(seg1, seg3) && pa->npoints == 3 )
1159  {
1160  double radius, d;
1161  POINT2D c;
1162  seg2 = getPoint2d_cp(pa, 1);
1163 
1164  /* Wait, it's just a point, so it can't contain anything */
1165  if ( lw_arc_is_pt(seg1, seg2, seg3) )
1166  return LW_OUTSIDE;
1167 
1168  /* See if the point is within the circle radius */
1169  radius = lw_arc_center(seg1, seg2, seg3, &c);
1170  d = distance2d_pt_pt(pt, &c);
1171  if ( FP_EQUALS(d, radius) )
1172  return LW_BOUNDARY; /* Boundary of circle */
1173  else if ( d < radius )
1174  return LW_INSIDE; /* Inside circle */
1175  else
1176  return LW_OUTSIDE; /* Outside circle */
1177  }
1178  else if ( p2d_same(seg1, pt) || p2d_same(seg3, pt) )
1179  {
1180  return LW_BOUNDARY; /* Boundary case */
1181  }
1182 
1183  /* Start on the ring */
1184  seg1 = getPoint2d_cp(pa, 0);
1185  for ( i=1; i < pa->npoints; i += 2 )
1186  {
1187  seg2 = getPoint2d_cp(pa, i);
1188  seg3 = getPoint2d_cp(pa, i+1);
1189 
1190  /* Catch an easy boundary case */
1191  if( p2d_same(seg3, pt) )
1192  return LW_BOUNDARY;
1193 
1194  /* Skip arcs that have no size */
1195  if ( lw_arc_is_pt(seg1, seg2, seg3) )
1196  {
1197  seg1 = seg3;
1198  continue;
1199  }
1200 
1201  /* Only test segments in our vertical range */
1202  lw_arc_calculate_gbox_cartesian_2d(seg1, seg2, seg3, &gbox);
1203  if ( pt->y > gbox.ymax || pt->y < gbox.ymin )
1204  {
1205  seg1 = seg3;
1206  continue;
1207  }
1208 
1209  /* Outside of horizontal range, and not between end points we also skip */
1210  if ( (pt->x > gbox.xmax || pt->x < gbox.xmin) &&
1211  (pt->y > FP_MAX(seg1->y, seg3->y) || pt->y < FP_MIN(seg1->y, seg3->y)) )
1212  {
1213  seg1 = seg3;
1214  continue;
1215  }
1216 
1217  side = lw_arc_side(seg1, seg2, seg3, pt);
1218 
1219  /* On the boundary */
1220  if ( (side == 0) && lw_pt_in_arc(pt, seg1, seg2, seg3) )
1221  {
1222  return LW_BOUNDARY;
1223  }
1224 
1225  /* Going "up"! Point to left of arc. */
1226  if ( side < 0 && (seg1->y <= pt->y) && (pt->y < seg3->y) )
1227  {
1228  wn++;
1229  }
1230 
1231  /* Going "down"! */
1232  if ( side > 0 && (seg3->y <= pt->y) && (pt->y < seg1->y) )
1233  {
1234  wn--;
1235  }
1236 
1237  /* Inside the arc! */
1238  if ( pt->x <= gbox.xmax && pt->x >= gbox.xmin )
1239  {
1240  POINT2D C;
1241  double radius = lw_arc_center(seg1, seg2, seg3, &C);
1242  double d = distance2d_pt_pt(pt, &C);
1243 
1244  /* On the boundary! */
1245  if ( d == radius )
1246  return LW_BOUNDARY;
1247 
1248  /* Within the arc! */
1249  if ( d < radius )
1250  {
1251  /* Left side, increment winding number */
1252  if ( side < 0 )
1253  wn++;
1254  /* Right side, decrement winding number */
1255  if ( side > 0 )
1256  wn--;
1257  }
1258  }
1259 
1260  seg1 = seg3;
1261  }
1262 
1263  /* Sent out the winding number for calls that are building on this as a primitive */
1264  if ( winding_number )
1265  *winding_number = wn;
1266 
1267  /* Outside */
1268  if (wn == 0)
1269  {
1270  return LW_OUTSIDE;
1271  }
1272 
1273  /* Inside */
1274  return LW_INSIDE;
1275 }
1276 
1282 double
1284 {
1285  const POINT2D *P1;
1286  const POINT2D *P2;
1287  const POINT2D *P3;
1288  double sum = 0.0;
1289  double x0, x, y1, y2;
1290  uint32_t i;
1291 
1292  if (! pa || pa->npoints < 3 )
1293  return 0.0;
1294 
1295  P1 = getPoint2d_cp(pa, 0);
1296  P2 = getPoint2d_cp(pa, 1);
1297  x0 = P1->x;
1298  for ( i = 2; i < pa->npoints; i++ )
1299  {
1300  P3 = getPoint2d_cp(pa, i);
1301  x = P2->x - x0;
1302  y1 = P3->y;
1303  y2 = P1->y;
1304  sum += x * (y2-y1);
1305 
1306  /* Move forwards! */
1307  P1 = P2;
1308  P2 = P3;
1309  }
1310  return sum / 2.0;
1311 }
1312 
1313 int
1315 {
1316  double area = 0;
1317  area = ptarray_signed_area(pa);
1318  if ( area > 0 ) return LW_FALSE;
1319  else return LW_TRUE;
1320 }
1321 
1322 POINTARRAY*
1323 ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm, double zval, double mval)
1324 {
1325  /* TODO handle zero-length point arrays */
1326  uint32_t i;
1327  int in_hasz = FLAGS_GET_Z(pa->flags);
1328  int in_hasm = FLAGS_GET_M(pa->flags);
1329  POINT4D pt;
1330  POINTARRAY *pa_out = ptarray_construct_empty(hasz, hasm, pa->npoints);
1331 
1332  for( i = 0; i < pa->npoints; i++ )
1333  {
1334  getPoint4d_p(pa, i, &pt);
1335  if( hasz && ! in_hasz )
1336  pt.z = zval;
1337  if( hasm && ! in_hasm )
1338  pt.m = mval;
1339  ptarray_append_point(pa_out, &pt, LW_TRUE);
1340  }
1341 
1342  return pa_out;
1343 }
1344 
1345 POINTARRAY *
1346 ptarray_substring(POINTARRAY *ipa, double from, double to, double tolerance)
1347 {
1348  POINTARRAY *dpa;
1349  POINT4D pt;
1350  POINT4D p1, p2;
1351  POINT4D *p1ptr=&p1; /* don't break strict-aliasing rule */
1352  POINT4D *p2ptr=&p2;
1353  int nsegs, i;
1354  double length, slength, tlength;
1355  int state = 0; /* 0=before, 1=inside */
1356 
1357  /*
1358  * Create a dynamic pointarray with an initial capacity
1359  * equal to full copy of input points
1360  */
1362 
1363  /* Compute total line length */
1364  length = ptarray_length_2d(ipa);
1365 
1366 
1367  LWDEBUGF(3, "Total length: %g", length);
1368 
1369 
1370  /* Get 'from' and 'to' lengths */
1371  from = length*from;
1372  to = length*to;
1373 
1374 
1375  LWDEBUGF(3, "From/To: %g/%g", from, to);
1376 
1377 
1378  tlength = 0;
1379  getPoint4d_p(ipa, 0, &p1);
1380  nsegs = ipa->npoints - 1;
1381  for ( i = 0; i < nsegs; i++ )
1382  {
1383  double dseg;
1384 
1385  getPoint4d_p(ipa, i+1, &p2);
1386 
1387 
1388  LWDEBUGF(3 ,"Segment %d: (%g,%g,%g,%g)-(%g,%g,%g,%g)",
1389  i, p1.x, p1.y, p1.z, p1.m, p2.x, p2.y, p2.z, p2.m);
1390 
1391 
1392  /* Find the length of this segment */
1393  slength = distance2d_pt_pt((POINT2D *)p1ptr, (POINT2D *)p2ptr);
1394 
1395  /*
1396  * We are before requested start.
1397  */
1398  if ( state == 0 ) /* before */
1399  {
1400 
1401  LWDEBUG(3, " Before start");
1402 
1403  if ( fabs ( from - ( tlength + slength ) ) <= tolerance )
1404  {
1405 
1406  LWDEBUG(3, " Second point is our start");
1407 
1408  /*
1409  * Second point is our start
1410  */
1411  ptarray_append_point(dpa, &p2, LW_FALSE);
1412  state=1; /* we're inside now */
1413  goto END;
1414  }
1415 
1416  else if ( fabs(from - tlength) <= tolerance )
1417  {
1418 
1419  LWDEBUG(3, " First point is our start");
1420 
1421  /*
1422  * First point is our start
1423  */
1424  ptarray_append_point(dpa, &p1, LW_FALSE);
1425 
1426  /*
1427  * We're inside now, but will check
1428  * 'to' point as well
1429  */
1430  state=1;
1431  }
1432 
1433  /*
1434  * Didn't reach the 'from' point,
1435  * nothing to do
1436  */
1437  else if ( from > tlength + slength ) goto END;
1438 
1439  else /* tlength < from < tlength+slength */
1440  {
1441 
1442  LWDEBUG(3, " Seg contains first point");
1443 
1444  /*
1445  * Our start is between first and
1446  * second point
1447  */
1448  dseg = (from - tlength) / slength;
1449 
1450  interpolate_point4d(&p1, &p2, &pt, dseg);
1451 
1452  ptarray_append_point(dpa, &pt, LW_FALSE);
1453 
1454  /*
1455  * We're inside now, but will check
1456  * 'to' point as well
1457  */
1458  state=1;
1459  }
1460  }
1461 
1462  if ( state == 1 ) /* inside */
1463  {
1464 
1465  LWDEBUG(3, " Inside");
1466 
1467  /*
1468  * 'to' point is our second point.
1469  */
1470  if ( fabs(to - ( tlength + slength ) ) <= tolerance )
1471  {
1472 
1473  LWDEBUG(3, " Second point is our end");
1474 
1475  ptarray_append_point(dpa, &p2, LW_FALSE);
1476  break; /* substring complete */
1477  }
1478 
1479  /*
1480  * 'to' point is our first point.
1481  * (should only happen if 'to' is 0)
1482  */
1483  else if ( fabs(to - tlength) <= tolerance )
1484  {
1485 
1486  LWDEBUG(3, " First point is our end");
1487 
1488  ptarray_append_point(dpa, &p1, LW_FALSE);
1489 
1490  break; /* substring complete */
1491  }
1492 
1493  /*
1494  * Didn't reach the 'end' point,
1495  * just copy second point
1496  */
1497  else if ( to > tlength + slength )
1498  {
1499  ptarray_append_point(dpa, &p2, LW_FALSE);
1500  goto END;
1501  }
1502 
1503  /*
1504  * 'to' point falls on this segment
1505  * Interpolate and break.
1506  */
1507  else if ( to < tlength + slength )
1508  {
1509 
1510  LWDEBUG(3, " Seg contains our end");
1511 
1512  dseg = (to - tlength) / slength;
1513  interpolate_point4d(&p1, &p2, &pt, dseg);
1514 
1515  ptarray_append_point(dpa, &pt, LW_FALSE);
1516 
1517  break;
1518  }
1519 
1520  else
1521  {
1522  LWDEBUG(3, "Unhandled case");
1523  }
1524  }
1525 
1526 
1527 END:
1528 
1529  tlength += slength;
1530  memcpy(&p1, &p2, sizeof(POINT4D));
1531  }
1532 
1533  LWDEBUGF(3, "Out of loop, ptarray has %d points", dpa->npoints);
1534 
1535  return dpa;
1536 }
1537 
1538 /*
1539  * Write into the *ret argument coordinates of the closes point on
1540  * the given segment to the reference input point.
1541  */
1542 void
1543 closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
1544 {
1545  double r;
1546 
1547  if ( FP_EQUALS(A->x, B->x) && FP_EQUALS(A->y, B->y) )
1548  {
1549  *ret = *A;
1550  return;
1551  }
1552 
1553  /*
1554  * We use comp.graphics.algorithms Frequently Asked Questions method
1555  *
1556  * (1) AC dot AB
1557  * r = ----------
1558  * ||AB||^2
1559  * r has the following meaning:
1560  * r=0 P = A
1561  * r=1 P = B
1562  * r<0 P is on the backward extension of AB
1563  * r>1 P is on the forward extension of AB
1564  * 0<r<1 P is interior to AB
1565  *
1566  */
1567  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) );
1568 
1569  if (r<=0)
1570  {
1571  *ret = *A;
1572  return;
1573  }
1574  if (r>=1)
1575  {
1576  *ret = *B;
1577  return;
1578  }
1579 
1580  ret->x = A->x + ( (B->x - A->x) * r );
1581  ret->y = A->y + ( (B->y - A->y) * r );
1582  ret->z = A->z + ( (B->z - A->z) * r );
1583  ret->m = A->m + ( (B->m - A->m) * r );
1584 }
1585 
1586 int
1587 ptarray_closest_segment_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
1588 {
1589  const POINT2D *start = getPoint2d_cp(pa, 0), *end = NULL;
1590  uint32_t t, seg=0;
1591  double mindist=DBL_MAX;
1592 
1593  /* Loop through pointarray looking for nearest segment */
1594  for (t=1; t<pa->npoints; t++)
1595  {
1596  double dist_sqr;
1597  end = getPoint2d_cp(pa, t);
1598  dist_sqr = distance2d_sqr_pt_seg(qp, start, end);
1599 
1600  if (dist_sqr < mindist)
1601  {
1602  mindist = dist_sqr;
1603  seg=t-1;
1604  if ( mindist == 0 )
1605  {
1606  LWDEBUG(3, "Breaking on mindist=0");
1607  break;
1608  }
1609  }
1610 
1611  start = end;
1612  }
1613 
1614  if ( dist ) *dist = sqrt(mindist);
1615  return seg;
1616 }
1617 
1618 
1619 int
1620 ptarray_closest_vertex_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
1621 {
1622  uint32_t t, pn=0;
1623  const POINT2D *p;
1624  double mindist = DBL_MAX;
1625 
1626  /* Loop through pointarray looking for nearest segment */
1627  for (t=0; t<pa->npoints; t++)
1628  {
1629  double dist_sqr;
1630  p = getPoint2d_cp(pa, t);
1631  dist_sqr = distance2d_sqr_pt_pt(p, qp);
1632 
1633  if (dist_sqr < mindist)
1634  {
1635  mindist = dist_sqr;
1636  pn = t;
1637  if ( mindist == 0 )
1638  {
1639  LWDEBUG(3, "Breaking on mindist=0");
1640  break;
1641  }
1642  }
1643  }
1644  if ( dist ) *dist = sqrt(mindist);
1645  return pn;
1646 }
1647 
1648 /*
1649  * Given a point, returns the location of closest point on pointarray
1650  * and, optionally, it's actual distance from the point array.
1651  */
1652 double
1653 ptarray_locate_point(const POINTARRAY *pa, const POINT4D *p4d, double *mindistout, POINT4D *proj4d)
1654 {
1655  double mindist=DBL_MAX;
1656  double tlen, plen;
1657  uint32_t t, seg=0;
1658  POINT4D start4d, end4d, projtmp;
1659  POINT2D proj, p;
1660  const POINT2D *start = NULL, *end = NULL;
1661 
1662  /* Initialize our 2D copy of the input parameter */
1663  p.x = p4d->x;
1664  p.y = p4d->y;
1665 
1666  if ( ! proj4d ) proj4d = &projtmp;
1667 
1668  /* Check for special cases (length 0 and 1) */
1669  if ( pa->npoints <= 1 )
1670  {
1671  if ( pa->npoints == 1 )
1672  {
1673  getPoint4d_p(pa, 0, proj4d);
1674  if ( mindistout )
1675  *mindistout = distance2d_pt_pt(&p, getPoint2d_cp(pa, 0));
1676  }
1677  return 0.0;
1678  }
1679 
1680  start = getPoint2d_cp(pa, 0);
1681  /* Loop through pointarray looking for nearest segment */
1682  for (t=1; t<pa->npoints; t++)
1683  {
1684  double dist_sqr;
1685  end = getPoint2d_cp(pa, t);
1686  dist_sqr = distance2d_sqr_pt_seg(&p, start, end);
1687 
1688  if (dist_sqr < mindist)
1689  {
1690  mindist = dist_sqr;
1691  seg=t-1;
1692  if ( mindist == 0 )
1693  {
1694  LWDEBUG(3, "Breaking on mindist=0");
1695  break;
1696  }
1697  }
1698 
1699  start = end;
1700  }
1701  mindist = sqrt(mindist);
1702 
1703  if ( mindistout ) *mindistout = mindist;
1704 
1705  LWDEBUGF(3, "Closest segment: %d", seg);
1706  LWDEBUGF(3, "mindist: %g", mindist);
1707 
1708  /*
1709  * We need to project the
1710  * point on the closest segment.
1711  */
1712  getPoint4d_p(pa, seg, &start4d);
1713  getPoint4d_p(pa, seg+1, &end4d);
1714  closest_point_on_segment(p4d, &start4d, &end4d, proj4d);
1715 
1716  /* Copy 4D values into 2D holder */
1717  proj.x = proj4d->x;
1718  proj.y = proj4d->y;
1719 
1720  LWDEBUGF(3, "Closest segment:%d, npoints:%d", seg, pa->npoints);
1721 
1722  /* For robustness, force 1 when closest point == endpoint */
1723  if ( (seg >= (pa->npoints-2)) && p2d_same(&proj, end) )
1724  {
1725  return 1.0;
1726  }
1727 
1728  LWDEBUGF(3, "Closest point on segment: %g,%g", proj.x, proj.y);
1729 
1730  tlen = ptarray_length_2d(pa);
1731 
1732  LWDEBUGF(3, "tlen %g", tlen);
1733 
1734  /* Location of any point on a zero-length line is 0 */
1735  /* See http://trac.osgeo.org/postgis/ticket/1772#comment:2 */
1736  if ( tlen == 0 ) return 0;
1737 
1738  plen=0;
1739  start = getPoint2d_cp(pa, 0);
1740  for (t=0; t<seg; t++, start=end)
1741  {
1742  end = getPoint2d_cp(pa, t+1);
1743  plen += distance2d_pt_pt(start, end);
1744 
1745  LWDEBUGF(4, "Segment %d made plen %g", t, plen);
1746  }
1747 
1748  plen+=distance2d_pt_pt(&proj, start);
1749 
1750  LWDEBUGF(3, "plen %g, tlen %g", plen, tlen);
1751 
1752  return plen/tlen;
1753 }
1754 
1764 void
1766 {
1767  uint32_t i;
1768  double x;
1769 
1770  for (i=0; i<pa->npoints; i++)
1771  {
1772  memcpy(&x, getPoint_internal(pa, i), sizeof(double));
1773  if ( x < 0 ) x+= 360;
1774  else if ( x > 180 ) x -= 360;
1775  memcpy(getPoint_internal(pa, i), &x, sizeof(double));
1776  }
1777 }
1778 
1779 
1780 /*
1781  * Returns a POINTARRAY with consecutive equal points
1782  * removed. Equality test on all dimensions of input.
1783  *
1784  * Always returns a newly allocated object.
1785  */
1786 static POINTARRAY *
1787 ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
1788 {
1789  POINTARRAY *out = ptarray_clone_deep(in);
1790  ptarray_remove_repeated_points_in_place(out, tolerance, minpoints);
1791  return out;
1792 }
1793 
1794 POINTARRAY *
1795 ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
1796 {
1797  return ptarray_remove_repeated_points_minpoints(in, tolerance, 2);
1798 }
1799 
1800 
1801 void
1802 ptarray_remove_repeated_points_in_place(POINTARRAY *pa, double tolerance, uint32_t min_points)
1803 {
1804  uint32_t i;
1805  double tolsq = tolerance * tolerance;
1806  const POINT2D *last = NULL;
1807  const POINT2D *pt;
1808  uint32_t n_points = pa->npoints;
1809  uint32_t n_points_out = 1;
1810  size_t pt_size = ptarray_point_size(pa);
1811 
1812  double dsq = FLT_MAX;
1813 
1814  /* No-op on short inputs */
1815  if ( n_points <= min_points ) return;
1816 
1817  last = getPoint2d_cp(pa, 0);
1818  void *p_to = ((char *)last) + pt_size;
1819  for (i = 1; i < n_points; i++)
1820  {
1821  int last_point = (i == n_points - 1);
1822 
1823  /* Look straight into the abyss */
1824  pt = getPoint2d_cp(pa, i);
1825 
1826  /* Don't drop points if we are running short of points */
1827  if (n_points + n_points_out > min_points + i)
1828  {
1829  if (tolerance > 0.0)
1830  {
1831  /* Only drop points that are within our tolerance */
1832  dsq = distance2d_sqr_pt_pt(last, pt);
1833  /* Allow any point but the last one to be dropped */
1834  if (!last_point && dsq <= tolsq)
1835  {
1836  continue;
1837  }
1838  }
1839  else
1840  {
1841  /* At tolerance zero, only skip exact dupes */
1842  if (memcmp((char*)pt, (char*)last, pt_size) == 0)
1843  continue;
1844  }
1845 
1846  /* Got to last point, and it's not very different from */
1847  /* the point that preceded it. We want to keep the last */
1848  /* point, not the second-to-last one, so we pull our write */
1849  /* index back one value */
1850  if (last_point && n_points_out > 1 && tolerance > 0.0 && dsq <= tolsq)
1851  {
1852  n_points_out--;
1853  p_to = (char*)p_to - pt_size;
1854  }
1855  }
1856 
1857  /* Compact all remaining values to front of array */
1858  memcpy(p_to, pt, pt_size);
1859  n_points_out++;
1860  p_to = (char*)p_to + pt_size;
1861  last = pt;
1862  }
1863  /* Adjust array length */
1864  pa->npoints = n_points_out;
1865  return;
1866 }
1867 
1868 /* Out of the points in pa [itfist .. itlast], finds the one that's farthest away from
1869  * the segment determined by pts[itfist] and pts[itlast].
1870  * Returns itfirst if no point was found futher away than max_distance_sqr
1871  */
1872 static uint32_t
1873 ptarray_dp_findsplit_in_place(const POINTARRAY *pts, uint32_t it_first, uint32_t it_last, double max_distance_sqr)
1874 {
1875  uint32_t split = it_first;
1876  if ((it_first - it_last) < 2)
1877  return it_first;
1878 
1879  const POINT2D *A = getPoint2d_cp(pts, it_first);
1880  const POINT2D *B = getPoint2d_cp(pts, it_last);
1881 
1882  if (distance2d_sqr_pt_pt(A, B) < DBL_EPSILON)
1883  {
1884  /* If p1 == p2, we can just calculate the distance from each point to A */
1885  for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1886  {
1887  const POINT2D *pk = getPoint2d_cp(pts, itk);
1888  double distance_sqr = distance2d_sqr_pt_pt(pk, A);
1889  if (distance_sqr > max_distance_sqr)
1890  {
1891  split = itk;
1892  max_distance_sqr = distance_sqr;
1893  }
1894  }
1895  return split;
1896  }
1897 
1898  /* This is based on distance2d_sqr_pt_seg, but heavily inlined here to avoid recalculations */
1899  double ba_x = (B->x - A->x);
1900  double ba_y = (B->y - A->y);
1901  double ab_length_sqr = (ba_x * ba_x + ba_y * ba_y);
1902  /* To avoid the division by ab_length_sqr in the 3rd path, we normalize here
1903  * and multiply in the first two paths [(dot_ac_ab < 0) and (> ab_length_sqr)] */
1904  max_distance_sqr *= ab_length_sqr;
1905  for (uint32_t itk = it_first + 1; itk < it_last; itk++)
1906  {
1907  const POINT2D *C = getPoint2d_cp(pts, itk);
1908  double distance_sqr;
1909  double ca_x = (C->x - A->x);
1910  double ca_y = (C->y - A->y);
1911  double dot_ac_ab = (ca_x * ba_x + ca_y * ba_y);
1912 
1913  if (dot_ac_ab <= 0.0)
1914  {
1915  distance_sqr = distance2d_sqr_pt_pt(C, A) * ab_length_sqr;
1916  }
1917  else if (dot_ac_ab >= ab_length_sqr)
1918  {
1919  distance_sqr = distance2d_sqr_pt_pt(C, B) * ab_length_sqr;
1920  }
1921  else
1922  {
1923  double s_numerator = ca_x * ba_y - ca_y * ba_x;
1924  distance_sqr = s_numerator * s_numerator; /* Missing division by ab_length_sqr on purpose */
1925  }
1926 
1927  if (distance_sqr > max_distance_sqr)
1928  {
1929  split = itk;
1930  max_distance_sqr = distance_sqr;
1931  }
1932  }
1933  return split;
1934 }
1935 
1936 /* O(N) simplification for tolearnce = 0 */
1937 static void
1939 {
1940  uint32_t kept_it = 0;
1941  uint32_t last_it = pa->npoints - 1;
1942  const POINT2D *kept_pt = getPoint2d_cp(pa, 0);
1943  const size_t pt_size = ptarray_point_size(pa);
1944 
1945  for (uint32_t i = 1; i < last_it; i++)
1946  {
1947  const POINT2D *curr_pt = getPoint2d_cp(pa, i);
1948  const POINT2D *next_pt = getPoint2d_cp(pa, i + 1);
1949 
1950  double ba_x = next_pt->x - kept_pt->x;
1951  double ba_y = next_pt->y - kept_pt->y;
1952  double ab_length_sqr = ba_x * ba_x + ba_y * ba_y;
1953 
1954  double ca_x = curr_pt->x - kept_pt->x;
1955  double ca_y = curr_pt->y - kept_pt->y;
1956  double dot_ac_ab = ca_x * ba_x + ca_y * ba_y;
1957  double s_numerator = ca_x * ba_y - ca_y * ba_x;
1958 
1959  if (dot_ac_ab < 0.0 || dot_ac_ab > ab_length_sqr || s_numerator != 0)
1960  {
1961  kept_it++;
1962  kept_pt = curr_pt;
1963  if (kept_it != i)
1964  memcpy(pa->serialized_pointlist + pt_size * kept_it,
1965  pa->serialized_pointlist + pt_size * i,
1966  pt_size);
1967  }
1968  }
1969 
1970  /* Append last point */
1971  kept_it++;
1972  if (kept_it != last_it)
1973  memcpy(pa->serialized_pointlist + pt_size * kept_it,
1974  pa->serialized_pointlist + pt_size * last_it,
1975  pt_size);
1976  pa->npoints = kept_it + 1;
1977 }
1978 
1979 void
1980 ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
1981 {
1982  /* Do not try to simplify really short things */
1983  if (pa->npoints < 3 || pa->npoints <= minpts)
1984  return;
1985 
1986  if (tolerance == 0 && minpts <= 2)
1987  {
1989  return;
1990  }
1991 
1992  /* We use this array to keep track of the points we are keeping, so
1993  * we store just TRUE / FALSE in their position */
1994  uint8_t *kept_points = lwalloc(sizeof(uint8_t) * pa->npoints);
1995  memset(kept_points, LW_FALSE, sizeof(uint8_t) * pa->npoints);
1996  kept_points[0] = LW_TRUE;
1997  kept_points[pa->npoints - 1] = LW_TRUE;
1998  uint32_t keptn = 2;
1999 
2000  /* We use this array as a stack to store the iterators that we are going to need
2001  * in the following steps.
2002  * This is ~10% faster than iterating over @kept_points looking for them
2003  */
2004  uint32_t *iterator_stack = lwalloc(sizeof(uint32_t) * pa->npoints);
2005  iterator_stack[0] = 0;
2006  uint32_t iterator_stack_size = 1;
2007 
2008  uint32_t it_first = 0;
2009  uint32_t it_last = pa->npoints - 1;
2010 
2011  const double tolerance_sqr = tolerance * tolerance;
2012  /* For the first @minpts points we ignore the tolerance */
2013  double it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
2014 
2015  while (iterator_stack_size)
2016  {
2017  uint32_t split = ptarray_dp_findsplit_in_place(pa, it_first, it_last, it_tol);
2018  if (split == it_first)
2019  {
2020  it_first = it_last;
2021  it_last = iterator_stack[--iterator_stack_size];
2022  }
2023  else
2024  {
2025  kept_points[split] = LW_TRUE;
2026  keptn++;
2027 
2028  iterator_stack[iterator_stack_size++] = it_last;
2029  it_last = split;
2030  it_tol = keptn >= minpts ? tolerance_sqr : -1.0;
2031  }
2032  }
2033 
2034  const size_t pt_size = ptarray_point_size(pa);
2035  /* The first point is already in place, so we don't need to copy it */
2036  size_t kept_it = 1;
2037  if (keptn == 2)
2038  {
2039  /* If there are 2 points remaining, it has to be first and last as
2040  * we added those at the start */
2041  memcpy(pa->serialized_pointlist + pt_size * kept_it,
2042  pa->serialized_pointlist + pt_size * (pa->npoints - 1),
2043  pt_size);
2044  }
2045  else if (pa->npoints != keptn) /* We don't need to move any points if we are keeping them all */
2046  {
2047  for (uint32_t i = 1; i < pa->npoints; i++)
2048  {
2049  if (kept_points[i])
2050  {
2051  memcpy(pa->serialized_pointlist + pt_size * kept_it,
2052  pa->serialized_pointlist + pt_size * i,
2053  pt_size);
2054  kept_it++;
2055  }
2056  }
2057  }
2058  pa->npoints = keptn;
2059 
2060  lwfree(kept_points);
2061  lwfree(iterator_stack);
2062 }
2063 
2064 /************************************************************************/
2065 
2071 double
2073 {
2074  double dist = 0.0;
2075  uint32_t i;
2076  const POINT2D *a1;
2077  const POINT2D *a2;
2078  const POINT2D *a3;
2079 
2080  if ( pts->npoints % 2 != 1 )
2081  lwerror("arc point array with even number of points");
2082 
2083  a1 = getPoint2d_cp(pts, 0);
2084 
2085  for ( i=2; i < pts->npoints; i += 2 )
2086  {
2087  a2 = getPoint2d_cp(pts, i-1);
2088  a3 = getPoint2d_cp(pts, i);
2089  dist += lw_arc_length(a1, a2, a3);
2090  a1 = a3;
2091  }
2092  return dist;
2093 }
2094 
2098 double
2100 {
2101  double dist = 0.0;
2102  uint32_t i;
2103  const POINT2D *frm;
2104  const POINT2D *to;
2105 
2106  if ( pts->npoints < 2 ) return 0.0;
2107 
2108  frm = getPoint2d_cp(pts, 0);
2109 
2110  for ( i=1; i < pts->npoints; i++ )
2111  {
2112  to = getPoint2d_cp(pts, i);
2113 
2114  dist += sqrt( ((frm->x - to->x)*(frm->x - to->x)) +
2115  ((frm->y - to->y)*(frm->y - to->y)) );
2116 
2117  frm = to;
2118  }
2119  return dist;
2120 }
2121 
2126 double
2128 {
2129  double dist = 0.0;
2130  uint32_t i;
2131  POINT3DZ frm;
2132  POINT3DZ to;
2133 
2134  if ( pts->npoints < 2 ) return 0.0;
2135 
2136  /* compute 2d length if 3d is not available */
2137  if ( ! FLAGS_GET_Z(pts->flags) ) return ptarray_length_2d(pts);
2138 
2139  getPoint3dz_p(pts, 0, &frm);
2140  for ( i=1; i < pts->npoints; i++ )
2141  {
2142  getPoint3dz_p(pts, i, &to);
2143  dist += sqrt( ((frm.x - to.x)*(frm.x - to.x)) +
2144  ((frm.y - to.y)*(frm.y - to.y)) +
2145  ((frm.z - to.z)*(frm.z - to.z)) );
2146  frm = to;
2147  }
2148  return dist;
2149 }
2150 
2151 
2152 
2156 void
2158 {
2159  if (FLAGS_GET_Z(pa->flags))
2160  {
2161  for (uint32_t i = 0; i < pa->npoints; i++)
2162  {
2163  POINT4D *p4d = (POINT4D *)(getPoint_internal(pa, i));
2164  double x = p4d->x;
2165  double y = p4d->y;
2166  double z = p4d->z;
2167  p4d->x = a->afac * x + a->bfac * y + a->cfac * z + a->xoff;
2168  p4d->y = a->dfac * x + a->efac * y + a->ffac * z + a->yoff;
2169  p4d->z = a->gfac * x + a->hfac * y + a->ifac * z + a->zoff;
2170  }
2171  }
2172  else
2173  {
2174  for (uint32_t i = 0; i < pa->npoints; i++)
2175  {
2176  POINT2D *pt = (POINT2D *)(getPoint_internal(pa, i));
2177  double x = pt->x;
2178  double y = pt->y;
2179  pt->x = a->afac * x + a->bfac * y + a->xoff;
2180  pt->y = a->dfac * x + a->efac * y + a->yoff;
2181  }
2182  }
2183 }
2184 
2189 #if 0
2190 static int gluInvertMatrix(const double *m, double *invOut)
2191 {
2192  double inv[16], det;
2193  int i;
2194 
2195  inv[0] = m[5] * m[10] * m[15] -
2196  m[5] * m[11] * m[14] -
2197  m[9] * m[6] * m[15] +
2198  m[9] * m[7] * m[14] +
2199  m[13] * m[6] * m[11] -
2200  m[13] * m[7] * m[10];
2201 
2202  inv[4] = -m[4] * m[10] * m[15] +
2203  m[4] * m[11] * m[14] +
2204  m[8] * m[6] * m[15] -
2205  m[8] * m[7] * m[14] -
2206  m[12] * m[6] * m[11] +
2207  m[12] * m[7] * m[10];
2208 
2209  inv[8] = m[4] * m[9] * m[15] -
2210  m[4] * m[11] * m[13] -
2211  m[8] * m[5] * m[15] +
2212  m[8] * m[7] * m[13] +
2213  m[12] * m[5] * m[11] -
2214  m[12] * m[7] * m[9];
2215 
2216  inv[12] = -m[4] * m[9] * m[14] +
2217  m[4] * m[10] * m[13] +
2218  m[8] * m[5] * m[14] -
2219  m[8] * m[6] * m[13] -
2220  m[12] * m[5] * m[10] +
2221  m[12] * m[6] * m[9];
2222 
2223  inv[1] = -m[1] * m[10] * m[15] +
2224  m[1] * m[11] * m[14] +
2225  m[9] * m[2] * m[15] -
2226  m[9] * m[3] * m[14] -
2227  m[13] * m[2] * m[11] +
2228  m[13] * m[3] * m[10];
2229 
2230  inv[5] = m[0] * m[10] * m[15] -
2231  m[0] * m[11] * m[14] -
2232  m[8] * m[2] * m[15] +
2233  m[8] * m[3] * m[14] +
2234  m[12] * m[2] * m[11] -
2235  m[12] * m[3] * m[10];
2236 
2237  inv[9] = -m[0] * m[9] * m[15] +
2238  m[0] * m[11] * m[13] +
2239  m[8] * m[1] * m[15] -
2240  m[8] * m[3] * m[13] -
2241  m[12] * m[1] * m[11] +
2242  m[12] * m[3] * m[9];
2243 
2244  inv[13] = m[0] * m[9] * m[14] -
2245  m[0] * m[10] * m[13] -
2246  m[8] * m[1] * m[14] +
2247  m[8] * m[2] * m[13] +
2248  m[12] * m[1] * m[10] -
2249  m[12] * m[2] * m[9];
2250 
2251  inv[2] = m[1] * m[6] * m[15] -
2252  m[1] * m[7] * m[14] -
2253  m[5] * m[2] * m[15] +
2254  m[5] * m[3] * m[14] +
2255  m[13] * m[2] * m[7] -
2256  m[13] * m[3] * m[6];
2257 
2258  inv[6] = -m[0] * m[6] * m[15] +
2259  m[0] * m[7] * m[14] +
2260  m[4] * m[2] * m[15] -
2261  m[4] * m[3] * m[14] -
2262  m[12] * m[2] * m[7] +
2263  m[12] * m[3] * m[6];
2264 
2265  inv[10] = m[0] * m[5] * m[15] -
2266  m[0] * m[7] * m[13] -
2267  m[4] * m[1] * m[15] +
2268  m[4] * m[3] * m[13] +
2269  m[12] * m[1] * m[7] -
2270  m[12] * m[3] * m[5];
2271 
2272  inv[14] = -m[0] * m[5] * m[14] +
2273  m[0] * m[6] * m[13] +
2274  m[4] * m[1] * m[14] -
2275  m[4] * m[2] * m[13] -
2276  m[12] * m[1] * m[6] +
2277  m[12] * m[2] * m[5];
2278 
2279  inv[3] = -m[1] * m[6] * m[11] +
2280  m[1] * m[7] * m[10] +
2281  m[5] * m[2] * m[11] -
2282  m[5] * m[3] * m[10] -
2283  m[9] * m[2] * m[7] +
2284  m[9] * m[3] * m[6];
2285 
2286  inv[7] = m[0] * m[6] * m[11] -
2287  m[0] * m[7] * m[10] -
2288  m[4] * m[2] * m[11] +
2289  m[4] * m[3] * m[10] +
2290  m[8] * m[2] * m[7] -
2291  m[8] * m[3] * m[6];
2292 
2293  inv[11] = -m[0] * m[5] * m[11] +
2294  m[0] * m[7] * m[9] +
2295  m[4] * m[1] * m[11] -
2296  m[4] * m[3] * m[9] -
2297  m[8] * m[1] * m[7] +
2298  m[8] * m[3] * m[5];
2299 
2300  inv[15] = m[0] * m[5] * m[10] -
2301  m[0] * m[6] * m[9] -
2302  m[4] * m[1] * m[10] +
2303  m[4] * m[2] * m[9] +
2304  m[8] * m[1] * m[6] -
2305  m[8] * m[2] * m[5];
2306 
2307  det = m[0] * inv[0] + m[1] * inv[4] + m[2] * inv[8] + m[3] * inv[12];
2308 
2309  if (det == 0)
2310  return LW_FALSE;
2311 
2312  det = 1.0 / det;
2313 
2314  for (i = 0; i < 16; i++)
2315  invOut[i] = inv[i] * det;
2316 
2317  return LW_TRUE;
2318 }
2319 #endif
2320 
2324 void
2326 {
2327  uint32_t i;
2328  POINT4D p4d;
2329  LWDEBUG(3, "ptarray_scale start");
2330  for (i=0; i<pa->npoints; i++)
2331  {
2332  getPoint4d_p(pa, i, &p4d);
2333  p4d.x *= fact->x;
2334  p4d.y *= fact->y;
2335  p4d.z *= fact->z;
2336  p4d.m *= fact->m;
2337  ptarray_set_point4d(pa, i, &p4d);
2338  }
2339  LWDEBUG(3, "ptarray_scale end");
2340 }
2341 
2342 int
2344 {
2345  return getPoint4d_p(pa, 0, pt);
2346 }
2347 
2348 
2349 /*
2350  * Stick an array of points to the given gridspec.
2351  * Return "gridded" points in *outpts and their number in *outptsn.
2352  *
2353  * Two consecutive points falling on the same grid cell are collapsed
2354  * into one single point.
2355  *
2356  */
2357 void
2359 {
2360  uint32_t j = 0;
2361  POINT4D *p, *p_out = NULL;
2362  double x, y, z = 0, m = 0;
2363  uint32_t ndims = FLAGS_NDIMS(pa->flags);
2364  uint32_t has_z = FLAGS_GET_Z(pa->flags);
2365  uint32_t has_m = FLAGS_GET_M(pa->flags);
2366 
2367  for (uint32_t i = 0; i < pa->npoints; i++)
2368  {
2369  /* Look straight into the abyss */
2370  p = (POINT4D *)(getPoint_internal(pa, i));
2371  x = p->x;
2372  y = p->y;
2373  if (ndims > 2)
2374  z = p->z;
2375  if (ndims > 3)
2376  m = p->m;
2377 
2378  if (grid->xsize > 0)
2379  x = rint((x - grid->ipx) / grid->xsize) * grid->xsize + grid->ipx;
2380 
2381  if (grid->ysize > 0)
2382  y = rint((y - grid->ipy) / grid->ysize) * grid->ysize + grid->ipy;
2383 
2384  /* Read and round this point */
2385  /* Z is always in third position */
2386  if (has_z && grid->zsize > 0)
2387  z = rint((z - grid->ipz) / grid->zsize) * grid->zsize + grid->ipz;
2388 
2389  /* M might be in 3rd or 4th position */
2390  if (has_m && grid->msize > 0)
2391  {
2392  /* In POINT ZM, M is in 4th position, in POINT M, M is in 3rd position which is Z in POINT4D */
2393  if (has_z)
2394  m = rint((m - grid->ipm) / grid->msize) * grid->msize + grid->ipm;
2395  else
2396  z = rint((z - grid->ipm) / grid->msize) * grid->msize + grid->ipm;
2397  }
2398 
2399  /* Skip duplicates */
2400  if (p_out && p_out->x == x && p_out->y == y && (ndims > 2 ? p_out->z == z : 1) &&
2401  (ndims > 3 ? p_out->m == m : 1))
2402  continue;
2403 
2404  /* Write rounded values into the next available point */
2405  p_out = (POINT4D *)(getPoint_internal(pa, j++));
2406  p_out->x = x;
2407  p_out->y = y;
2408  if (ndims > 2)
2409  p_out->z = z;
2410  if (ndims > 3)
2411  p_out->m = m;
2412  }
2413 
2414  /* Update output ptarray length */
2415  pa->npoints = j;
2416  return;
2417 }
2418 
2419 int
2420 ptarray_npoints_in_rect(const POINTARRAY *pa, const GBOX *gbox)
2421 {
2422  const POINT2D *pt;
2423  int n = 0;
2424  uint32_t i;
2425  for ( i = 0; i < pa->npoints; i++ )
2426  {
2427  pt = getPoint2d_cp(pa, i);
2428  if ( gbox_contains_point2d(gbox, pt) )
2429  n++;
2430  }
2431  return n;
2432 }
2433 
2434 
2435 /*
2436  * Reorder the vertices of a closed pointarray so that the
2437  * given point is the first/last one.
2438  *
2439  * Error out if pointarray is not closed or it does not
2440  * contain the given point.
2441  */
2442 int
2444 {
2445  POINTARRAY *tmp;
2446  int found;
2447  uint32_t it;
2448  int ptsize;
2449 
2450  if ( ! ptarray_is_closed_2d(pa) )
2451  {
2452  lwerror("ptarray_scroll_in_place: input POINTARRAY is not closed");
2453  return LW_FAILURE;
2454  }
2455 
2456  ptsize = ptarray_point_size(pa);
2457 
2458  /* Find the point in the array */
2459  found = 0;
2460  for ( it = 0; it < pa->npoints; ++it )
2461  {
2462  if ( ! memcmp(getPoint_internal(pa, it), pt, ptsize) )
2463  {
2464  found = 1;
2465  break;
2466  }
2467  }
2468 
2469  if ( ! found )
2470  {
2471  lwerror("ptarray_scroll_in_place: input POINTARRAY does not contain the given point");
2472  return LW_FAILURE;
2473  }
2474 
2475  if ( 0 == it )
2476  {
2477  /* Point is already the start/end point, just clone the input */
2478  return LW_SUCCESS;
2479  }
2480 
2481  /* TODO: reduce allocations */
2483 
2484  bzero(getPoint_internal(tmp, 0), ptsize * pa->npoints);
2485  /* Copy the block from found point to last point into the output array */
2486  memcpy(
2487  getPoint_internal(tmp, 0),
2488  getPoint_internal(pa, it),
2489  ptsize * ( pa->npoints - it )
2490  );
2491 
2492  /* Copy the block from second point to the found point into the last portion of the
2493  * return */
2494  memcpy(
2495  getPoint_internal(tmp, pa->npoints - it),
2496  getPoint_internal(pa, 1),
2497  ptsize * ( it )
2498  );
2499 
2500  /* Copy the resulting pointarray back to source one */
2501  memcpy(
2502  getPoint_internal(pa, 0),
2503  getPoint_internal(tmp, 0),
2504  ptsize * ( pa->npoints )
2505  );
2506 
2507  ptarray_free(tmp);
2508 
2509  return LW_SUCCESS;
2510 }
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
double distance2d_pt_pt(const POINT2D *p1, const POINT2D *p2)
Definition: measures.c:2306
#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:650
#define LW_SUCCESS
Definition: liblwgeom.h:111
int getPoint2d_p(const POINTARRAY *pa, uint32_t n, POINT2D *point)
Definition: lwgeom_api.c:343
#define FLAGS_GET_READONLY(flags)
Definition: liblwgeom.h:183
#define FLAGS_GET_Z(flags)
Definition: liblwgeom.h:179
int getPoint3dz_p(const POINTARRAY *pa, uint32_t n, POINT3DZ *point)
Definition: lwgeom_api.c:216
void * lwrealloc(void *mem, size_t size)
Definition: lwutil.c:235
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:126
double distance2d_sqr_pt_seg(const POINT2D *p, const POINT2D *A, const POINT2D *B)
Definition: measures.c:2316
#define FLAGS_GET_ZM(flags)
Definition: liblwgeom.h:194
void * lwalloc(size_t size)
Definition: lwutil.c:227
#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:370
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.
Definition: lwalgorithm.c:127
int lw_pt_on_segment(const POINT2D *A1, const POINT2D *A2, const POINT2D *P)
Definition: lwalgorithm.c:103
#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.
Definition: lwalgorithm.c:237
int lw_arc_side(const POINT2D *A1, const POINT2D *A2, const POINT2D *A3, const POINT2D *Q)
Definition: lwalgorithm.c:187
#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:63
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) .
Definition: lwalgorithm.c:114
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:84
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
Definition: lwin_wkt_lex.c:334
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:101
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:58
static uint8_t * getPoint_internal(const POINTARRAY *pa, uint32_t n)
Definition: lwinline.h:77
int ptarray_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
Definition: ptarray.c:1017
POINTARRAY * ptarray_force_dims(const POINTARRAY *pa, int hasz, int hasm, double zval, double mval)
Definition: ptarray.c:1323
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:1873
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_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:1765
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_closest_segment_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
Definition: ptarray.c:1587
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 ptarray_closest_vertex_2d(const POINTARRAY *pa, const POINT2D *qp, double *dist)
Definition: ptarray.c:1620
int ptarrayarc_contains_point(const POINTARRAY *pa, const POINT2D *pt)
For POINTARRAYs representing CIRCULARSTRINGS.
Definition: ptarray.c:1110
int ptarray_is_closed_2d(const POINTARRAY *in)
Definition: ptarray.c:701
int ptarray_is_closed_z(const POINTARRAY *in)
Definition: ptarray.c:727
int ptarray_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
Definition: ptarray.c:818
double ptarray_length(const POINTARRAY *pts)
Find the 3d/2d length of the given POINTARRAY (depending on its dimensionality)
Definition: ptarray.c:2127
double ptarray_signed_area(const POINTARRAY *pa)
Returns the area in cartesian units.
Definition: ptarray.c:1283
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 closest_point_on_segment(const POINT4D *p, const POINT4D *A, const POINT4D *B, POINT4D *ret)
Definition: ptarray.c:1543
int ptarray_startpoint(const POINTARRAY *pa, POINT4D *pt)
Definition: ptarray.c:2343
int ptarray_is_closed(const POINTARRAY *in)
Check for ring closure using whatever dimensionality is declared on the pointarray.
Definition: ptarray.c:687
POINTARRAY * ptarray_clone_deep(const POINTARRAY *in)
Deep clone a pointarray (also clones serialized pointlist)
Definition: ptarray.c:634
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
void ptarray_grid_in_place(POINTARRAY *pa, const gridspec *grid)
Snap to grid.
Definition: ptarray.c:2358
double ptarray_length_2d(const POINTARRAY *pts)
Find the 2d length of the given POINTARRAY (even if it's 3d)
Definition: ptarray.c:2099
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:1802
int ptarrayarc_raycast_intersections(const POINTARRAY *pa, const POINT2D *p, int *on_boundary)
Definition: ptarray.c:935
POINTARRAY * ptarray_removePoint(POINTARRAY *pa, uint32_t which)
Remove a point from a pointarray.
Definition: ptarray.c:561
int ptarray_scroll_in_place(POINTARRAY *pa, const POINT4D *pt)
Definition: ptarray.c:2443
void ptarray_simplify_in_place(POINTARRAY *pa, double tolerance, uint32_t minpts)
Definition: ptarray.c:1980
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:2072
POINTARRAY * ptarray_remove_repeated_points(const POINTARRAY *in, double tolerance)
Definition: ptarray.c:1795
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:894
void ptarray_affine(POINTARRAY *pa, const AFFINE *a)
Affine transform a pointarray.
Definition: ptarray.c:2157
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:746
int ptarrayarc_contains_point_partial(const POINTARRAY *pa, const POINT2D *pt, int check_closed, int *winding_number)
Definition: ptarray.c:1125
POINTARRAY * ptarray_construct_empty(char hasz, char hasm, uint32_t maxpoints)
Create a new POINTARRAY with no points.
Definition: ptarray.c:59
void ptarray_free(POINTARRAY *pa)
Definition: ptarray.c:327
static void ptarray_simplify_in_place_tolerance0(POINTARRAY *pa)
Definition: ptarray.c:1938
int ptarray_isccw(const POINTARRAY *pa)
Definition: ptarray.c:1314
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
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:2325
static POINTARRAY * ptarray_remove_repeated_points_minpoints(const POINTARRAY *in, double tolerance, int minpoints)
Definition: ptarray.c:1787
void ptarray_swap_ordinates(POINTARRAY *pa, LWORD o1, LWORD o2)
Swap ordinate values o1 and o2 on a given POINTARRAY.
Definition: ptarray.c:387
int ptarray_npoints_in_rect(const POINTARRAY *pa, const GBOX *gbox)
Definition: ptarray.c:2420
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:1346
POINTARRAY * ptarray_merge(POINTARRAY *pa1, POINTARRAY *pa2)
Merge two given POINTARRAY and returns a pointer on the new aggregate one.
Definition: ptarray.c:603
POINTARRAY * ptarray_flip_coordinates(POINTARRAY *pa)
Reverse X and Y axis on a given POINTARRAY.
Definition: ptarray.c:368
int ptarray_has_m(const POINTARRAY *pa)
Definition: ptarray.c:44
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
double ptarray_locate_point(const POINTARRAY *pa, const POINT4D *p4d, double *mindistout, POINT4D *proj4d)
Definition: ptarray.c:1653
double gfac
Definition: liblwgeom.h:346
double zoff
Definition: liblwgeom.h:346
double bfac
Definition: liblwgeom.h:346
double ifac
Definition: liblwgeom.h:346
double xoff
Definition: liblwgeom.h:346
double dfac
Definition: liblwgeom.h:346
double afac
Definition: liblwgeom.h:346
double ffac
Definition: liblwgeom.h:346
double cfac
Definition: liblwgeom.h:346
double hfac
Definition: liblwgeom.h:346
double efac
Definition: liblwgeom.h:346
double yoff
Definition: liblwgeom.h:346
double ymax
Definition: liblwgeom.h:371
double xmax
Definition: liblwgeom.h:369
double ymin
Definition: liblwgeom.h:370
double xmin
Definition: liblwgeom.h:368
double y
Definition: liblwgeom.h:404
double x
Definition: liblwgeom.h:404
double z
Definition: liblwgeom.h:410
double x
Definition: liblwgeom.h:410
double y
Definition: liblwgeom.h:410
double m
Definition: liblwgeom.h:428
double x
Definition: liblwgeom.h:428
double z
Definition: liblwgeom.h:428
double y
Definition: liblwgeom.h:428
lwflags_t flags
Definition: liblwgeom.h:445
uint32_t maxpoints
Definition: liblwgeom.h:442
uint32_t npoints
Definition: liblwgeom.h:441
uint8_t * serialized_pointlist
Definition: liblwgeom.h:448
double ipm
Definition: liblwgeom.h:1386
double zsize
Definition: liblwgeom.h:1389
double ysize
Definition: liblwgeom.h:1388
double xsize
Definition: liblwgeom.h:1387
double ipx
Definition: liblwgeom.h:1383
double msize
Definition: liblwgeom.h:1390
double ipy
Definition: liblwgeom.h:1384
double ipz
Definition: liblwgeom.h:1385
Snap-to-grid.
Definition: liblwgeom.h:1382