PostGIS  2.3.7dev-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 vertically */
719  if (FLT_EQ(line1[X2] - line1[X1], 0.) && FLT_EQ(line2[X2] - line2[X1], 0.))
720  byHeight = 0;
721  /* parallel */
722  else if (FLT_EQ(((line1[Y2] - line1[Y1]) / (line1[X2] - line1[X1])), ((line2[Y2] - line2[Y1]) / (line2[X2] - line2[X1]))))
723  byHeight = 0;
724 
725  if (byHeight)
726  dimValue = height2;
727  else
728  dimValue = width2;
729  RASTER_DEBUGF(4, "byHeight: %d, dimValue: %d", byHeight, dimValue);
730 
731  /* 3 x 3 search */
732  for (coloffset = 0; coloffset < 3; coloffset++) {
733  for (rowoffset = 0; rowoffset < 3; rowoffset++) {
734  /* smaller raster */
735  for (col = coloffset; col <= width1; col += 3) {
736 
738  rast1,
739  col, 0,
740  &(line1[X1]), &(line1[Y1]),
741  gt1
742  );
743 
745  rast1,
746  col, height1,
747  &(line1[X2]), &(line1[Y2]),
748  gt1
749  );
750 
751  /* larger raster */
752  for (row = rowoffset; row <= dimValue; row += 3) {
753 
754  if (byHeight) {
756  rast2,
757  0, row,
758  &(line2[X1]), &(line2[Y1]),
759  gt2
760  );
761 
763  rast2,
764  width2, row,
765  &(line2[X2]), &(line2[Y2]),
766  gt2
767  );
768  }
769  else {
771  rast2,
772  row, 0,
773  &(line2[X1]), &(line2[Y1]),
774  gt2
775  );
776 
778  rast2,
779  row, height2,
780  &(line2[X2]), &(line2[Y2]),
781  gt2
782  );
783  }
784 
785  RASTER_DEBUGF(4, "(col, row) = (%d, %d)", col, row);
786  RASTER_DEBUGF(4, "line1(x1, y1, x2, y2) = (%f, %f, %f, %f)",
787  line1[X1], line1[Y1], line1[X2], line1[Y2]);
788  RASTER_DEBUGF(4, "line2(x1, y1, x2, y2) = (%f, %f, %f, %f)",
789  line2[X1], line2[Y1], line2[X2], line2[Y2]);
790 
791  /* intersection */
792  /* http://en.wikipedia.org/wiki/Line-line_intersection */
793  d = ((line1[X1] - line1[X2]) * (line2[Y1] - line2[Y2])) - ((line1[Y1] - line1[Y2]) * (line2[X1] - line2[X2]));
794  /* no intersection */
795  if (FLT_EQ(d, 0.)) {
796  continue;
797  }
798 
799  P[pX] = (((line1[X1] * line1[Y2]) - (line1[Y1] * line1[X2])) * (line2[X1] - line2[X2])) -
800  ((line1[X1] - line1[X2]) * ((line2[X1] * line2[Y2]) - (line2[Y1] * line2[X2])));
801  P[pX] = P[pX] / d;
802 
803  P[pY] = (((line1[X1] * line1[Y2]) - (line1[Y1] * line1[X2])) * (line2[Y1] - line2[Y2])) -
804  ((line1[Y1] - line1[Y2]) * ((line2[X1] * line2[Y2]) - (line2[Y1] * line2[X2])));
805  P[pY] = P[pY] / d;
806 
807  RASTER_DEBUGF(4, "P(x, y) = (%f, %f)", P[pX], P[pY]);
808 
809  /* intersection within bounds */
810  if ((
811  (FLT_EQ(P[pX], line1[X1]) || FLT_EQ(P[pX], line1[X2])) ||
812  (P[pX] > fmin(line1[X1], line1[X2]) && (P[pX] < fmax(line1[X1], line1[X2])))
813  ) && (
814  (FLT_EQ(P[pY], line1[Y1]) || FLT_EQ(P[pY], line1[Y2])) ||
815  (P[pY] > fmin(line1[Y1], line1[Y2]) && (P[pY] < fmax(line1[Y1], line1[Y2])))
816  ) && (
817  (FLT_EQ(P[pX], line2[X1]) || FLT_EQ(P[pX], line2[X2])) ||
818  (P[pX] > fmin(line2[X1], line2[X2]) && (P[pX] < fmax(line2[X1], line2[X2])))
819  ) && (
820  (FLT_EQ(P[pY], line2[Y1]) || FLT_EQ(P[pY], line2[Y2])) ||
821  (P[pY] > fmin(line2[Y1], line2[Y2]) && (P[pY] < fmax(line2[Y1], line2[Y2])))
822  )) {
823  RASTER_DEBUG(4, "within bounds");
824 
825  for (i = 0; i < 8; i++) adjacent[i] = 0;
826 
827  /* test points around intersection */
828  for (i = 0; i < 8; i++) {
829  switch (i) {
830  case 7:
831  Qw[pX] = P[pX] - xscale;
832  Qw[pY] = P[pY] + yscale;
833  break;
834  /* 270 degrees = 09:00 */
835  case 6:
836  Qw[pX] = P[pX] - xscale;
837  Qw[pY] = P[pY];
838  break;
839  case 5:
840  Qw[pX] = P[pX] - xscale;
841  Qw[pY] = P[pY] - yscale;
842  break;
843  /* 180 degrees = 06:00 */
844  case 4:
845  Qw[pX] = P[pX];
846  Qw[pY] = P[pY] - yscale;
847  break;
848  case 3:
849  Qw[pX] = P[pX] + xscale;
850  Qw[pY] = P[pY] - yscale;
851  break;
852  /* 90 degrees = 03:00 */
853  case 2:
854  Qw[pX] = P[pX] + xscale;
855  Qw[pY] = P[pY];
856  break;
857  /* 45 degrees */
858  case 1:
859  Qw[pX] = P[pX] + xscale;
860  Qw[pY] = P[pY] + yscale;
861  break;
862  /* 0 degrees = 00:00 */
863  case 0:
864  Qw[pX] = P[pX];
865  Qw[pY] = P[pY] + yscale;
866  break;
867  }
868 
869  /* unable to convert point to cell */
870  noval1 = 0;
872  rast1,
873  Qw[pX], Qw[pY],
874  &(Qr[pX]), &(Qr[pY]),
875  igt1
876  ) != ES_NONE) {
877  noval1 = 1;
878  }
879  /* cell is outside bounds of grid */
880  else if (
881  (Qr[pX] < 0 || Qr[pX] > width1 || FLT_EQ(Qr[pX], width1)) ||
882  (Qr[pY] < 0 || Qr[pY] > height1 || FLT_EQ(Qr[pY], height1))
883  ) {
884  noval1 = 1;
885  }
886  else if (hasnodata1 == FALSE)
887  val1 = 1;
888  /* unable to get value at cell */
889  else if (rt_band_get_pixel(band1, Qr[pX], Qr[pY], &val1, &isnodata1) != ES_NONE)
890  noval1 = 1;
891 
892  /* unable to convert point to cell */
893  noval2 = 0;
895  rast2,
896  Qw[pX], Qw[pY],
897  &(Qr[pX]), &(Qr[pY]),
898  igt2
899  ) != ES_NONE) {
900  noval2 = 1;
901  }
902  /* cell is outside bounds of grid */
903  else if (
904  (Qr[pX] < 0 || Qr[pX] > width2 || FLT_EQ(Qr[pX], width2)) ||
905  (Qr[pY] < 0 || Qr[pY] > height2 || FLT_EQ(Qr[pY], height2))
906  ) {
907  noval2 = 1;
908  }
909  else if (hasnodata2 == FALSE)
910  val2 = 1;
911  /* unable to get value at cell */
912  else if (rt_band_get_pixel(band2, Qr[pX], Qr[pY], &val2, &isnodata2) != ES_NONE)
913  noval2 = 1;
914 
915  if (!noval1) {
916  RASTER_DEBUGF(4, "val1 = %f", val1);
917  }
918  if (!noval2) {
919  RASTER_DEBUGF(4, "val2 = %f", val2);
920  }
921 
922  /* pixels touch */
923  if (!noval1 && (
924  (hasnodata1 == FALSE) || !isnodata1
925  )) {
926  adjacent[i]++;
927  }
928  if (!noval2 && (
929  (hasnodata2 == FALSE) || !isnodata2
930  )) {
931  adjacent[i] += 3;
932  }
933 
934  /* two pixel values not present */
935  if (noval1 || noval2) {
936  RASTER_DEBUGF(4, "noval1 = %d, noval2 = %d", noval1, noval2);
937  continue;
938  }
939 
940  /* pixels valid, so intersect */
941  if (
942  ((hasnodata1 == FALSE) || !isnodata1) &&
943  ((hasnodata2 == FALSE) || !isnodata2)
944  ) {
945  RASTER_DEBUG(3, "The two rasters do intersect");
946 
947  return 1;
948  }
949  }
950 
951  /* pixels touch */
952  for (i = 0; i < 4; i++) {
953  RASTER_DEBUGF(4, "adjacent[%d] = %d, adjacent[%d] = %d"
954  , i, adjacent[i], i + 4, adjacent[i + 4]);
955  if (adjacent[i] == 0) continue;
956 
957  if (adjacent[i] + adjacent[i + 4] == 4) {
958  RASTER_DEBUG(3, "The two rasters touch");
959 
960  return 1;
961  }
962  }
963  }
964  else {
965  RASTER_DEBUG(4, "outside of bounds");
966  }
967  }
968  }
969  }
970  }
971 
972  return 0;
973 }
974 
993  rt_raster rast1, int nband1,
994  rt_raster rast2, int nband2,
995  int *intersects
996 ) {
997  int i;
998  int j;
999  int within = 0;
1000 
1001  LWGEOM *hull[2] = {NULL};
1002  GEOSGeometry *ghull[2] = {NULL};
1003 
1004  uint16_t width1;
1005  uint16_t height1;
1006  uint16_t width2;
1007  uint16_t height2;
1008  double area1;
1009  double area2;
1010  double pixarea1;
1011  double pixarea2;
1012  rt_raster rastS = NULL;
1013  rt_raster rastL = NULL;
1014  uint16_t *widthS = NULL;
1015  uint16_t *heightS = NULL;
1016  uint16_t *widthL = NULL;
1017  uint16_t *heightL = NULL;
1018  int nbandS;
1019  int nbandL;
1020  rt_band bandS = NULL;
1021  rt_band bandL = NULL;
1022  int hasnodataS = FALSE;
1023  int hasnodataL = FALSE;
1024  double nodataS = 0;
1025  double nodataL = 0;
1026  int isnodataS = 0;
1027  int isnodataL = 0;
1028  double gtS[6] = {0};
1029  double igtL[6] = {0};
1030 
1031  uint32_t row;
1032  uint32_t rowoffset;
1033  uint32_t col;
1034  uint32_t coloffset;
1035 
1036  enum line_points {X1, Y1, X2, Y2};
1037  enum point {pX, pY};
1038  double lineS[4];
1039  double Qr[2];
1040  double valS;
1041  double valL;
1042 
1043  RASTER_DEBUG(3, "Starting");
1044 
1045  assert(NULL != rast1);
1046  assert(NULL != rast2);
1047  assert(NULL != intersects);
1048 
1049  if (nband1 < 0 && nband2 < 0) {
1050  nband1 = -1;
1051  nband2 = -1;
1052  }
1053  else {
1054  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
1055  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
1056  }
1057 
1058  /* same srid */
1059  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
1060  rterror("rt_raster_intersects: The two rasters provided have different SRIDs");
1061  *intersects = 0;
1062  return ES_ERROR;
1063  }
1064 
1065  /* raster extents need to intersect */
1066  do {
1067  int rtn;
1068 
1069  initGEOS(rtinfo, lwgeom_geos_error);
1070 
1071  rtn = 1;
1072  for (i = 0; i < 2; i++) {
1073  if ((rt_raster_get_convex_hull(i < 1 ? rast1 : rast2, &(hull[i])) != ES_NONE) || NULL == hull[i]) {
1074  for (j = 0; j < i; j++) {
1075  GEOSGeom_destroy(ghull[j]);
1076  lwgeom_free(hull[j]);
1077  }
1078  rtn = 0;
1079  break;
1080  }
1081  ghull[i] = (GEOSGeometry *) LWGEOM2GEOS(hull[i], 0);
1082  if (NULL == ghull[i]) {
1083  for (j = 0; j < i; j++) {
1084  GEOSGeom_destroy(ghull[j]);
1085  lwgeom_free(hull[j]);
1086  }
1087  lwgeom_free(hull[i]);
1088  rtn = 0;
1089  break;
1090  }
1091  }
1092  if (!rtn) break;
1093 
1094  /* test to see if raster within the other */
1095  within = 0;
1096  if (GEOSWithin(ghull[0], ghull[1]) == 1)
1097  within = -1;
1098  else if (GEOSWithin(ghull[1], ghull[0]) == 1)
1099  within = 1;
1100 
1101  if (within != 0)
1102  rtn = 1;
1103  else
1104  rtn = GEOSIntersects(ghull[0], ghull[1]);
1105 
1106  for (i = 0; i < 2; i++) {
1107  GEOSGeom_destroy(ghull[i]);
1108  lwgeom_free(hull[i]);
1109  }
1110 
1111  if (rtn != 2) {
1112  RASTER_DEBUGF(4, "convex hulls of rasters do %sintersect", rtn != 1 ? "NOT " : "");
1113  if (rtn != 1) {
1114  *intersects = 0;
1115  return ES_NONE;
1116  }
1117  /* band isn't specified */
1118  else if (nband1 < 0) {
1119  *intersects = 1;
1120  return ES_NONE;
1121  }
1122  }
1123  else {
1124  RASTER_DEBUG(4, "GEOSIntersects() returned a 2!!!!");
1125  }
1126  }
1127  while (0);
1128 
1129  /* smaller raster by area or width */
1130  width1 = rt_raster_get_width(rast1);
1131  height1 = rt_raster_get_height(rast1);
1132  width2 = rt_raster_get_width(rast2);
1133  height2 = rt_raster_get_height(rast2);
1134  pixarea1 = fabs(rt_raster_get_x_scale(rast1) * rt_raster_get_y_scale(rast1));
1135  pixarea2 = fabs(rt_raster_get_x_scale(rast2) * rt_raster_get_y_scale(rast2));
1136  area1 = fabs(width1 * height1 * pixarea1);
1137  area2 = fabs(width2 * height2 * pixarea2);
1138  RASTER_DEBUGF(4, "pixarea1, pixarea2, area1, area2 = %f, %f, %f, %f",
1139  pixarea1, pixarea2, area1, area2);
1140  if (
1141  (within <= 0) ||
1142  (area1 < area2) ||
1143  FLT_EQ(area1, area2) ||
1144  (area1 < pixarea2) || /* area of rast1 smaller than pixel area of rast2 */
1145  FLT_EQ(area1, pixarea2)
1146  ) {
1147  rastS = rast1;
1148  nbandS = nband1;
1149  widthS = &width1;
1150  heightS = &height1;
1151 
1152  rastL = rast2;
1153  nbandL = nband2;
1154  widthL = &width2;
1155  heightL = &height2;
1156  }
1157  else {
1158  rastS = rast2;
1159  nbandS = nband2;
1160  widthS = &width2;
1161  heightS = &height2;
1162 
1163  rastL = rast1;
1164  nbandL = nband1;
1165  widthL = &width1;
1166  heightL = &height1;
1167  }
1168 
1169  /* no band to use, set band to zero */
1170  if (nband1 < 0) {
1171  nbandS = 0;
1172  nbandL = 0;
1173  }
1174 
1175  RASTER_DEBUGF(4, "rast1 @ %p", rast1);
1176  RASTER_DEBUGF(4, "rast2 @ %p", rast2);
1177  RASTER_DEBUGF(4, "rastS @ %p", rastS);
1178  RASTER_DEBUGF(4, "rastL @ %p", rastL);
1179 
1180  /* load band of smaller raster */
1181  bandS = rt_raster_get_band(rastS, nbandS);
1182  if (NULL == bandS) {
1183  rterror("rt_raster_intersects: Could not get band %d of the first raster", nbandS);
1184  *intersects = 0;
1185  return ES_ERROR;
1186  }
1187 
1188  hasnodataS = rt_band_get_hasnodata_flag(bandS);
1189  if (hasnodataS != FALSE)
1190  rt_band_get_nodata(bandS, &nodataS);
1191 
1192  /* load band of larger raster */
1193  bandL = rt_raster_get_band(rastL, nbandL);
1194  if (NULL == bandL) {
1195  rterror("rt_raster_intersects: Could not get band %d of the first raster", nbandL);
1196  *intersects = 0;
1197  return ES_ERROR;
1198  }
1199 
1200  hasnodataL = rt_band_get_hasnodata_flag(bandL);
1201  if (hasnodataL != FALSE)
1202  rt_band_get_nodata(bandL, &nodataL);
1203 
1204  /* no band to use, ignore nodata */
1205  if (nband1 < 0) {
1206  hasnodataS = FALSE;
1207  hasnodataL = FALSE;
1208  }
1209 
1210  /* hasnodata(S|L) = TRUE and one of the two bands is isnodata */
1211  if (
1212  (hasnodataS && rt_band_get_isnodata_flag(bandS)) ||
1213  (hasnodataL && rt_band_get_isnodata_flag(bandL))
1214  ) {
1215  RASTER_DEBUG(3, "One of the two raster bands is NODATA. The two rasters do not intersect");
1216  *intersects = 0;
1217  return ES_NONE;
1218  }
1219 
1220  /* special case where a raster can fit inside another raster's pixel */
1221  if (within != 0 && ((pixarea1 > area2) || (pixarea2 > area1))) {
1222  RASTER_DEBUG(4, "Using special case of raster fitting into another raster's pixel");
1223  /* 3 x 3 search */
1224  for (coloffset = 0; coloffset < 3; coloffset++) {
1225  for (rowoffset = 0; rowoffset < 3; rowoffset++) {
1226  for (col = coloffset; col < *widthS; col += 3) {
1227  for (row = rowoffset; row < *heightS; row += 3) {
1228  if (hasnodataS == FALSE)
1229  valS = 1;
1230  else if (rt_band_get_pixel(bandS, col, row, &valS, &isnodataS) != ES_NONE)
1231  continue;
1232 
1233  if ((hasnodataS == FALSE) || !isnodataS) {
1235  rastS,
1236  col, row,
1237  &(lineS[X1]), &(lineS[Y1]),
1238  gtS
1239  );
1240 
1242  rastL,
1243  lineS[X1], lineS[Y1],
1244  &(Qr[pX]), &(Qr[pY]),
1245  igtL
1246  ) != ES_NONE) {
1247  continue;
1248  }
1249 
1250  if (
1251  (Qr[pX] < 0 || Qr[pX] > *widthL || FLT_EQ(Qr[pX], *widthL)) ||
1252  (Qr[pY] < 0 || Qr[pY] > *heightL || FLT_EQ(Qr[pY], *heightL))
1253  ) {
1254  continue;
1255  }
1256 
1257  if (hasnodataS == FALSE)
1258  valL = 1;
1259  else if (rt_band_get_pixel(bandL, Qr[pX], Qr[pY], &valL, &isnodataL) != ES_NONE)
1260  continue;
1261 
1262  if ((hasnodataL == FALSE) || !isnodataL) {
1263  RASTER_DEBUG(3, "The two rasters do intersect");
1264  *intersects = 1;
1265  return ES_NONE;
1266  }
1267  }
1268  }
1269  }
1270  }
1271  }
1272  RASTER_DEBUG(4, "Smaller raster not in the other raster's pixel. Continuing");
1273  }
1274 
1275  RASTER_DEBUG(4, "Testing smaller raster vs larger raster");
1276  *intersects = rt_raster_intersects_algorithm(
1277  rastS, rastL,
1278  bandS, bandL,
1279  hasnodataS, hasnodataL,
1280  nodataS, nodataL
1281  );
1282 
1283  if (*intersects) return ES_NONE;
1284 
1285  RASTER_DEBUG(4, "Testing larger raster vs smaller raster");
1286  *intersects = rt_raster_intersects_algorithm(
1287  rastL, rastS,
1288  bandL, bandS,
1289  hasnodataL, hasnodataS,
1290  nodataL, nodataS
1291  );
1292 
1293  if (*intersects) return ES_NONE;
1294 
1295  RASTER_DEBUG(3, "The two rasters do not intersect");
1296 
1297  *intersects = 0;
1298  return ES_NONE;
1299 }
1300 
Datum coveredby(PG_FUNCTION_ARGS)
double skewY
Definition: librtcore.h:2261
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:1063
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:2197
void rterror(const char *fmt,...)
Wrappers used for reporting errors and info.
Definition: rt_context.c:199
rt_errorstate
Enum definitions.
Definition: librtcore.h:191
double ipY
Definition: librtcore.h:2259
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:229
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's convex hull.
Definition: rt_geometry.c:803
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:2196
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:2260
#define RASTER_DEBUGF(level, msg,...)
Definition: librtcore.h:311
int32_t rt_raster_get_srid(rt_raster raster)
Get raster'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:307
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:2258
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:2256
LWGEOM * lwmpoly_as_lwgeom(const LWMPOLY *obj)
Definition: lwgeom.c:217
double scaleY
Definition: librtcore.h:2257