Skip to content

Commit

Permalink
Fix datetime.timedelta casting bug in coding.times.infer_datetime_uni…
Browse files Browse the repository at this point in the history
…ts (#2128)

* Fix #2127

* Fix typo in time-series.rst

* Use pd.to_timedelta to convert to np.timedelta64 objects

* Install cftime through netcdf4 through pip

* box=False
  • Loading branch information
spencerkclark authored and shoyer committed May 14, 2018
1 parent d1b669e commit 188141f
Show file tree
Hide file tree
Showing 4 changed files with 10 additions and 4 deletions.
3 changes: 1 addition & 2 deletions ci/requirements-py27-windows.yml
Original file line number Diff line number Diff line change
Expand Up @@ -8,7 +8,6 @@ dependencies:
- h5py
- h5netcdf
- matplotlib
- netcdf4
- pathlib2
- pytest
- flake8
Expand All @@ -21,4 +20,4 @@ dependencies:
- rasterio
- zarr
- pip:
- cftime
- netcdf4
2 changes: 1 addition & 1 deletion doc/time-series.rst
Original file line number Diff line number Diff line change
Expand Up @@ -73,7 +73,7 @@ native representation of dates to those that fall between the years 1678 and
returned as arrays of ``cftime.datetime`` objects and a ``CFTimeIndex``
can be used for indexing. The ``CFTimeIndex`` enables only a subset of
the indexing functionality of a ``pandas.DatetimeIndex`` and is only enabled
when using standalone version of ``cftime`` (not the version packaged with
when using the standalone version of ``cftime`` (not the version packaged with
earlier versions ``netCDF4``). See :ref:`CFTimeIndex` for more information.

Datetime indexing
Expand Down
6 changes: 5 additions & 1 deletion xarray/coding/times.py
Original file line number Diff line number Diff line change
Expand Up @@ -253,7 +253,11 @@ def infer_datetime_units(dates):
else:
reference_date = dates[0] if len(dates) > 0 else '1970-01-01'
reference_date = format_cftime_datetime(reference_date)
unique_timedeltas = np.unique(np.diff(dates)).astype('timedelta64[ns]')
unique_timedeltas = np.unique(np.diff(dates))
if unique_timedeltas.dtype == np.dtype('O'):
# Convert to np.timedelta64 objects using pandas to work around a
# NumPy casting bug: https://github.com/numpy/numpy/issues/11096
unique_timedeltas = pd.to_timedelta(unique_timedeltas, box=False)
units = _infer_time_units_from_diff(unique_timedeltas)
return '%s since %s' % (units, reference_date)

Expand Down
3 changes: 3 additions & 0 deletions xarray/tests/test_coding_times.py
Original file line number Diff line number Diff line change
Expand Up @@ -587,6 +587,9 @@ def test_infer_cftime_datetime_units():
'seconds since 1900-01-01 00:00:00.000000'),
([date_type(1900, 1, 1),
date_type(1900, 1, 2, 0, 0, 0, 5)],
'days since 1900-01-01 00:00:00.000000'),
([date_type(1900, 1, 1), date_type(1900, 1, 8),
date_type(1900, 1, 16)],
'days since 1900-01-01 00:00:00.000000')]:
assert expected == coding.times.infer_datetime_units(dates)

Expand Down

0 comments on commit 188141f

Please sign in to comment.