Skip to content

Commit 3a96c8a

Browse files
committed
Add FutureWarning for change in decode_timedelta behavior
1 parent 714e17d commit 3a96c8a

File tree

6 files changed

+86
-16
lines changed

6 files changed

+86
-16
lines changed

doc/whats-new.rst

Lines changed: 6 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -74,6 +74,12 @@ Breaking changes
7474

7575
Deprecations
7676
~~~~~~~~~~~~
77+
- In a future version of xarray decoding of variables into
78+
:py:class:`numpy.timedelta64` values will be disabled by default. To silence
79+
warnings associated with this, set ``decode_timedelta`` to ``True``,
80+
``False``, or a :py:class:`coders.CFTimedeltaCoder` instance when opening
81+
data (:issue:`1621`, :pull:`9966`). By `Spencer Clark
82+
<https://github.com/spencerkclark>`_.
7783

7884

7985
Bug fixes

xarray/coding/times.py

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1356,6 +1356,7 @@ def __init__(
13561356
time_unit: PDDatetimeUnitOptions = "ns",
13571357
) -> None:
13581358
self.time_unit = time_unit
1359+
self._emit_decode_timedelta_future_warning = False
13591360

13601361
def encode(self, variable: Variable, name: T_Name = None) -> Variable:
13611362
if np.issubdtype(variable.data.dtype, np.timedelta64):
@@ -1373,6 +1374,14 @@ def encode(self, variable: Variable, name: T_Name = None) -> Variable:
13731374
def decode(self, variable: Variable, name: T_Name = None) -> Variable:
13741375
units = variable.attrs.get("units", None)
13751376
if isinstance(units, str) and units in TIME_UNITS:
1377+
if self._emit_decode_timedelta_future_warning:
1378+
emit_user_level_warning(
1379+
"In a future version of xarray decode_timedelta will "
1380+
"default to False rather than None. To silence this "
1381+
"warning, set decode_timedelta to True, False, or a "
1382+
"'CFTimedeltaCoder' instance.",
1383+
FutureWarning,
1384+
)
13761385
dims, data, attrs, encoding = unpack_for_decoding(variable)
13771386

13781387
units = pop_to(attrs, encoding, "units")

xarray/conventions.py

Lines changed: 9 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
from __future__ import annotations
22

