Skip to content

Make cubic resample prefilter explicit so chunk seams stay sub-eps#1478

Merged
brendancol merged 3 commits intomainfrom
issue-1464
May 5, 2026
Merged

Make cubic resample prefilter explicit so chunk seams stay sub-eps#1478
brendancol merged 3 commits intomainfrom
issue-1464

Conversation

@brendancol
Copy link
Copy Markdown
Contributor

Closes #1464.

Summary

scipy.ndimage.map_coordinates(order >= 2, prefilter=True) does three things behind the scenes for mode='nearest':

  1. Edge-pad the input by 12 pixels.
  2. Run spline_filter on the padded array.
  3. Shift the sample coordinates by 12 to account for the pad.

The whole step happens inside each map_coordinates call. In the NaN-aware path it runs twice per block (once for filled, once for weights). In the dask path it runs once per chunk on the overlapped block, with the prefilter's IIR boundary transient anchored to the chunk's overlapped edge rather than the global edge. With depth=10 overlap the transient leaked into chunk-interior samples at ~1e-6 near chunk seams.

This PR runs the prepad + spline_filter step ourselves, then calls map_coordinates(prefilter=False). Boundary handling is now the same between numpy and dask paths, and the per-block prefilter on filled and weights runs once each instead of twice.

The IIR transient still depends on overlap depth. The previous _INTERP_DEPTH['cubic'] = 10 left ~2e-6 residual, which is what surfaced this issue. Bumped to 16, which leaves ~7e-10 residual, below float32 epsilon, so the output rounds bit-identical between eager and chunked paths.

Same change in the cupy path with cupyx.scipy.ndimage.spline_filter.

Test plan

  • Existing 62 resample tests still pass
  • New TestCubicPrefilterParity suite (6 tests) asserts atol=1e-10 between eager numpy and dask+numpy for cubic resample of:
    • degree-2 polynomial (60x60, chunks=13)
    • random data (50x50, chunks=11)
    • random data at scale factors 0.5, 2.0, 0.7 (48x48, chunks=11)
    • random data with NaN (50x50, chunks=11)
  • Confirmed new tests fail without the fix: max diff was ~1e-6 on the random-data cases
  • Confirmed new tests pass with the fix at the 1e-10 bound

…1464)

scipy.ndimage.map_coordinates with order >= 2 silently runs spline_filter
on a 12-pixel edge-padded copy of the input.  Calling it implicitly per
block (and twice per block in the NaN-aware case, once for filled and
once for weights) couples the IIR boundary transient to the chunk
geometry; the eager and chunked paths could disagree by ~1e-6 near
seams.

Run the prepad+spline_filter step ourselves, shift the sample
coordinates by the pad, and pass prefilter=False to map_coordinates.
Mirror the change in the cupy path.  Bump cubic overlap depth from 10
to 16 so the residual transient (0.268^16 ~ 7e-10) rounds to zero in
the float32 output.

Add a TestCubicPrefilterParity suite that asserts atol=1e-10 between
eager numpy and dask+numpy paths on polynomial, random, and NaN inputs.
@github-actions github-actions Bot added the performance PR touches performance-sensitive code label May 4, 2026
@brendancol
Copy link
Copy Markdown
Contributor Author

@copilot resolve the merge conflicts in this pull request

… prefilter + updated threshold

Co-authored-by: brendancol <433221+brendancol@users.noreply.github.com>
Copy link
Copy Markdown
Contributor

Copilot AI commented May 5, 2026

@copilot resolve the merge conflicts in this pull request

Merge conflicts resolved in commit 5d358a4.

The three conflicts in resample.py were resolved by keeping the explicit prefilter path from this PR and adopting the majority-weight gate threshold of 0.5 (updated in main). The single conflict in test_resample.py was resolved by keeping both TestCubicPrefilterParity (from this PR) and all new test classes added in main (TestDaskCuPyParity, TestIntegerInput, TestTargetResolutionTuple). All 111 tests pass.

@brendancol
Copy link
Copy Markdown
Contributor Author

@copilot resolve the merge conflicts in this pull request

…eep explicit prefilter

Co-authored-by: brendancol <433221+brendancol@users.noreply.github.com>
Copy link
Copy Markdown
Contributor

Copilot AI commented May 5, 2026

@copilot resolve the merge conflicts in this pull request

Merge conflicts resolved again in commit 2d29663. The new commit f6a7759 on main ("Preserve input float dtype") changed weights = (~mask).astype(np.float64) to weights = (~mask).astype(block.dtype) in both the numpy and cupy NaN-aware paths. Both conflicts were resolved by keeping the dtype preservation from main and retaining the explicit prefilter path from this PR. All 133 tests pass.

@brendancol brendancol merged commit 401256f into main May 5, 2026
1 check passed
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

performance PR touches performance-sensitive code

Projects

None yet

Development

Successfully merging this pull request may close these issues.

Make cubic resample prefilter explicit so chunk seams stay sub-eps

2 participants