diff --git a/docs/releases/development.rst b/docs/releases/development.rst index 53b6188e..6189f3c9 100644 --- a/docs/releases/development.rst +++ b/docs/releases/development.rst @@ -20,3 +20,6 @@ Next release (in development) Only the first frame is used to generate the clim to avoid loading more data than required. The new clim parameter allows users to specify the data limits explicitly if this is insufficient (:pr:`191`). +* Detect time coordinates with correct dtype and some standard attributes even without a `units` attribute. + Variables with `units` are preferred + (:issue:`190`, :pr:`193`). diff --git a/src/emsarray/conventions/_base.py b/src/emsarray/conventions/_base.py index e35c38f1..ccb00b2d 100644 --- a/src/emsarray/conventions/_base.py +++ b/src/emsarray/conventions/_base.py @@ -296,17 +296,31 @@ def time_coordinate(self) -> xarray.DataArray: ---------- .. [1] `CF Conventions v1.10, 4.4 Time Coordinate `_ """ + # First look for a datetime64 variable with a 'units' field in the encoding for name in self.dataset.variables.keys(): variable = self.dataset[name] - # xarray will automatically decode all time variables - # and move the 'units' attribute over to encoding to store this change. - if 'units' in variable.encoding: - units = variable.encoding['units'] - # A time variable must have units of the form ' since ' - if 'since' in units: - # The variable must now be a numpy datetime - if variable.dtype.type == numpy.datetime64: + # The variable must be a numpy datetime + if variable.dtype.type == numpy.datetime64: + # xarray will automatically decode all time variables + # and move the 'units' attribute over to encoding to store this change. + if 'units' in variable.encoding: + units = variable.encoding['units'] + # A time variable must have units of the form ' since ' + if 'since' in units: return variable + + # Next, look for any datetime64 variable with an appropriate attribute + for name in self.dataset.variables.keys(): + variable = self.dataset[name] + if variable.dtype.type == numpy.datetime64: + possible_attributes = { + 'coordinate_type': 'time', + 'standard_name': 'time', + 'axis': 'T', + } + if any(variable.attrs.get(name, None) == value for name, value in possible_attributes.items()): + return variable + raise NoSuchCoordinateError("Could not find time coordinate in dataset") @cached_property diff --git a/tests/conventions/test_base.py b/tests/conventions/test_base.py index 28200dd8..a459220f 100644 --- a/tests/conventions/test_base.py +++ b/tests/conventions/test_base.py @@ -146,7 +146,7 @@ def apply_clip_mask(self, clip_mask: xarray.Dataset, work_dir: Pathish) -> xarra return masking.mask_grid_dataset(self.dataset, clip_mask, work_dir) -def test_time_coordinate(datasets: pathlib.Path) -> None: +def test_time_coordinate_with_units(datasets: pathlib.Path) -> None: dataset = xarray.open_dataset(datasets / 'times.nc') SimpleConvention(dataset).bind() time_coordinate = dataset.ems.time_coordinate @@ -154,6 +154,39 @@ def test_time_coordinate(datasets: pathlib.Path) -> None: xarray.testing.assert_equal(time_coordinate, dataset['time']) +def test_time_coordinate_no_units(datasets: pathlib.Path) -> None: + dataset = xarray.Dataset({ + 'time': xarray.DataArray( + data=pandas.date_range('2025-09-08', '2025-10-08'), + attrs={'coordinate_type': 'time', 'standard_name': 'time'}, + ), + }) + + SimpleConvention(dataset).bind() + time_coordinate = dataset.ems.time_coordinate + assert time_coordinate.name == 'time' + xarray.testing.assert_equal(time_coordinate, dataset['time']) + + +def test_time_coordinate_with_units_first(datasets: pathlib.Path) -> None: + dataset = xarray.Dataset({ + 'time_no_units': xarray.DataArray( + data=pandas.date_range('2025-09-08', '2025-10-08'), + attrs={'coordinate_type': 'time', 'standard_name': 'time'}, + ), + 'time_with_units': xarray.DataArray( + data=pandas.date_range('2025-09-08', '2025-10-08'), + attrs={'coordinate_type': 'time', 'standard_name': 'time'}, + ), + }) + + dataset['time_with_units'].encoding['units'] = 'days since 1970-01-01' + + SimpleConvention(dataset).bind() + time_coordinate = dataset.ems.time_coordinate + assert time_coordinate.name == 'time_with_units' + + def test_time_coordinate_missing() -> None: dataset = xarray.Dataset() SimpleConvention(dataset).bind()