33
import itertools
4+
import warnings
45
from collections import defaultdict
56
from collections.abc import Hashable, Iterable, Mapping, MutableMapping
67
from typing import TYPE_CHECKING, Any, Literal, TypeVar, Union
@@ -172,6 +173,7 @@ def decode_cf_variable(
172173

173174
original_dtype = var.dtype
174175

176+
decode_timedelta_was_none = decode_timedelta is None
175177
if decode_timedelta is None:
176178
if isinstance(decode_times, CFDatetimeCoder):
177179
decode_timedelta = CFTimedeltaCoder(time_unit=decode_times.time_unit)
@@ -200,6 +202,9 @@ def decode_cf_variable(
200202
if decode_timedelta:
201203
if not isinstance(decode_timedelta, CFTimedeltaCoder):
202204
decode_timedelta = CFTimedeltaCoder()
205+
decode_timedelta._emit_decode_timedelta_future_warning = (
206+
decode_timedelta_was_none
207+
)
203208
var = decode_timedelta.decode(var, name=name)
204209
if decode_times:
205210
# remove checks after end of deprecation cycle
@@ -352,6 +357,10 @@ def decode_cf_variables(
352357
353358
See: decode_cf_variable
354359
"""
360+
# Only emit once instance of the decode_timedelta default change
361+
# FutureWarning. This can be removed once this change is made.
362+
warnings.filterwarnings("once", "decode_timedelta", FutureWarning)
363+
355364
dimensions_used_by = defaultdict(list)
356365
for v in variables.values():
357366
for d in v.dims:

xarray/tests/test_backends.py

Lines changed: 11 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -49,7 +49,7 @@
4949
from xarray.backends.pydap_ import PydapDataStore
5050
from xarray.backends.scipy_ import ScipyBackendEntrypoint
5151
from xarray.backends.zarr import ZarrStore
52-
from xarray.coders import CFDatetimeCoder
52+
from xarray.coders import CFDatetimeCoder, CFTimedeltaCoder
5353
from xarray.coding.cftime_offsets import cftime_range
5454
from xarray.coding.strings import check_vlen_dtype, create_vlen_dtype
5555
from xarray.coding.variables import SerializationWarning
@@ -639,7 +639,9 @@ def test_roundtrip_timedelta_data(self) -> None:
639639
# to support large ranges
640640
time_deltas = pd.to_timedelta(["1h", "2h", "NaT"]).as_unit("s") # type: ignore[arg-type, unused-ignore]
641641
expected = Dataset({"td": ("td", time_deltas), "td0": time_deltas[0]})
642-
with self.roundtrip(expected) as actual:
642+
with self.roundtrip(
643+
expected, open_kwargs={"decode_timedelta": CFTimedeltaCoder(time_unit="ns")}
644+
) as actual:
643645
assert_identical(expected, actual)
644646

645647
def test_roundtrip_float64_data(self) -> None:
@@ -3267,7 +3269,13 @@ def test_attributes(self, obj) -> None:
32673269
def test_chunked_datetime64_or_timedelta64(self, dtype) -> None:
32683270
# Generalized from @malmans2's test in PR #8253
32693271
original = create_test_data().astype(dtype).chunk(1)
3270-
with self.roundtrip(original, open_kwargs={"chunks": {}}) as actual:
3272+
with self.roundtrip(
3273+
original,
3274+
open_kwargs={
3275+
"chunks": {},
3276+
"decode_timedelta": CFTimedeltaCoder(time_unit="ns"),
3277+
},
3278+
) as actual:
32713279
for name, actual_var in actual.variables.items():
32723280
assert original[name].chunks == actual_var.chunks
32733281
assert original.chunks == actual.chunks

xarray/tests/test_coding_times.py

Lines changed: 35 additions & 11 deletions
Original file line numberDiff line numberDiff line change
@@ -1542,7 +1542,10 @@ def test_roundtrip_timedelta64_nanosecond_precision(
15421542

15431543
encoded_var = conventions.encode_cf_variable(var)
15441544
decoded_var = conventions.decode_cf_variable(
1545-
"foo", encoded_var, decode_times=CFDatetimeCoder(time_unit=time_unit)
1545+
"foo",
1546+
encoded_var,
1547+
decode_times=CFDatetimeCoder(time_unit=time_unit),
1548+
decode_timedelta=CFTimedeltaCoder(time_unit=time_unit),
15461549
)
15471550

15481551
assert_identical(var, decoded_var)
@@ -1569,7 +1572,9 @@ def test_roundtrip_timedelta64_nanosecond_precision_warning() -> None:
15691572
assert encoded_var.dtype == np.int64
15701573
assert encoded_var.attrs["units"] == needed_units
15711574
assert encoded_var.attrs["_FillValue"] == 20
1572-
decoded_var = conventions.decode_cf_variable("foo", encoded_var)
1575+
decoded_var = conventions.decode_cf_variable(
1576+
"foo", encoded_var, decode_timedelta=CFTimedeltaCoder(time_unit="ns")
1577+
)
15731578
assert_identical(var, decoded_var)
15741579
assert decoded_var.encoding["dtype"] == np.int64
15751580

@@ -1617,7 +1622,9 @@ def test_roundtrip_float_times(fill_value, times, units, encoded_values) -> None
16171622
assert encoded_var.attrs["units"] == units
16181623
assert encoded_var.attrs["_FillValue"] == fill_value
16191624

1620-
decoded_var = conventions.decode_cf_variable("foo", encoded_var)
1625+
decoded_var = conventions.decode_cf_variable(
1626+
"foo", encoded_var, decode_timedelta=CFTimedeltaCoder(time_unit="ns")
1627+
)
16211628
assert_identical(var, decoded_var)
16221629
assert decoded_var.encoding["units"] == units
16231630
assert decoded_var.encoding["_FillValue"] == fill_value
@@ -1808,7 +1815,9 @@ def test_encode_cf_timedelta_casting_value_error(use_dask) -> None:
18081815
with pytest.warns(UserWarning, match="Timedeltas can't be serialized"):
18091816
encoded = conventions.encode_cf_variable(variable)
18101817
assert encoded.attrs["units"] == "hours"
1811-
decoded = conventions.decode_cf_variable("name", encoded)
1818+
decoded = conventions.decode_cf_variable(
1819+
"name", encoded, decode_timedelta=CFTimedeltaCoder(time_unit="ns")
1820+
)
18121821
assert_equal(variable, decoded)
18131822
else:
18141823
with pytest.raises(ValueError, match="Not possible"):
@@ -1832,43 +1841,58 @@ def test_encode_cf_timedelta_casting_overflow_error(use_dask, dtype) -> None:
18321841

18331842

18341843
_DECODE_TIMEDELTA_TESTS = {
1835-
"default": (True, None, np.dtype("timedelta64[ns]")),
1836-
"decode_timdelta=False": (True, False, np.dtype("int64")),
1844+
"default": (True, None, np.dtype("timedelta64[ns]"), True),
1845+
"decode_timdelta=False": (True, False, np.dtype("int64"), False),
18371846
"inherit-time_unit-from-decode_times": (
18381847
CFDatetimeCoder(time_unit="s"),
18391848
None,
18401849
np.dtype("timedelta64[s]"),
1850+
True,
18411851
),
18421852
"set-time_unit-via-CFTimedeltaCoder-decode_times=True": (
18431853
True,
18441854
CFTimedeltaCoder(time_unit="s"),
18451855
np.dtype("timedelta64[s]"),
1856+
False,
18461857
),
18471858
"set-time_unit-via-CFTimedeltaCoder-decode_times=False": (
18481859
False,
18491860
CFTimedeltaCoder(time_unit="s"),
18501861
np.dtype("timedelta64[s]"),
1862+
False,
18511863
),
18521864
"override-time_unit-from-decode_times": (
18531865
CFDatetimeCoder(time_unit="ns"),
18541866
CFTimedeltaCoder(time_unit="s"),
18551867
np.dtype("timedelta64[s]"),
1868+
False,
18561869
),
18571870
}
18581871

18591872

18601873
@pytest.mark.parametrize(
1861-
("decode_times", "decode_timedelta", "expected_dtype"),
1874+
("decode_times", "decode_timedelta", "expected_dtype", "warns"),
18621875
list(_DECODE_TIMEDELTA_TESTS.values()),
18631876
ids=list(_DECODE_TIMEDELTA_TESTS.keys()),
18641877
)
1865-
def test_decode_timedelta(decode_times, decode_timedelta, expected_dtype) -> None:
1878+
def test_decode_timedelta(
1879+
decode_times, decode_timedelta, expected_dtype, warns
1880+
) -> None:
18661881
timedeltas = pd.timedelta_range(0, freq="d", periods=3)
18671882
var = Variable(["time"], timedeltas)
18681883
encoded = conventions.encode_cf_variable(var)
1869-
decoded = conventions.decode_cf_variable(
1870-
"foo", encoded, decode_times=decode_times, decode_timedelta=decode_timedelta
1871-
)
1884+
if warns:
1885+
with pytest.warns(FutureWarning, match="decode_timedelta"):
1886+
decoded = conventions.decode_cf_variable(
1887+
"foo",
1888+
encoded,
1889+
decode_times=decode_times,
1890+
decode_timedelta=decode_timedelta,
1891+
)
1892+
else:
1893+
decoded = conventions.decode_cf_variable(
1894+
"foo", encoded, decode_times=decode_times, decode_timedelta=decode_timedelta
1895+
)
18721896
if decode_timedelta is False:
18731897
assert_equal(encoded, decoded)
18741898
else:

xarray/tests/test_conventions.py

Lines changed: 16 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -18,7 +18,7 @@
1818
)
1919
from xarray.backends.common import WritableCFDataStore
2020
from xarray.backends.memory import InMemoryDataStore
21-
from xarray.coders import CFDatetimeCoder
21+
from xarray.coders import CFDatetimeCoder, CFTimedeltaCoder
2222
from xarray.conventions import decode_cf
2323
from xarray.testing import assert_identical
2424
from xarray.tests import (
@@ -536,7 +536,9 @@ def test_decode_cf_time_kwargs(self, time_unit) -> None:
536536
)
537537

538538
dsc = conventions.decode_cf(
539-
ds, decode_times=CFDatetimeCoder(time_unit=time_unit)
539+
ds,
540+
decode_times=CFDatetimeCoder(time_unit=time_unit),
541+
decode_timedelta=CFTimedeltaCoder(time_unit=time_unit),
540542
)
541543
assert dsc.timedelta.dtype == np.dtype(f"m8[{time_unit}]")
542544
assert dsc.time.dtype == np.dtype(f"M8[{time_unit}]")
@@ -655,3 +657,15 @@ def test_encode_cf_variable_with_vlen_dtype() -> None:
655657
encoded_v = conventions.encode_cf_variable(v)
656658
assert encoded_v.data.dtype.kind == "O"
657659
assert coding.strings.check_vlen_dtype(encoded_v.data.dtype) is str
660+
661+
662+
def test_decode_cf_variables_decode_timedelta_warning() -> None:
663+
v = Variable(["time"], [1, 2], attrs={"units": "seconds"})
664+
variables = {"a": v}
665+
666+
with warnings.catch_warnings():
667+
warnings.filterwarnings("error", "decode_timedelta", FutureWarning)
668+
conventions.decode_cf_variables(variables, {}, decode_timedelta=True)
669+
670+
with pytest.warns(FutureWarning, match="decode_timedelta"):
671+
conventions.decode_cf_variables(variables, {})

0 commit comments

Comments
 (0)