Skip to content
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

Use the same function to floatize coords in polyfit and polyval #9691

Merged
merged 7 commits into from
Oct 29, 2024
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 @@ -37,6 +37,8 @@ Bug fixes

- Fix inadvertent deep-copying of child data in DataTree.
By `Stephan Hoyer <https://github.com/shoyer>`_.
- Fix regression in the interoperability of :py:meth:`DataArray.polyfit` and :py:meth:`xr.polyval` for date-time coordinates. (:pull:`9691`).
By `Pascal Bourgault <https://github.com/aulemahal>`_.

Documentation
~~~~~~~~~~~~~
Expand Down
15 changes: 3 additions & 12 deletions xarray/core/dataset.py
Original file line number Diff line number Diff line change
Expand Up @@ -60,7 +60,7 @@
_contains_datetime_like_objects,
get_chunksizes,
)
from xarray.core.computation import unify_chunks
from xarray.core.computation import _ensure_numeric, unify_chunks
from xarray.core.coordinates import (
Coordinates,
DatasetCoordinates,
Expand All @@ -87,7 +87,6 @@
merge_coordinates_without_align,
merge_core,
)
from xarray.core.missing import _floatize_x
from xarray.core.options import OPTIONS, _get_keep_attrs
from xarray.core.types import (
Bins,
Expand Down Expand Up @@ -9066,22 +9065,14 @@ def polyfit(
variables = {}
skipna_da = skipna

x: Any = self.coords[dim].variable
x = _floatize_x((x,), (x,))[0][0]

try:
x = x.values.astype(np.float64)
except TypeError as e:
raise TypeError(
f"Dim {dim!r} must be castable to float64, got {type(x).__name__}."
) from e
x = np.asarray(_ensure_numeric(self.coords[dim]).astype(np.float64))

xname = f"{self[dim].name}_"
order = int(deg) + 1
lhs = np.vander(x, order)

if rcond is None:
rcond = x.shape[0] * np.finfo(x.dtype).eps
rcond = x.shape[0] * np.finfo(x.dtype).eps # type: ignore[assignment]

# Weights:
if w is not None:
Expand Down
30 changes: 29 additions & 1 deletion xarray/tests/test_dataset.py
Original file line number Diff line number Diff line change
Expand Up @@ -6720,6 +6720,34 @@ def test_polyfit_warnings(self) -> None:
ds.var1.polyfit("dim2", 10, full=True)
assert len(ws) == 1

def test_polyfit_polyval(self) -> None:
da = xr.DataArray(
np.arange(1, 10).astype(np.float64), dims=["x"], coords=dict(x=np.arange(9))
)

out = da.polyfit("x", 3, full=False)
da_fitval = xr.polyval(da.x, out.polyfit_coefficients)
# polyval introduces very small errors (1e-16 here)
xr.testing.assert_allclose(da_fitval, da)

da = da.assign_coords(x=xr.date_range("2001-01-01", periods=9, freq="YS"))
out = da.polyfit("x", 3, full=False)
da_fitval = xr.polyval(da.x, out.polyfit_coefficients)
xr.testing.assert_allclose(da_fitval, da, rtol=1e-3)

@requires_cftime
def test_polyfit_polyval_cftime(self) -> None:
da = xr.DataArray(
np.arange(1, 10).astype(np.float64),
dims=["x"],
coords=dict(
x=xr.date_range("2001-01-01", periods=9, freq="YS", calendar="noleap")
),
)
out = da.polyfit("x", 3, full=False)
da_fitval = xr.polyval(da.x, out.polyfit_coefficients)
np.testing.assert_allclose(da_fitval, da)

@staticmethod
def _test_data_var_interior(
original_data_var, padded_data_var, padded_dim_name, expected_pad_values
Expand Down Expand Up @@ -7230,7 +7258,7 @@ def test_differentiate_datetime(dask) -> None:
assert np.allclose(actual, 1.0)


@pytest.mark.skipif(not has_cftime, reason="Test requires cftime.")
@requires_cftime
@pytest.mark.parametrize("dask", [True, False])
def test_differentiate_cftime(dask) -> None:
rs = np.random.RandomState(42)
Expand Down
Loading