Skip to content

Guard _spline TPS against unbounded memory allocations (#1372)#1374

Merged
brendancol merged 1 commit into
mainfrom
issue-1372
Apr 29, 2026
Merged

Guard _spline TPS against unbounded memory allocations (#1372)#1374
brendancol merged 1 commit into
mainfrom
issue-1372

Conversation

@brendancol
Copy link
Copy Markdown
Contributor

Summary

Closes #1372.

spline() (Thin Plate Spline) builds N x N kernel matrices and an (N+3) x (N+3) augmented system before any guard runs. With ~10000 points the K matrix alone is ~800 MB; with ~40000 points it is ~13 GB. A user with a moderately large point set hits this silently — the process either OOMs or numpy raises a raw error without naming the input that drove the allocation.

This adds _check_spline_memory(n_points, grid_pixels) that estimates the peak allocation and raises MemoryError against 80% of _available_memory_bytes(). Same pattern as the kriging guard added in #1309 (issue #1307).

Estimate is the larger of:

  • kernel build: 4 * N**2 * 8 bytes (dx, dy, r2, K each N x N float64)
  • solve peak: 2 * (N+3)**2 * 8 bytes (A plus its LU factorization copy)

The TPS evaluator iterates the grid cell-by-cell inside a numba kernel, so prediction does not materialize a grid x N matrix. The guard only depends on N.

The 0.8 fraction matches the kriging guard so users get the same behaviour across interpolators.

Test plan

  • pytest xrspatial/tests/test_interpolation.py — 38 passed (was 34, four new in TestSplineMemoryGuard)
  • kernel-block trigger: N=600 with 8 MB available raises MemoryError matching "kernel block K"
  • grid size irrelevant: 1000x1000 template with N=4 passes the guard at 8 MB available
  • small input passes at 16 MB available
  • helper-level test confirms the GB-formatted message names the culprit

Add `_check_spline_memory(n_points, grid_pixels)` mirroring the
kriging guard from #1309. The TPS solver builds N x N kernel
matrices and an (N+3) x (N+3) augmented system before any guard
runs, so a moderately large point set takes the process down with
no indication of which input was at fault.

Estimate is the larger of:
  - kernel build: 4 * N**2 * 8 bytes (dx, dy, r2, K)
  - solve peak:   2 * (N+3)**2 * 8 bytes (A plus its LU copy)

Raises MemoryError when the estimate exceeds 80% of available
memory, naming the offending allocation in the message.
@github-actions github-actions Bot added the performance PR touches performance-sensitive code label Apr 29, 2026
@brendancol brendancol merged commit 0804296 into main Apr 29, 2026
10 of 11 checks passed
@brendancol brendancol deleted the issue-1372 branch May 4, 2026 19:50
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.

spline TPS: unbounded N\xc3\x97N system has no memory guard

1 participant