PostGIS  2.5.7dev-r@@SVN_REVISION@@

◆ rt_raster_intersects()

rt_errorstate rt_raster_intersects ( rt_raster  rast1,
int  nband1,
rt_raster  rast2,
int  nband2,
int *  intersects 
)

Return zero if error occurred in function.

Return ES_ERROR if error occurred in function.

Parameter intersects returns non-zero if two rasters intersect

Parameters
rast1: the first raster whose band will be tested
nband1: the 0-based band of raster rast1 to use if value is less than zero, bands are ignored. if nband1 gte zero, nband2 must be gte zero
rast2: the second raster whose band will be tested
nband2: the 0-based band of raster rast2 to use if value is less than zero, bands are ignored if nband2 gte zero, nband1 must be gte zero
intersects: non-zero value if the two rasters' bands intersects
Returns
ES_NONE if success, ES_ERROR if error

Definition at line 990 of file rt_spatial_relationship.c.

994  {
995  int i;
996  int within = 0;
997 
998  LWGEOM *hull[2] = {NULL};
999  GEOSGeometry *ghull[2] = {NULL};
1000 
1001  uint16_t width1;
1002  uint16_t height1;
1003  uint16_t width2;
1004  uint16_t height2;
1005  double area1;
1006  double area2;
1007  double pixarea1;
1008  double pixarea2;
1009  rt_raster rastS = NULL;
1010  rt_raster rastL = NULL;
1011  uint16_t *widthS = NULL;
1012  uint16_t *heightS = NULL;
1013  uint16_t *widthL = NULL;
1014  uint16_t *heightL = NULL;
1015  int nbandS;
1016  int nbandL;
1017  rt_band bandS = NULL;
1018  rt_band bandL = NULL;
1019  int hasnodataS = FALSE;
1020  int hasnodataL = FALSE;
1021  double nodataS = 0;
1022  double nodataL = 0;
1023  int isnodataS = 0;
1024  int isnodataL = 0;
1025  double gtS[6] = {0};
1026  double igtL[6] = {0};
1027 
1028  uint32_t row;
1029  uint32_t rowoffset;
1030  uint32_t col;
1031  uint32_t coloffset;
1032 
1033  enum line_points {X1, Y1, X2, Y2};
1034  enum point {pX, pY};
1035  double lineS[4];
1036  double Qr[2];
1037  double valS;
1038  double valL;
1039 
1040  RASTER_DEBUG(3, "Starting");
1041 
1042  assert(NULL != rast1);
1043  assert(NULL != rast2);
1044  assert(NULL != intersects);
1045 
1046  if (nband1 < 0 && nband2 < 0) {
1047  nband1 = -1;
1048  nband2 = -1;
1049  }
1050  else {
1051  assert(nband1 >= 0 && nband1 < rt_raster_get_num_bands(rast1));
1052  assert(nband2 >= 0 && nband2 < rt_raster_get_num_bands(rast2));
1053  }
1054 
1055  /* same srid */
1056  if (rt_raster_get_srid(rast1) != rt_raster_get_srid(rast2)) {
1057  rterror("rt_raster_intersects: The two rasters provided have different SRIDs");
1058  *intersects = 0;
1059  return ES_ERROR;
1060  }
1061 
1062  /* raster extents need to intersect */
1063  do {
1064  int rtn;
1065 
1066  initGEOS(rtinfo, lwgeom_geos_error);
1067 
1068  rtn = 1;
1069 
1070  if ((rt_raster_get_convex_hull(rast1, &(hull[0])) != ES_NONE) || !hull[0]) {
1071  break;
1072  }
1073  ghull[0] = (GEOSGeometry *) LWGEOM2GEOS(hull[0], 0);
1074  if (!ghull[0]) {
1075  lwgeom_free(hull[0]);
1076  break;
1077  }
1078 
1079  if ((rt_raster_get_convex_hull(rast2, &(hull[1])) != ES_NONE) || !hull[1]) {
1080  GEOSGeom_destroy(ghull[0]);
1081  lwgeom_free(hull[0]);
1082  break;
1083  }
1084  ghull[1] = (GEOSGeometry *) LWGEOM2GEOS(hull[1], 0);
1085  if (!ghull[0]) {
1086  GEOSGeom_destroy(ghull[0]);
1087  lwgeom_free(hull[1]);
1088  lwgeom_free(hull[0]);
1089  break;
1090  }
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");
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");
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 }
#define FALSE
Definition: dbfopen.c:168
GEOSGeometry * LWGEOM2GEOS(const LWGEOM *lwgeom, uint8_t autofix)
void lwgeom_geos_error(const char *fmt,...)
void lwgeom_free(LWGEOM *geom)
Definition: lwgeom.c:1144
void rterror(const char *fmt,...)
Wrappers used for reporting errors and info.
Definition: rt_context.c:199
#define RASTER_DEBUG(level, msg)
Definition: librtcore.h:295
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
int32_t rt_raster_get_srid(rt_raster raster)
Get raster's SRID.
Definition: rt_raster.c:356
#define RASTER_DEBUGF(level, msg,...)
Definition: librtcore.h:299
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:674
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
rt_errorstate rt_band_get_pixel(rt_band band, int x, int y, double *value, int *nodata)
Get pixel value.
Definition: rt_band.c:1221
int rt_band_get_isnodata_flag(rt_band band)
Get isnodata flag value.
Definition: rt_band.c:714
#define FLT_EQ(x, y)
Definition: librtcore.h:2234
double rt_raster_get_x_scale(rt_raster raster)
Get scale X in projection units.
Definition: rt_raster.c:150
@ ES_NONE
Definition: librtcore.h:180
@ ES_ERROR
Definition: librtcore.h:181
uint16_t rt_raster_get_num_bands(rt_raster raster)
Definition: rt_raster.c:372
rt_errorstate rt_raster_get_convex_hull(rt_raster raster, LWGEOM **hull)
Get raster's convex hull.
Definition: rt_geometry.c:803
uint16_t rt_raster_get_height(rt_raster raster)
Definition: rt_raster.c:129
rt_errorstate rt_band_get_nodata(rt_band band, double *nodata)
Get NODATA value.
Definition: rt_band.c:1730
uint16_t rt_raster_get_width(rt_raster raster)
Definition: rt_raster.c:121
double rt_raster_get_y_scale(rt_raster raster)
Get scale Y in projection units.
Definition: rt_raster.c:159
rt_band rt_raster_get_band(rt_raster raster, int bandNum)
Return Nth band, or NULL if unavailable.
Definition: rt_raster.c:381
Datum intersects(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)
unsigned int uint32_t
Definition: uthash.h:78

References ES_ERROR, ES_NONE, FALSE, FLT_EQ, intersects(), LWGEOM2GEOS(), lwgeom_free(), lwgeom_geos_error(), RASTER_DEBUG, RASTER_DEBUGF, rt_band_get_hasnodata_flag(), rt_band_get_isnodata_flag(), rt_band_get_nodata(), rt_band_get_pixel(), rt_raster_cell_to_geopoint(), rt_raster_geopoint_to_cell(), rt_raster_get_band(), rt_raster_get_convex_hull(), rt_raster_get_height(), rt_raster_get_num_bands(), rt_raster_get_srid(), rt_raster_get_width(), rt_raster_get_x_scale(), rt_raster_get_y_scale(), rt_raster_intersects_algorithm(), rterror(), and rtinfo().

Referenced by RASTER_intersects(), and test_raster_intersects().

Here is the call graph for this function:
Here is the caller graph for this function: