Skip to content

Fix bounds_error=True ignored with 1D interpolation #4855

New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Merged
merged 6 commits into from
Feb 10, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
2 changes: 2 additions & 0 deletions doc/whats-new.rst
Original file line number Diff line number Diff line change
Expand Up @@ -98,6 +98,8 @@ Bug fixes
:py:meth:`Dataset.to_zarr` (:issue:`4783`, :pull:`4795`).
By `Julien Seguinot <https://github.com/juseg>`_.
- Add :py:meth:`Dataset.drop_isel` and :py:meth:`DataArray.drop_isel` (:issue:`4658`, :pull:`4819`). By `Daniel Mesejo <https://github.com/mesejo>`_.
- Ensure that :py:meth:`Dataset.interp` raises ``ValueError`` when interpolating outside coordinate range and ``bounds_error=True`` (:issue:`4854`, :pull:`4855`).
By `Leif Denby <https://github.com/leifdenby>`_.
- Fix time encoding bug associated with using cftime versions greater than
1.4.0 with xarray (:issue:`4870`, :pull:`4871`). By `Spencer Clark <https://github.com/spencerkclark>`_.

Expand Down
2 changes: 1 addition & 1 deletion xarray/core/missing.py
Original file line number Diff line number Diff line change
Expand Up @@ -154,7 +154,7 @@ def __init__(
yi,
kind=self.method,
fill_value=fill_value,
bounds_error=False,
bounds_error=bounds_error,
assume_sorted=assume_sorted,
copy=copy,
**self.cons_kwargs,
Expand Down
15 changes: 15 additions & 0 deletions xarray/tests/test_interp.py
Original file line number Diff line number Diff line change
Expand Up @@ -866,3 +866,18 @@ def test_interpolate_chunk_advanced(method):
z = z.chunk(3)
actual = da.interp(t=0.5, x=x, y=y, z=z, kwargs=kwargs, method=method)
assert_identical(actual, expected)


@requires_scipy
def test_interp1d_bounds_error():
"""Ensure exception on bounds error is raised if requested"""
da = xr.DataArray(
np.sin(0.3 * np.arange(4)),
[("time", np.arange(4))],
)

with pytest.raises(ValueError):
da.interp(time=3.5, kwargs=dict(bounds_error=True))

# default is to fill with nans, so this should pass
da.interp(time=3.5)