Skip to content

Commit dc3379f

Browse files
authored
Fix resample interpolation coordinate mapping (#1202) (#1204)
The interpolation path used scipy.ndimage.zoom's edge-aligned formula which maps output pixel o to input pixel o*(N_in-1)/(N_out-1). This pinned the first output pixel to the first input pixel and the last to the last, creating a systematic half-pixel spatial shift vs. the block-centered output coordinate metadata. Switch all four backend paths (numpy, cupy, dask+numpy, dask+cupy) to use map_coordinates with the block-centered formula: input_pixel = (out_pixel + 0.5) * (N_in / N_out) - 0.5 This brings the interpolation path into alignment with both the aggregation path and the output coordinate labels.
1 parent 3339ea5 commit dc3379f

File tree

2 files changed

+120
-43
lines changed

2 files changed

+120
-43
lines changed

xrspatial/resample.py

Lines changed: 65 additions & 43 deletions
Original file line numberDiff line numberDiff line change
@@ -10,7 +10,6 @@
1010
import numpy as np
1111
import xarray as xr
1212
from scipy.ndimage import map_coordinates as _scipy_map_coords
13-
from scipy.ndimage import zoom as _scipy_zoom
1413

1514
try:
1615
import dask.array as da
@@ -19,10 +18,8 @@
1918

2019
try:
2120
import cupy
22-
import cupyx.scipy.ndimage as _cupy_ndimage
2321
except ImportError:
2422
cupy = None
25-
_cupy_ndimage = None
2623

2724
from xrspatial.utils import (
2825
ArrayTypeFunctionMapping,
@@ -64,52 +61,88 @@ def _output_chunks(in_chunks, scale):
6461
for i in range(len(in_chunks)))
6562

6663

67-
# -- NaN-aware zoom (NumPy) --------------------------------------------------
64+
# -- Block-centered coordinate mapping ---------------------------------------
6865

69-
def _nan_aware_zoom_np(data, zoom_yx, order):
70-
"""``scipy.ndimage.zoom`` with NaN-aware weighting.
66+
def _block_centered_coords(n_in, n_out):
67+
"""Return input coordinates for each output pixel using block-centered mapping.
68+
69+
Maps output pixel ``o`` to input pixel ``(o + 0.5) * (n_in / n_out) - 0.5``.
70+
This places each output pixel at the center of its spatial footprint,
71+
matching the convention used by ``_new_coords`` for output coordinate
72+
metadata.
73+
"""
74+
o = np.arange(n_out, dtype=np.float64)
75+
return (o + 0.5) * (n_in / n_out) - 0.5
76+
77+
78+
# -- NaN-aware interpolation (NumPy) ----------------------------------------
79+
80+
def _nan_aware_interp_np(data, out_h, out_w, order):
81+
"""Interpolate *data* to *(out_h, out_w)* with NaN-aware weighting.
82+
83+
Uses ``scipy.ndimage.map_coordinates`` with block-centered coordinate
84+
mapping so that sample positions match the output coordinate metadata.
7185
7286
For *order* 0 (nearest-neighbour) NaN propagates naturally.
7387
For higher orders the zero-fill / weight-mask trick is used so that
7488
NaN pixels do not corrupt their neighbours.
7589
"""
90+
iy = _block_centered_coords(data.shape[0], out_h)
91+
ix = _block_centered_coords(data.shape[1], out_w)
92+
yy, xx = np.meshgrid(iy, ix, indexing='ij')
93+
coords = np.array([yy.ravel(), xx.ravel()])
94+
7695
if order == 0:
77-
return _scipy_zoom(data, zoom_yx, order=0, mode='nearest')
96+
result = _scipy_map_coords(data, coords, order=0, mode='nearest')
97+
return result.reshape(out_h, out_w)
7898

7999
mask = np.isnan(data)
80100
if not mask.any():
81-
return _scipy_zoom(data, zoom_yx, order=order, mode='nearest')
101+
result = _scipy_map_coords(data, coords, order=order, mode='nearest')
102+
return result.reshape(out_h, out_w)
82103

83104
filled = np.where(mask, 0.0, data)
84105
weights = (~mask).astype(data.dtype)
85106

86-
z_data = _scipy_zoom(filled, zoom_yx, order=order, mode='nearest')
87-
z_wt = _scipy_zoom(weights, zoom_yx, order=order, mode='nearest')
107+
z_data = _scipy_map_coords(filled, coords, order=order, mode='nearest')
108+
z_wt = _scipy_map_coords(weights, coords, order=order, mode='nearest')
109+
110+
result = np.where(z_wt > 0.01,
111+
z_data / np.maximum(z_wt, 1e-10),
112+
np.nan)
113+
return result.reshape(out_h, out_w)
88114

89-
return np.where(z_wt > 0.01,
90-
z_data / np.maximum(z_wt, 1e-10),
91-
np.nan)
92115

116+
# -- NaN-aware interpolation (CuPy) -----------------------------------------
93117

94-
# -- NaN-aware zoom (CuPy) ---------------------------------------------------
118+
def _nan_aware_interp_cupy(data, out_h, out_w, order):
119+
"""CuPy variant of :func:`_nan_aware_interp_np`."""
120+
from cupyx.scipy.ndimage import map_coordinates as _cupy_map_coords
121+
122+
iy = cupy.asarray(_block_centered_coords(data.shape[0], out_h))
123+
ix = cupy.asarray(_block_centered_coords(data.shape[1], out_w))
124+
yy, xx = cupy.meshgrid(iy, ix, indexing='ij')
125+
coords = cupy.array([yy.ravel(), xx.ravel()])
95126

96-
def _nan_aware_zoom_cupy(data, zoom_yx, order):
97127
if order == 0:
98-
return _cupy_ndimage.zoom(data, zoom_yx, order=0, mode='nearest')
128+
result = _cupy_map_coords(data, coords, order=0, mode='nearest')
129+
return result.reshape(out_h, out_w)
99130

100131
mask = cupy.isnan(data)
101132
if not mask.any():
102-
return _cupy_ndimage.zoom(data, zoom_yx, order=order, mode='nearest')
133+
result = _cupy_map_coords(data, coords, order=order, mode='nearest')
134+
return result.reshape(out_h, out_w)
103135

104136
filled = cupy.where(mask, 0.0, data)
105137
weights = (~mask).astype(data.dtype)
106138

107-
z_data = _cupy_ndimage.zoom(filled, zoom_yx, order=order, mode='nearest')
108-
z_wt = _cupy_ndimage.zoom(weights, zoom_yx, order=order, mode='nearest')
139+
z_data = _cupy_map_coords(filled, coords, order=order, mode='nearest')
140+
z_wt = _cupy_map_coords(weights, coords, order=order, mode='nearest')
109141

110-
return cupy.where(z_wt > 0.01,
111-
z_data / cupy.maximum(z_wt, 1e-10),
112-
cupy.nan)
142+
result = cupy.where(z_wt > 0.01,
143+
z_data / cupy.maximum(z_wt, 1e-10),
144+
cupy.nan)
145+
return result.reshape(out_h, out_w)
113146

114147

115148
# -- Block-aggregation kernels (NumPy, numba) --------------------------------
@@ -283,13 +316,9 @@ def _interp_block_np(block, global_in_h, global_in_w,
283316
oy = np.arange(cum_out_y[yi], cum_out_y[yi + 1], dtype=np.float64)
284317
ox = np.arange(cum_out_x[xi], cum_out_x[xi + 1], dtype=np.float64)
285318

286-
# Map to global input coordinates (same formula scipy.ndimage.zoom uses)
287-
iy = (oy * (global_in_h - 1) / (global_out_h - 1)
288-
if global_out_h > 1
289-
else np.full(len(oy), (global_in_h - 1) / 2.0))
290-
ix = (ox * (global_in_w - 1) / (global_out_w - 1)
291-
if global_out_w > 1
292-
else np.full(len(ox), (global_in_w - 1) / 2.0))
319+
# Map to global input coordinates using block-centered formula
320+
iy = (oy + 0.5) * (global_in_h / global_out_h) - 0.5
321+
ix = (ox + 0.5) * (global_in_w / global_out_w) - 0.5
293322

294323
# Convert to local block coordinates (overlap shifts the origin)
295324
iy_local = iy - (cum_in_y[yi] - depth)
@@ -331,12 +360,9 @@ def _interp_block_cupy(block, global_in_h, global_in_w,
331360
ox = cupy.arange(int(cum_out_x[xi]), int(cum_out_x[xi + 1]),
332361
dtype=cupy.float64)
333362

334-
iy = (oy * (global_in_h - 1) / (global_out_h - 1)
335-
if global_out_h > 1
336-
else cupy.full(len(oy), (global_in_h - 1) / 2.0))
337-
ix = (ox * (global_in_w - 1) / (global_out_w - 1)
338-
if global_out_w > 1
339-
else cupy.full(len(ox), (global_in_w - 1) / 2.0))
363+
# Map to global input coordinates using block-centered formula
364+
iy = (oy + 0.5) * (global_in_h / global_out_h) - 0.5
365+
ix = (ox + 0.5) * (global_in_w / global_out_w) - 0.5
340366

341367
iy_local = iy - float(cum_in_y[yi] - depth)
342368
ix_local = ix - float(cum_in_x[xi] - depth)
@@ -410,10 +436,8 @@ def _run_numpy(data, scale_y, scale_x, method):
410436
out_h, out_w = _output_shape(*data.shape, scale_y, scale_x)
411437

412438
if method in INTERP_METHODS:
413-
zy = out_h / data.shape[0]
414-
zx = out_w / data.shape[1]
415-
return _nan_aware_zoom_np(data, (zy, zx),
416-
INTERP_METHODS[method]).astype(np.float32)
439+
return _nan_aware_interp_np(data, out_h, out_w,
440+
INTERP_METHODS[method]).astype(np.float32)
417441

418442
return _AGG_FUNCS[method](data, out_h, out_w).astype(np.float32)
419443

@@ -423,10 +447,8 @@ def _run_cupy(data, scale_y, scale_x, method):
423447
out_h, out_w = _output_shape(*data.shape, scale_y, scale_x)
424448

425449
if method in INTERP_METHODS:
426-
zy = out_h / data.shape[0]
427-
zx = out_w / data.shape[1]
428-
return _nan_aware_zoom_cupy(data, (zy, zx),
429-
INTERP_METHODS[method]).astype(cupy.float32)
450+
return _nan_aware_interp_cupy(data, out_h, out_w,
451+
INTERP_METHODS[method]).astype(cupy.float32)
430452

431453
# Aggregate: GPU reshape+reduce for integer factors, CPU fallback otherwise
432454
fy, fx = data.shape[0] / out_h, data.shape[1] / out_w

xrspatial/tests/test_resample.py

Lines changed: 55 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -182,6 +182,48 @@ def test_bilinear_upsample_smooth(self, grid_8x8):
182182
# Verify interior is within tolerance of the linear gradient
183183
assert np.all(np.isfinite(out.values))
184184

185+
@pytest.mark.parametrize('method', ['nearest', 'bilinear', 'cubic'])
186+
def test_interp_coordinate_alignment_downsample(self, method):
187+
"""Interpolated values should match output coordinate labels (#1202).
188+
189+
On a linear gradient where value == x-coordinate, a correct
190+
block-centered resample produces output values equal to the output
191+
coordinate labels (within floating-point tolerance).
192+
"""
193+
data = np.tile(np.arange(8, dtype=np.float32), (8, 1))
194+
agg = create_test_raster(data, attrs={'res': (1.0, 1.0)},
195+
dims=['y', 'x'])
196+
out = resample(agg, scale_factor=0.5, method=method)
197+
198+
# Output x-coords are block-centered: 0.5, 2.5, 4.5, 6.5
199+
# Values should match because the input is a linear gradient
200+
np.testing.assert_allclose(
201+
out.values[0], out.x.values, atol=0.6,
202+
err_msg=f"{method}: values should be close to x-coordinates"
203+
)
204+
# For bilinear on a linear gradient, the match should be exact
205+
if method == 'bilinear':
206+
np.testing.assert_allclose(
207+
out.values[0], out.x.values, atol=1e-5,
208+
err_msg="bilinear on linear gradient must be exact"
209+
)
210+
211+
def test_bilinear_coordinate_alignment_upsample(self):
212+
"""Upsampled interior pixels should match coordinates on a gradient."""
213+
data = np.tile(np.arange(8, dtype=np.float32), (8, 1))
214+
agg = create_test_raster(data, attrs={'res': (1.0, 1.0)},
215+
dims=['y', 'x'])
216+
out = resample(agg, scale_factor=2.0, method='bilinear')
217+
218+
# Interior pixels (away from boundary clamping) should be exact
219+
# for bilinear on a linear gradient. Skip first and last pixel
220+
# which may be clamped by mode='nearest' boundary handling.
221+
interior = slice(1, -1)
222+
np.testing.assert_allclose(
223+
out.values[0, interior], out.x.values[interior], atol=1e-4,
224+
err_msg="bilinear: interior values should match x-coordinates"
225+
)
226+
185227

186228
# ---------------------------------------------------------------------------
187229
# NaN handling
@@ -281,6 +323,19 @@ def test_interp_parity(self, numpy_and_dask_rasters, method, sf):
281323
np.testing.assert_allclose(dk_out.values, np_out.values,
282324
atol=1e-5, equal_nan=True)
283325

326+
@pytest.mark.parametrize('method', ['bilinear'])
327+
def test_dask_coordinate_alignment(self, method):
328+
"""Dask bilinear on a linear gradient should match coordinates (#1202)."""
329+
data = np.tile(np.arange(20, dtype=np.float32), (20, 1))
330+
dk_agg = create_test_raster(data, backend='dask+numpy',
331+
attrs={'res': (1.0, 1.0)},
332+
chunks=(8, 8))
333+
out = resample(dk_agg, scale_factor=0.5, method=method)
334+
np.testing.assert_allclose(
335+
out.values[0], out.x.values, atol=1e-4,
336+
err_msg="dask bilinear values should match x-coordinates"
337+
)
338+
284339
@pytest.mark.parametrize('method', ['average', 'min', 'max', 'median', 'mode'])
285340
def test_aggregate_parity(self, numpy_and_dask_rasters, method):
286341
np_agg, dk_agg = numpy_and_dask_rasters

0 commit comments

Comments
 (0)