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

[Bug]: get_bounds() can't get bounds for non-CF axes without var_key argument #707

Closed
tomvothecoder opened this issue Oct 7, 2024 · 0 comments · Fixed by #708
Closed
Labels
type: bug Inconsistencies or issues which will cause an issue or problem for users or implementors.

Comments

@tomvothecoder
Copy link
Collaborator

tomvothecoder commented Oct 7, 2024

What happened?

Notice in the get_bounds() conditional below, we retrieve bounds for the Dataset object using _get_bounds_keys() (which does not work with non-CF axes). On the other hand, if we pass var_key, we get the bounds using the "bounds" attribute of non-CF axes (which works fine).

xcdat/xcdat/bounds.py

Lines 233 to 250 in e3dda64

if var_key is None:
# Get all bounds keys in the Dataset for this axis.
bounds_keys = self._get_bounds_keys(axis)
else:
# Get the obj in the Dataset using the key.
obj = _get_data_var(self._dataset, key=var_key)
# Check if the object is a data variable or a coordinate variable.
# If it is a data variable, derive the axis coordinate variable.
if obj.name in list(self._dataset.data_vars):
coord = get_dim_coords(obj, axis)
elif obj.name in list(self._dataset.coords):
coord = obj
try:
bounds_keys = [coord.attrs["bounds"]]
except KeyError:
bounds_keys = []

xcdat/xcdat/bounds.py

Lines 477 to 508 in e3dda64

def _get_bounds_keys(self, axis: CFAxisKey) -> List[str]:
"""Get bounds keys for an axis's coordinate variables in the dataset.
This function attempts to map bounds to an axis using ``cf_xarray``
and its interpretation of the CF "bounds" attribute.
Parameters
----------
axis : CFAxisKey
The CF axis key ("X", "Y", "T", or "Z").
Returns
-------
List[str]
The axis bounds key(s).
"""
cf_method = self._dataset.cf.bounds
cf_attrs = CF_ATTR_MAP[axis]
keys: List[str] = []
try:
keys = keys + cf_method[cf_attrs["axis"]]
except KeyError:
pass
try:
keys = cf_method[cf_attrs["coordinate"]]
except KeyError:
pass
return list(set(keys))

What did you expect to happen? Are there are possible answers you came across?

We should be able to get bounds for non-CF compliant axes that xCDAT can map to via dimension name (e.g., "lev", "lat", "lon"). We need to update _get_bounds_keys() to consider the "bounds" attribute.

Minimal Complete Verifiable Example (MVCE)

# %%
import xcdat as xc
import xarray as xr

# Create a dummy dataset
data = xr.DataArray([[1, 2, 3], [4, 5, 6]], dims=("lat", "lon"))
ds = xr.Dataset({"data": data})

# Add a "lev" Z axis with bounds
lev = xr.DataArray([1000, 2000], dims=("lev"))
ds["lev"] = lev
ds = ds.bounds.add_missing_bounds(axes="Z")

# %%
# Notice how we can map to "lev" without CF attributes because "lev" is a standard
# Z axis name.
xc.get_dim_coords(ds, "Z")

# %%
# However, we can't map to "lev_bnds".
ds.bounds.get_bounds(axis="Z")
# KeyError: "No bounds data variables were found for the 'Z' axis. Make sure the dataset has bound data vars and their names match the 'bounds' attributes found on their related time coordinate variables. Alternatively, you can add bounds with `ds.bounds.add_missing_bounds()` or `ds.bounds.add_bounds()`."

# %%
# Workaround 1 - add CF attribute(s)
ds_cf = ds.copy()
ds_cf["lev"].attrs["axis"] = "Z"
ds_cf.bounds.get_bounds(axis="Z")

# %%
# Workaround 2 - pass the var_key arg to get bounds via "bounds" key
ds.bounds.get_bounds(axis="Z", var_key="lev")

Relevant log output

No response

Anything else we need to know?

Found this bug through e3sm_diags GH issue 863

Environment

Latest xCDAT and main

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
type: bug Inconsistencies or issues which will cause an issue or problem for users or implementors.
Projects
Status: Done
Development

Successfully merging a pull request may close this issue.

1 participant