PostGIS  2.4.9dev-r@@SVN_REVISION@@
rt_spatial_relationship.c
Go to the documentation of this file.
1 /*
2  *
3  * WKTRaster - Raster Types for PostGIS
4  * http://trac.osgeo.org/postgis/wiki/WKTRaster
5  *
6  * Copyright (C) 2011-2013 Regents of the University of California
7  * <bkpark@ucdavis.edu>
8  * Copyright (C) 2010-2011 Jorge Arevalo <jorge.arevalo@deimos-space.com>
9  * Copyright (C) 2010-2011 David Zwarg <dzwarg@azavea.com>
10  * Copyright (C) 2009-2011 Pierre Racine <pierre.racine@sbf.ulaval.ca>
11  * Copyright (C) 2009-2011 Mateusz Loskot <mateusz@loskot.net>
12  * Copyright (C) 2008-2009 Sandro Santilli <strk@kbt.io>
13  *
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * This program is distributed in the hope that it will be useful,
20  * but WITHOUT ANY WARRANTY; without even the implied warranty of
21  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
22  * GNU General Public License for more details.
23  *
24  * You should have received a copy of the GNU General Public License
25  * along with this program; if not, write to the Free Software Foundation,
26  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
27  *
28  */
29 
30 #include "librtcore.h"
31 #include "librtcore_internal.h"
32 
33 /*
34  * Return ES_ERROR if error occurred in function.
35  * Paramter aligned returns non-zero if two rasters are aligned
36  *
37  * @param rast1 : the first raster for alignment test
38  * @param rast2 : the second raster for alignment test
39  * @param *aligned : non-zero value if the two rasters are aligned
40  * @param *reason : reason why rasters are not aligned
41  *
42  * @return ES_NONE if success, ES_ERROR if error
43  */
46  rt_raster rast1,
47  rt_raster rast2,
48  int *aligned, char **reason
49 ) {
50  double xr;
51  double yr;
52  double xw;
53  double yw;
54  int err = 0;
55 
56  assert(NULL != rast1);
57  assert(NULL != rast2);
58  assert(NULL != aligned);
59 
60  err = 0;
61  /* same srid */
62  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
63  if (reason != NULL) *reason = "The rasters have different SRIDs";
64  err = 1;
65  }
66  /* scales must match */
67  else if (FLT_NEQ(fabs(rast1->scaleX), fabs(rast2->scaleX))) {
68  if (reason != NULL) *reason = "The rasters have different scales on the X axis";
69  err = 1;
70  }
71  else if (FLT_NEQ(fabs(rast1->scaleY), fabs(rast2->scaleY))) {
72  if (reason != NULL) *reason = "The rasters have different scales on the Y axis";
73  err = 1;
74  }
75  /* skews must match */
76  else if (FLT_NEQ(rast1->skewX, rast2->skewX)) {
77  if (reason != NULL) *reason = "The rasters have different skews on the X axis";
78  err = 1;
79  }
80  else if (FLT_NEQ(rast1->skewY, rast2->skewY)) {
81  if (reason != NULL) *reason = "The rasters have different skews on the Y axis";
82  err = 1;
83  }
84 
85  if (err) {
86  *aligned = 0;
87  return ES_NONE;
88  }
89 
90  /* raster coordinates in context of second raster of first raster's upper-left corner */
92  rast2,
93  rast1->ipX, rast1->ipY,
94  &xr, &yr,
95  NULL
96  ) != ES_NONE) {
97  rterror("rt_raster_same_alignment: Could not get raster coordinates of second raster from first raster's spatial coordinates");
98  *aligned = 0;
99  return ES_ERROR;
100  }
101 
102  /* spatial coordinates of raster coordinates from above */
104  rast2,
105  xr, yr,
106  &xw, &yw,
107  NULL
108  ) != ES_NONE) {
109  rterror("rt_raster_same_alignment: Could not get spatial coordinates of second raster from raster coordinates");
110  *aligned = 0;
111  return ES_ERROR;
112  }
113 
114  RASTER_DEBUGF(4, "rast1(ipX, ipxY) = (%f, %f)", rast1->ipX, rast1->ipY);
115  RASTER_DEBUGF(4, "rast2(xr, yr) = (%f, %f)", xr, yr);
116  RASTER_DEBUGF(4, "rast2(xw, yw) = (%f, %f)", xw, yw);
117 
118  /* spatial coordinates are identical to that of first raster's upper-left corner */
119  if (FLT_EQ(xw, rast1->ipX) && FLT_EQ(yw, rast1->ipY)) {
120  if (reason != NULL) *reason = "The rasters are aligned";
121  *aligned = 1;
122  return ES_NONE;
123  }
124 
125  /* no alignment */
126  if (reason != NULL) *reason = "The rasters (pixel corner coordinates) are not aligned";
127 
128  *aligned = 0;
129  return ES_NONE;
130 }
131 
132 /******************************************************************************
133 * GEOS-based spatial relationship tests
134 ******************************************************************************/
135 
136 static
138  rt_raster rast1, int nband1,
139  rt_raster rast2, int nband2,
140  rt_geos_spatial_test testtype,
141  int *testresult
142 ) {
143  LWMPOLY *surface1 = NULL;
144  LWMPOLY *surface2 = NULL;
145  GEOSGeometry *geom1 = NULL;
146  GEOSGeometry *geom2 = NULL;
147  int rtn = 0;
148  int flag = 0;
149 
150  RASTER_DEBUG(3, "Starting");
151 
152  assert(NULL != rast1);
153  assert(NULL != rast2);
154  assert(NULL != testresult);
155 
156  if (nband1 < 0 && nband2 < 0) {
157  nband1 = -1;
158  nband2 = -1;
159  }
160  else {
161  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
162  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
163  }
164 
165  /* initialize to zero, false result of spatial relationship test */
166  *testresult = 0;
167 
168  /* same srid */
169  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
170  rterror("rt_raster_geos_spatial_relationship: The two rasters provided have different SRIDs");
171  return ES_ERROR;
172  }
173 
174  initGEOS(rtinfo, lwgeom_geos_error);
175 
176  /* get LWMPOLY of each band */
177  if (rt_raster_surface(rast1, nband1, &surface1) != ES_NONE) {
178  rterror("rt_raster_geos_spatial_relationship: Could not get surface of the specified band from the first raster");
179  return ES_ERROR;
180  }
181  if (rt_raster_surface(rast2, nband2, &surface2) != ES_NONE) {
182  rterror("rt_raster_geos_spatial_relationship: Could not get surface of the specified band from the second raster");
183  lwmpoly_free(surface1);
184  return ES_ERROR;
185  }
186 
187  /* either surface is NULL, spatial relationship test is false */
188  if (surface1 == NULL || surface2 == NULL) {
189  if (surface1 != NULL) lwmpoly_free(surface1);
190  if (surface2 != NULL) lwmpoly_free(surface2);
191  return ES_NONE;
192  }
193 
194  /* convert LWMPOLY to GEOSGeometry */
195  geom1 = LWGEOM2GEOS(lwmpoly_as_lwgeom(surface1), 0);
196  lwmpoly_free(surface1);
197  if (geom1 == NULL) {
198  rterror("rt_raster_geos_spatial_relationship: Could not convert surface of the specified band from the first raster to a GEOSGeometry");
199  lwmpoly_free(surface2);
200  return ES_ERROR;
201  }
202 
203  geom2 = LWGEOM2GEOS(lwmpoly_as_lwgeom(surface2), 0);
204  lwmpoly_free(surface2);
205  if (geom2 == NULL) {
206  rterror("rt_raster_geos_spatial_relationship: Could not convert surface of the specified band from the second raster to a GEOSGeometry");
207  return ES_ERROR;
208  }
209 
210  flag = 0;
211  switch (testtype) {
212  case GSR_OVERLAPS:
213  rtn = GEOSOverlaps(geom1, geom2);
214  break;
215  case GSR_TOUCHES:
216  rtn = GEOSTouches(geom1, geom2);
217  break;
218  case GSR_CONTAINS:
219  rtn = GEOSContains(geom1, geom2);
220  break;
222  rtn = GEOSRelatePattern(geom1, geom2, "T**FF*FF*");
223  break;
224  case GSR_COVERS:
225  rtn = GEOSRelatePattern(geom1, geom2, "******FF*");
226  break;
227  case GSR_COVEREDBY:
228  rtn = GEOSRelatePattern(geom1, geom2, "**F**F***");
229  break;
230  default:
231  rterror("rt_raster_geos_spatial_relationship: Unknown or unsupported GEOS spatial relationship test");
232  flag = -1;
233  break;
234  }
235  GEOSGeom_destroy(geom1);
236  GEOSGeom_destroy(geom2);
237 
238  /* something happened in the spatial relationship test */
239  if (rtn == 2) {
240  rterror("rt_raster_geos_spatial_relationship: Could not run the appropriate GEOS spatial relationship test");
241  flag = ES_ERROR;
242  }
243  /* spatial relationship test ran fine */
244  else if (flag >= 0) {
245  if (rtn != 0)
246  *testresult = 1;
247  flag = ES_NONE;
248  }
249  /* flag < 0 for when testtype is unknown */
250  else
251  flag = ES_ERROR;
252 
253  return flag;
254 }
255 
273  rt_raster rast1, int nband1,
274  rt_raster rast2, int nband2,
275  int *overlaps
276 ) {
277  RASTER_DEBUG(3, "Starting");
278 
280  rast1, nband1,
281  rast2, nband2,
282  GSR_OVERLAPS,
283  overlaps
284  );
285 }
286 
304  rt_raster rast1, int nband1,
305  rt_raster rast2, int nband2,
306  int *touches
307 ) {
308  RASTER_DEBUG(3, "Starting");
309 
311  rast1, nband1,
312  rast2, nband2,
313  GSR_TOUCHES,
314  touches
315  );
316 }
317 
335  rt_raster rast1, int nband1,
336  rt_raster rast2, int nband2,
337  int *contains
338 ) {
339  RASTER_DEBUG(3, "Starting");
340 
342  rast1, nband1,
343  rast2, nband2,
344  GSR_CONTAINS,
345  contains
346  );
347 }
348 
366  rt_raster rast1, int nband1,
367  rt_raster rast2, int nband2,
368  int *contains
369 ) {
370  RASTER_DEBUG(3, "Starting");
371 
373  rast1, nband1,
374  rast2, nband2,
376  contains
377  );
378 }
379 
397  rt_raster rast1, int nband1,
398  rt_raster rast2, int nband2,
399  int *covers
400 ) {
401  RASTER_DEBUG(3, "Starting");
402 
404  rast1, nband1,
405  rast2, nband2,
406  GSR_COVERS,
407  covers
408  );
409 }
410 
428  rt_raster rast1, int nband1,
429  rt_raster rast2, int nband2,
430  int *coveredby
431 ) {
432  RASTER_DEBUG(3, "Starting");
433 
435  rast1, nband1,
436  rast2, nband2,
438  coveredby
439  );
440 }
441 
461  rt_raster rast1, int nband1,
462  rt_raster rast2, int nband2,
463  double distance,
464  int *dwithin
465 ) {
466  LWMPOLY *surface = NULL;
467  LWGEOM *surface1 = NULL;
468  LWGEOM *surface2 = NULL;
469  double mindist = 0;
470 
471  RASTER_DEBUG(3, "Starting");
472 
473  assert(NULL != rast1);
474  assert(NULL != rast2);
475  assert(NULL != dwithin);
476 
477  if (nband1 < 0 && nband2 < 0) {
478  nband1 = -1;
479  nband2 = -1;
480  }
481  else {
482  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
483  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
484  }
485 
486  /* initialize to zero, false result */
487  *dwithin = 0;
488 
489  /* same srid */
490  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
491  rterror("rt_raster_distance_within: The two rasters provided have different SRIDs");
492  return ES_ERROR;
493  }
494 
495  /* distance cannot be less than zero */
496  if (distance < 0) {
497  rterror("rt_raster_distance_within: Distance cannot be less than zero");
498  return ES_ERROR;
499  }
500 
501  /* get LWMPOLY of each band */
502  if (rt_raster_surface(rast1, nband1, &surface) != ES_NONE) {
503  rterror("rt_raster_distance_within: Could not get surface of the specified band from the first raster");
504  return ES_ERROR;
505  }
506  surface1 = lwmpoly_as_lwgeom(surface);
507 
508  if (rt_raster_surface(rast2, nband2, &surface) != ES_NONE) {
509  rterror("rt_raster_distance_within: Could not get surface of the specified band from the second raster");
510  lwgeom_free(surface1);
511  return ES_ERROR;
512  }
513  surface2 = lwmpoly_as_lwgeom(surface);
514 
515  /* either surface is NULL, test is false */
516  if (surface1 == NULL || surface2 == NULL) {
517  if (surface1 != NULL) lwgeom_free(surface1);
518  if (surface2 != NULL) lwgeom_free(surface2);
519  return ES_NONE;
520  }
521 
522  /* get the min distance between the two surfaces */
523  mindist = lwgeom_mindistance2d_tolerance(surface1, surface2, distance);
524 
525  lwgeom_free(surface1);
526  lwgeom_free(surface2);
527 
528  /* if distance >= mindist, true */
529  if (FLT_EQ(mindist, distance) || distance > mindist)
530  *dwithin = 1;
531 
532  RASTER_DEBUGF(3, "(mindist, distance) = (%f, %f, %d)", mindist, distance, *dwithin);
533 
534  return ES_NONE;
535 }
536 
556  rt_raster rast1, int nband1,
557  rt_raster rast2, int nband2,
558  double distance,
559  int *dfwithin
560 ) {
561  LWMPOLY *surface = NULL;
562  LWGEOM *surface1 = NULL;
563  LWGEOM *surface2 = NULL;
564  double maxdist = 0;
565 
566  RASTER_DEBUG(3, "Starting");
567 
568  assert(NULL != rast1);
569  assert(NULL != rast2);
570  assert(NULL != dfwithin);
571 
572  if (nband1 < 0 && nband2 < 0) {
573  nband1 = -1;
574  nband2 = -1;
575  }
576  else {
577  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
578  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
579  }
580 
581  /* initialize to zero, false result */
582  *dfwithin = 0;
583 
584  /* same srid */
585  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
586  rterror("rt_raster_fully_within_distance: The two rasters provided have different SRIDs");
587  return ES_ERROR;
588  }
589 
590  /* distance cannot be less than zero */
591  if (distance < 0) {
592  rterror("rt_raster_fully_within_distance: Distance cannot be less than zero");
593  return ES_ERROR;
594  }
595 
596  /* get LWMPOLY of each band */
597  if (rt_raster_surface(rast1, nband1, &surface) != ES_NONE) {
598  rterror("rt_raster_fully_within_distance: Could not get surface of the specified band from the first raster");
599  return ES_ERROR;
600  }
601  surface1 = lwmpoly_as_lwgeom(surface);
602 
603  if (rt_raster_surface(rast2, nband2, &surface) != ES_NONE) {
604  rterror("rt_raster_fully_within_distance: Could not get surface of the specified band from the second raster");
605  lwgeom_free(surface1);
606  return ES_ERROR;
607  }
608  surface2 = lwmpoly_as_lwgeom(surface);
609 
610  /* either surface is NULL, test is false */
611  if (surface1 == NULL || surface2 == NULL) {
612  if (surface1 != NULL) lwgeom_free(surface1);
613  if (surface2 != NULL) lwgeom_free(surface2);
614  return ES_NONE;
615  }
616 
617  /* get the maximum distance between the two surfaces */
618  maxdist = lwgeom_maxdistance2d_tolerance(surface1, surface2, distance);
619 
620  lwgeom_free(surface1);
621  lwgeom_free(surface2);
622 
623  /* if distance >= maxdist, true */
624  if (FLT_EQ(maxdist, distance) || distance > maxdist)
625  *dfwithin = 1;
626 
627  RASTER_DEBUGF(3, "(maxdist, distance, dfwithin) = (%f, %f, %d)", maxdist, distance, *dfwithin);
628 
629  return ES_NONE;
630 }
631 
632 /******************************************************************************
633 * rt_raster_intersects()
634 ******************************************************************************/
635 
636 static
638  rt_raster rast1, rt_raster rast2,
639  rt_band band1, rt_band band2,
640  int hasnodata1, int hasnodata2,
641  double nodata1, double nodata2
642 ) {
643  int i;
644  int byHeight = 1;
645  uint32_t dimValue;
646 
647  uint32_t row;
648  uint32_t rowoffset;
649  uint32_t col;
650  uint32_t coloffset;
651 
652  enum line_points {X1, Y1, X2, Y2};
653  enum point {pX, pY};
654  double line1[4] = {0.};
655  double line2[4] = {0.};
656  double P[2] = {0.};
657  double Qw[2] = {0.};
658  double Qr[2] = {0.};
659  double gt1[6] = {0.};
660  double gt2[6] = {0.};
661  double igt1[6] = {0};
662  double igt2[6] = {0};
663  double d;
664  double val1;
665  int noval1;
666  int isnodata1;
667  double val2;
668  int noval2;
669  int isnodata2;
670  uint32_t adjacent[8] = {0};
671 
672  double xscale;
673  double yscale;
674 
675  uint16_t width1;
676  uint16_t height1;
677  uint16_t width2;
678  uint16_t height2;
679 
680  width1 = rt_raster_get_width(rast1);
681  height1 = rt_raster_get_height(rast1);
682  width2 = rt_raster_get_width(rast2);
683  height2 = rt_raster_get_height(rast2);
684 
685  /* sampling scale */
686  xscale = fmin(rt_raster_get_x_scale(rast1), rt_raster_get_x_scale(rast2)) / 10.;
687  yscale = fmin(rt_raster_get_y_scale(rast1), rt_raster_get_y_scale(rast2)) / 10.;
688 
689  /* see if skew made rast2's rows are parallel to rast1's cols */
691  rast1,
692  0, 0,
693  &(line1[X1]), &(line1[Y1]),
694  gt1
695  );
696 
698  rast1,
699  0, height1,
700  &(line1[X2]), &(line1[Y2]),
701  gt1
702  );
703 
705  rast2,
706  0, 0,
707  &(line2[X1]), &(line2[Y1]),
708  gt2
709  );
710 
712  rast2,
713  width2, 0,
714  &(line2[X2]), &(line2[Y2]),
715  gt2
716  );
717 
718  /* Parallel lines */
719  if (FLT_EQ(((line1[X2] - line1[X1]) * (line2[Y2] - line2[Y1])),
720  ((line2[X2] - line2[X1]) * (line1[Y2] - line1[Y1]))))
721  byHeight = 0;
722 
723  if (byHeight)
724  dimValue = height2;
725  else
726  dimValue = width2;
727  RASTER_DEBUGF(4, "byHeight: %d, dimValue: %d", byHeight, dimValue);
728 
729  /* 3 x 3 search */
730  for (coloffset = 0; coloffset < 3; coloffset++) {
731  for (rowoffset = 0; rowoffset < 3; rowoffset++) {
732  /* smaller raster */
733  for (col = coloffset; col <= width1; col += 3) {
734 
736  rast1,
737  col, 0,
738  &(line1[X1]), &(line1[Y1]),
739  gt1
740  );
741 
743  rast1,
744  col, height1,
745  &(line1[X2]), &(line1[Y2]),
746  gt1
747  );
748 
749  /* larger raster */
750  for (row = rowoffset; row <= dimValue; row += 3) {
751 
752  if (byHeight) {
754  rast2,
755  0, row,
756  &(line2[X1]), &(line2[Y1]),
757  gt2
758  );
759 
761  rast2,
762  width2, row,
763  &(line2[X2]), &(line2[Y2]),
764  gt2
765  );
766  }
767  else {
769  rast2,
770  row, 0,
771  &(line2[X1]), &(line2[Y1]),
772  gt2
773  );
774 
776  rast2,
777  row, height2,
778  &(line2[X2]), &(line2[Y2]),
779  gt2
780  );
781  }
782 
783  RASTER_DEBUGF(4, "(col, row) = (%d, %d)", col, row);
784  RASTER_DEBUGF(4, "line1(x1, y1, x2, y2) = (%f, %f, %f, %f)",
785  line1[X1], line1[Y1], line1[X2], line1[Y2]);
786  RASTER_DEBUGF(4, "line2(x1, y1, x2, y2) = (%f, %f, %f, %f)",
787  line2[X1], line2[Y1], line2[X2], line2[Y2]);
788 
789  /* intersection */
790  /* http://en.wikipedia.org/wiki/Line-line_intersection */
791  d = ((line1[X1] - line1[X2]) * (line2[Y1] - line2[Y2])) - ((line1[Y1] - line1[Y2]) * (line2[X1] - line2[X2]));
792  /* no intersection */
793  if (FLT_EQ(d, 0.)) {
794  continue;
795  }
796 
797  P[pX] = (((line1[X1] * line1[Y2]) - (line1[Y1] * line1[X2])) * (line2[X1] - line2[X2])) -
798  ((line1[X1] - line1[X2]) * ((line2[X1] * line2[Y2]) - (line2[Y1] * line2[X2])));
799  P[pX] = P[pX] / d;
800 
801  P[pY] = (((line1[X1] * line1[Y2]) - (line1[Y1] * line1[X2])) * (line2[Y1] - line2[Y2])) -
802  ((line1[Y1] - line1[Y2]) * ((line2[X1] * line2[Y2]) - (line2[Y1] * line2[X2])));
803  P[pY] = P[pY] / d;
804 
805  RASTER_DEBUGF(4, "P(x, y) = (%f, %f)", P[pX], P[pY]);
806 
807  /* intersection within bounds */
808  if ((
809  (FLT_EQ(P[pX], line1[X1]) || FLT_EQ(P[pX], line1[X2])) ||
810  (P[pX] > fmin(line1[X1], line1[X2]) && (P[pX] < fmax(line1[X1], line1[X2])))
811  ) && (
812  (FLT_EQ(P[pY], line1[Y1]) || FLT_EQ(P[pY], line1[Y2])) ||
813  (P[pY] > fmin(line1[Y1], line1[Y2]) && (P[pY] < fmax(line1[Y1], line1[Y2])))
814  ) && (
815  (FLT_EQ(P[pX], line2[X1]) || FLT_EQ(P[pX], line2[X2])) ||
816  (P[pX] > fmin(line2[X1], line2[X2]) && (P[pX] < fmax(line2[X1], line2[X2])))
817  ) && (
818  (FLT_EQ(P[pY], line2[Y1]) || FLT_EQ(P[pY], line2[Y2])) ||
819  (P[pY] > fmin(line2[Y1], line2[Y2]) && (P[pY] < fmax(line2[Y1], line2[Y2])))
820  )) {
821  RASTER_DEBUG(4, "within bounds");
822 
823  for (i = 0; i < 8; i++) adjacent[i] = 0;
824 
825  /* test points around intersection */
826  for (i = 0; i < 8; i++) {
827  switch (i) {
828  case 7:
829  Qw[pX] = P[pX] - xscale;
830  Qw[pY] = P[pY] + yscale;
831  break;
832  /* 270 degrees = 09:00 */
833  case 6:
834  Qw[pX] = P[pX] - xscale;
835  Qw[pY] = P[pY];
836  break;
837  case 5:
838  Qw[pX] = P[pX] - xscale;
839  Qw[pY] = P[pY] - yscale;
840  break;
841  /* 180 degrees = 06:00 */
842  case 4:
843  Qw[pX] = P[pX];
844  Qw[pY] = P[pY] - yscale;
845  break;
846  case 3:
847  Qw[pX] = P[pX] + xscale;
848  Qw[pY] = P[pY] - yscale;
849  break;
850  /* 90 degrees = 03:00 */
851  case 2:
852  Qw[pX] = P[pX] + xscale;
853  Qw[pY] = P[pY];
854  break;
855  /* 45 degrees */
856  case 1:
857  Qw[pX] = P[pX] + xscale;
858  Qw[pY] = P[pY] + yscale;
859  break;
860  /* 0 degrees = 00:00 */
861  case 0:
862  Qw[pX] = P[pX];
863  Qw[pY] = P[pY] + yscale;
864  break;
865  }
866 
867  /* unable to convert point to cell */
868  noval1 = 0;
870  rast1,
871  Qw[pX], Qw[pY],
872  &(Qr[pX]), &(Qr[pY]),
873  igt1
874  ) != ES_NONE) {
875  noval1 = 1;
876  }
877  /* cell is outside bounds of grid */
878  else if (
879  (Qr[pX] < 0 || Qr[pX] > width1 || FLT_EQ(Qr[pX], width1)) ||
880  (Qr[pY] < 0 || Qr[pY] > height1 || FLT_EQ(Qr[pY], height1))
881  ) {
882  noval1 = 1;
883  }
884  else if (hasnodata1 == FALSE)
885  val1 = 1;
886  /* unable to get value at cell */
887  else if (rt_band_get_pixel(band1, Qr[pX], Qr[pY], &val1, &isnodata1) != ES_NONE)
888  noval1 = 1;
889 
890  /* unable to convert point to cell */
891  noval2 = 0;
893  rast2,
894  Qw[pX], Qw[pY],
895  &(Qr[pX]), &(Qr[pY]),
896  igt2
897  ) != ES_NONE) {
898  noval2 = 1;
899  }
900  /* cell is outside bounds of grid */
901  else if (
902  (Qr[pX] < 0 || Qr[pX] > width2 || FLT_EQ(Qr[pX], width2)) ||
903  (Qr[pY] < 0 || Qr[pY] > height2 || FLT_EQ(Qr[pY], height2))
904  ) {
905  noval2 = 1;
906  }
907  else if (hasnodata2 == FALSE)
908  val2 = 1;
909  /* unable to get value at cell */
910  else if (rt_band_get_pixel(band2, Qr[pX], Qr[pY], &val2, &isnodata2) != ES_NONE)
911  noval2 = 1;
912 
913  if (!noval1) {
914  RASTER_DEBUGF(4, "val1 = %f", val1);
915  }
916  if (!noval2) {
917  RASTER_DEBUGF(4, "val2 = %f", val2);
918  }
919 
920  /* pixels touch */
921  if (!noval1 && (
922  (hasnodata1 == FALSE) || !isnodata1
923  )) {
924  adjacent[i]++;
925  }
926  if (!noval2 && (
927  (hasnodata2 == FALSE) || !isnodata2
928  )) {
929  adjacent[i] += 3;
930  }
931 
932  /* two pixel values not present */
933  if (noval1 || noval2) {
934  RASTER_DEBUGF(4, "noval1 = %d, noval2 = %d", noval1, noval2);
935  continue;
936  }
937 
938  /* pixels valid, so intersect */
939  if (
940  ((hasnodata1 == FALSE) || !isnodata1) &&
941  ((hasnodata2 == FALSE) || !isnodata2)
942  ) {
943  RASTER_DEBUG(3, "The two rasters do intersect");
944 
945  return 1;
946  }
947  }
948 
949  /* pixels touch */
950  for (i = 0; i < 4; i++) {
951  RASTER_DEBUGF(4, "adjacent[%d] = %d, adjacent[%d] = %d"
952  , i, adjacent[i], i + 4, adjacent[i + 4]);
953  if (adjacent[i] == 0) continue;
954 
955  if (adjacent[i] + adjacent[i + 4] == 4) {
956  RASTER_DEBUG(3, "The two rasters touch");
957 
958  return 1;
959  }
960  }
961  }
962  else {
963  RASTER_DEBUG(4, "outside of bounds");
964  }
965  }
966  }
967  }
968  }
969 
970  return 0;
971 }
972 
991  rt_raster rast1, int nband1,
992  rt_raster rast2, int nband2,
993  int *intersects
994 ) {
995  int i;
996  int j;
997  int within = 0;
998 
999  LWGEOM *hull[2] = {NULL};
1000  GEOSGeometry *ghull[2] = {NULL};
1001 
1002  uint16_t width1;
1003  uint16_t height1;
1004  uint16_t width2;
1005  uint16_t height2;
1006  double area1;
1007  double area2;
1008  double pixarea1;
1009  double pixarea2;
1010  rt_raster rastS = NULL;
1011  rt_raster rastL = NULL;
1012  uint16_t *widthS = NULL;
1013  uint16_t *heightS = NULL;
1014  uint16_t *widthL = NULL;
1015  uint16_t *heightL = NULL;
1016  int nbandS;
1017  int nbandL;
1018  rt_band bandS = NULL;
1019  rt_band bandL = NULL;
1020  int hasnodataS = FALSE;
1021  int hasnodataL = FALSE;
1022  double nodataS = 0;
1023  double nodataL = 0;
1024  int isnodataS = 0;
1025  int isnodataL = 0;
1026  double gtS[6] = {0};
1027  double igtL[6] = {0};
1028 
1029  uint32_t row;
1030  uint32_t rowoffset;
1031  uint32_t col;
1032  uint32_t coloffset;
1033 
1034  enum line_points {X1, Y1, X2, Y2};
1035  enum point {pX, pY};
1036  double lineS[4];
1037  double Qr[2];
1038  double valS;
1039  double valL;
1040 
1041  RASTER_DEBUG(3, "Starting");
1042 
1043  assert(NULL != rast1);
1044  assert(NULL != rast2);
1045  assert(NULL != intersects);
1046 
1047  if (nband1 < 0 && nband2 < 0) {
1048  nband1 = -1;
1049  nband2 = -1;
1050  }
1051  else {
1052  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
1053  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
1054  }
1055 
1056  /* same srid */
1057  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
1058  rterror("rt_raster_intersects: The two rasters provided have different SRIDs");
1059  *intersects = 0;
1060  return ES_ERROR;
1061  }
1062 
1063  /* raster extents need to intersect */
1064  do {
1065  int rtn;
1066 
1067  initGEOS(rtinfo, lwgeom_geos_error);
1068 
1069  rtn = 1;
1070  for (i = 0; i < 2; i++) {
1071  if ((rt_raster_get_convex_hull(i < 1 ? rast1 : rast2, &(hull[i])) != ES_NONE) || NULL == hull[i]) {
1072  for (j = 0; j < i; j++) {
1073  GEOSGeom_destroy(ghull[j]);
1074  lwgeom_free(hull[j]);
1075  }
1076  rtn = 0;
1077  break;
1078  }
1079  ghull[i] = (GEOSGeometry *) LWGEOM2GEOS(hull[i], 0);
1080  if (NULL == ghull[i]) {
1081  for (j = 0; j < i; j++) {
1082  GEOSGeom_destroy(ghull[j]);
1083  lwgeom_free(hull[j]);
1084  }
1085  lwgeom_free(hull[i]);
1086  rtn = 0;
1087  break;
1088  }
1089  }
1090  if (!rtn) break;
1091 
1092  /* test to see if raster within the other */
1093  within = 0;
1094  if (GEOSWithin(ghull[0], ghull[1]) == 1)
1095  within = -1;
1096  else if (GEOSWithin(ghull[1], ghull[0]) == 1)
1097  within = 1;
1098 
1099  if (within != 0)
1100  rtn = 1;
1101  else
1102  rtn = GEOSIntersects(ghull[0], ghull[1]);
1103 
1104  for (i = 0; i < 2; i++) {
1105  GEOSGeom_destroy(ghull[i]);
1106  lwgeom_free(hull[i]);
1107  }
1108 
1109  if (rtn != 2) {
1110  RASTER_DEBUGF(4, "convex hulls of rasters do %sintersect", rtn != 1 ? "NOT " : "");
1111  if (rtn != 1) {
1112  *intersects = 0;
1113  return ES_NONE;
1114  }
1115  /* band isn't specified */
1116  else if (nband1 < 0) {
1117  *intersects = 1;
1118  return ES_NONE;
1119  }
1120  }
1121  else {
1122  RASTER_DEBUG(4, "GEOSIntersects() returned a 2!!!!");
1123  }
1124  }
1125  while (0);
1126 
1127  /* smaller raster by area or width */
1128  width1 = rt_raster_get_width(rast1);
1129  height1 = rt_raster_get_height(rast1);
1130  width2 = rt_raster_get_width(rast2);
1131  height2 = rt_raster_get_height(rast2);
1132  pixarea1 = fabs(rt_raster_get_x_scale(rast1) * rt_raster_get_y_scale(rast1));
1133  pixarea2 = fabs(rt_raster_get_x_scale(rast2) * rt_raster_get_y_scale(rast2));
1134  area1 = fabs(width1 * height1 * pixarea1);
1135  area2 = fabs(width2 * height2 * pixarea2);
1136  RASTER_DEBUGF(4, "pixarea1, pixarea2, area1, area2 = %f, %f, %f, %f",
1137  pixarea1, pixarea2, area1, area2);
1138  if (
1139  (within <= 0) ||
1140  (area1 < area2) ||
1141  FLT_EQ(area1, area2) ||
1142  (area1 < pixarea2) || /* area of rast1 smaller than pixel area of rast2 */
1143  FLT_EQ(area1, pixarea2)
1144  ) {
1145  rastS = rast1;
1146  nbandS = nband1;
1147  widthS = &width1;
1148  heightS = &height1;
1149 
1150  rastL = rast2;
1151  nbandL = nband2;
1152  widthL = &width2;
1153  heightL = &height2;
1154  }
1155  else {
1156  rastS = rast2;
1157  nbandS = nband2;
1158  widthS = &width2;
1159  heightS = &height2;
1160 
1161  rastL = rast1;
1162  nbandL = nband1;
1163  widthL = &width1;
1164  heightL = &height1;
1165  }
1166 
1167  /* no band to use, set band to zero */
1168  if (nband1 < 0) {
1169  nbandS = 0;
1170  nbandL = 0;
1171  }
1172 
1173  RASTER_DEBUGF(4, "rast1 @ %p", rast1);
1174  RASTER_DEBUGF(4, "rast2 @ %p", rast2);
1175  RASTER_DEBUGF(4, "rastS @ %p", rastS);
1176  RASTER_DEBUGF(4, "rastL @ %p", rastL);
1177 
1178  /* load band of smaller raster */
1179  bandS = rt_raster_get_band(rastS, nbandS);
1180  if (NULL == bandS) {
1181  rterror("rt_raster_intersects: Could not get band %d of the first raster", nbandS);
1182  *intersects = 0;
1183  return ES_ERROR;
1184  }
1185 
1186  hasnodataS = rt_band_get_hasnodata_flag(bandS);
1187  if (hasnodataS != FALSE)
1188  rt_band_get_nodata(bandS, &nodataS);
1189 
1190  /* load band of larger raster */
1191  bandL = rt_raster_get_band(rastL, nbandL);
1192  if (NULL == bandL) {
1193  rterror("rt_raster_intersects: Could not get band %d of the first raster", nbandL);
1194  *intersects = 0;
1195  return ES_ERROR;
1196  }
1197 
1198  hasnodataL = rt_band_get_hasnodata_flag(bandL);
1199  if (hasnodataL != FALSE)
1200  rt_band_get_nodata(bandL, &nodataL);
1201 
1202  /* no band to use, ignore nodata */
1203  if (nband1 < 0) {
1204  hasnodataS = FALSE;
1205  hasnodataL = FALSE;
1206  }
1207 
1208  /* hasnodata(S|L) = TRUE and one of the two bands is isnodata */
1209  if (
1210  (hasnodataS && rt_band_get_isnodata_flag(bandS)) ||
1211  (hasnodataL && rt_band_get_isnodata_flag(bandL))
1212  ) {
1213  RASTER_DEBUG(3, "One of the two raster bands is NODATA. The two rasters do not intersect");
1214  *intersects = 0;
1215  return ES_NONE;
1216  }
1217 
1218  /* special case where a raster can fit inside another raster's pixel */
1219  if (within != 0 && ((pixarea1 > area2) || (pixarea2 > area1))) {
1220  RASTER_DEBUG(4, "Using special case of raster fitting into another raster's pixel");
1221  /* 3 x 3 search */
1222  for (coloffset = 0; coloffset < 3; coloffset++) {
1223  for (rowoffset = 0; rowoffset < 3; rowoffset++) {
1224  for (col = coloffset; col < *widthS; col += 3) {
1225  for (row = rowoffset; row < *heightS; row += 3) {
1226  if (hasnodataS == FALSE)
1227  valS = 1;
1228  else if (rt_band_get_pixel(bandS, col, row, &valS, &isnodataS) != ES_NONE)
1229  continue;
1230 
1231  if ((hasnodataS == FALSE) || !isnodataS) {
1233  rastS,
1234  col, row,
1235  &(lineS[X1]), &(lineS[Y1]),
1236  gtS
1237  );
1238 
1240  rastL,
1241  lineS[X1], lineS[Y1],
1242  &(Qr[pX]), &(Qr[pY]),
1243  igtL
1244  ) != ES_NONE) {
1245  continue;
1246  }
1247 
1248  if (
1249  (Qr[pX] < 0 || Qr[pX] > *widthL || FLT_EQ(Qr[pX], *widthL)) ||
1250  (Qr[pY] < 0 || Qr[pY] > *heightL || FLT_EQ(Qr[pY], *heightL))
1251  ) {
1252  continue;
1253  }
1254 
1255  if (hasnodataS == FALSE)
1256  valL = 1;
1257  else if (rt_band_get_pixel(bandL, Qr[pX], Qr[pY], &valL, &isnodataL) != ES_NONE)
1258  continue;
1259 
1260  if ((hasnodataL == FALSE) || !isnodataL) {
1261  RASTER_DEBUG(3, "The two rasters do intersect");
1262  *intersects = 1;
1263  return ES_NONE;
1264  }
1265  }
1266  }
1267  }
1268  }
1269  }
1270  RASTER_DEBUG(4, "Smaller raster not in the other raster's pixel. Continuing");
1271  }
1272 
1273  RASTER_DEBUG(4, "Testing smaller raster vs larger raster");
1274  *intersects = rt_raster_intersects_algorithm(
1275  rastS, rastL,
1276  bandS, bandL,
1277  hasnodataS, hasnodataL,
1278  nodataS, nodataL
1279  );
1280 
1281  if (*intersects) return ES_NONE;
1282 
1283  RASTER_DEBUG(4, "Testing larger raster vs smaller raster");
1284  *intersects = rt_raster_intersects_algorithm(
1285  rastL, rastS,
1286  bandL, bandS,
1287  hasnodataL, hasnodataS,
1288  nodataL, nodataS
1289  );
1290 
1291  if (*intersects) return ES_NONE;
1292 
1293  RASTER_DEBUG(3, "The two rasters do not intersect");
1294 
1295  *intersects = 0;
1296  return ES_NONE;
1297 }
1298 
Datum coveredby(PG_FUNCTION_ARGS)
double skewY
Definition: librtcore.h:2249
int rt_raster_get_num_bands(rt_raster raster)
Definition: rt_raster.c:372
Datum covers(PG_FUNCTION_ARGS)
rt_errorstate rt_raster_fully_within_distance(rt_raster rast1, int nband1, rt_raster rast2, int nband2, double distance, int *dfwithin)
Return ES_ERROR if error occurred in function.
rt_errorstate rt_raster_same_alignment(rt_raster rast1, rt_raster rast2, int *aligned, char **reason)
double lwgeom_maxdistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
Function handling max distance calculations and dfyllywithin calculations.
Definition: measures.c:181
double lwgeom_mindistance2d_tolerance(const LWGEOM *lw1, const LWGEOM *lw2, double tolerance)
Function handling min distance calculations and dwithin calculations.
Definition: measures.c:213
rt_errorstate rt_raster_covers(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *covers)
Return ES_ERROR if error occurred in function.
void lwgeom_free(LWGEOM *geom)
Definition: lwgeom.c:1099
Datum contains(PG_FUNCTION_ARGS)
static int rt_raster_intersects_algorithm(rt_raster rast1, rt_raster rast2, rt_band band1, rt_band band2, int hasnodata1, int hasnodata2, double nodata1, double nodata2)
#define FLT_EQ(x, y)
Definition: librtcore.h:2185
void rterror(const char *fmt,...)
Wrappers used for reporting errors and info.
Definition: rt_context.c:199
rt_errorstate
Enum definitions.
Definition: librtcore.h:179
double ipY
Definition: librtcore.h:2247
rt_errorstate rt_raster_overlaps(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *overlaps)
Return ES_ERROR if error occurred in function.
Datum overlaps(PG_FUNCTION_ARGS)
rt_errorstate rt_raster_coveredby(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *coveredby)
Return ES_ERROR if error occurred in function.
static rt_errorstate rt_raster_geos_spatial_relationship(rt_raster rast1, int nband1, rt_raster rast2, int nband2, rt_geos_spatial_test testtype, int *testresult)
Datum touches(PG_FUNCTION_ARGS)
rt_geos_spatial_test
GEOS spatial relationship tests available.
Definition: librtcore.h:217
rt_errorstate rt_band_get_nodata(rt_band band, double *nodata)
Get NODATA value.
Definition: rt_band.c:1597
rt_errorstate rt_raster_get_convex_hull(rt_raster raster, LWGEOM **hull)
Get raster&#39;s convex hull.
Definition: rt_geometry.c:803
unsigned int uint32_t
Definition: uthash.h:78
rt_errorstate rt_raster_cell_to_geopoint(rt_raster raster, double xr, double yr, double *xw, double *yw, double *gt)
Convert an xr, yr raster point to an xw, yw point on map.
Definition: rt_raster.c:755
rt_errorstate rt_raster_surface(rt_raster raster, int nband, LWMPOLY **surface)
Get a raster as a surface (multipolygon).
Definition: rt_geometry.c:355
void lwmpoly_free(LWMPOLY *mpoly)
Definition: lwmpoly.c:53
rt_errorstate rt_raster_within_distance(rt_raster rast1, int nband1, rt_raster rast2, int nband2, double distance, int *dwithin)
Return ES_ERROR if error occurred in function.
rt_errorstate rt_raster_contains_properly(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *contains)
Return ES_ERROR if error occurred in function.
void lwgeom_geos_error(const char *fmt,...)
Datum intersects(PG_FUNCTION_ARGS)
rt_errorstate rt_band_get_pixel(rt_band band, int x, int y, double *value, int *nodata)
Get pixel value.
Definition: rt_band.c:1088
rt_band rt_raster_get_band(rt_raster raster, int bandNum)
Return Nth band, or NULL if unavailable.
Definition: rt_raster.c:381
#define FLT_NEQ(x, y)
Definition: librtcore.h:2184
void rtinfo(const char *fmt,...)
Definition: rt_context.c:211
int rt_band_get_hasnodata_flag(rt_band band)
Get hasnodata flag value.
Definition: rt_band.c:541
double rt_raster_get_x_scale(rt_raster raster)
Get scale X in projection units.
Definition: rt_raster.c:150
double skewX
Definition: librtcore.h:2248
#define RASTER_DEBUGF(level, msg,...)
Definition: librtcore.h:299
int32_t rt_raster_get_srid(rt_raster raster)
Get raster&#39;s SRID.
Definition: rt_raster.c:356
rt_errorstate rt_raster_geopoint_to_cell(rt_raster raster, double xw, double yw, double *xr, double *yr, double *igt)
Convert an xw, yw map point to a xr, yr raster point.
Definition: rt_raster.c:806
GEOSGeometry * LWGEOM2GEOS(const LWGEOM *lwgeom, int autofix)
double rt_raster_get_y_scale(rt_raster raster)
Get scale Y in projection units.
Definition: rt_raster.c:159
Datum distance(PG_FUNCTION_ARGS)
rt_errorstate rt_raster_touches(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *touches)
Return ES_ERROR if error occurred in function.
uint16_t rt_raster_get_width(rt_raster raster)
Definition: rt_raster.c:121
rt_errorstate rt_raster_intersects(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *intersects)
Return zero if error occurred in function.
#define FALSE
Definition: dbfopen.c:168
#define RASTER_DEBUG(level, msg)
Definition: librtcore.h:295
This library is the generic raster handling section of PostGIS.
rt_errorstate rt_raster_contains(rt_raster rast1, int nband1, rt_raster rast2, int nband2, int *contains)
Return ES_ERROR if error occurred in function.
double ipX
Definition: librtcore.h:2246
int rt_band_get_isnodata_flag(rt_band band)
Get isnodata flag value.
Definition: rt_band.c:581
uint16_t rt_raster_get_height(rt_raster raster)
Definition: rt_raster.c:129
double scaleX
Definition: librtcore.h:2244
LWGEOM * lwmpoly_as_lwgeom(const LWMPOLY *obj)
Definition: lwgeom.c:253
double scaleY
Definition: librtcore.h:2245