General Dataset Utilities#

Authors: Tom Vo & Stephen Po-Chedley

Updated: 04/16/24 [xcdat v0.6.1]

Overview#

This notebook demonstrates the use of general utility methods available in xcdat, including the reorientation of the longitude axis, centering of time coordinates using time bounds, and adding and getting bounds.

The data used in this example is pulled directly from the Earth System Grid Federation (ESGF).

Notebook Kernel Setup#

Users can install their own instance of xcdat and follow these examples using their own environment (e.g., with VS Code, Jupyter, Spyder, iPython) or enable xcdat with existing JupyterHub instances.

First, create the conda environment:

conda create -n xcdat_notebook_0.7.0 -c conda-forge xcdat=0.7.0 xesmf matplotlib ipython ipykernel cartopy nc-time-axis gsw-xarray jupyter

Then install the kernel from the xcdat_notebook_0.7.0 environment using ipykernel and name the kernel with the display name (e.g., xcdat_notebook_0.7.0):

python -m ipykernel install --user --name xcdat_notebook_0.7.0 --display-name xcdat_notebook_0.7.0

Then to select the kernel xcdat_notebook_0.7.0 in Jupyter to use this kernel.

[1]:
import xcdat as xc

Open a dataset#

Datasets can be opened and read using open_dataset() or open_mfdataset() (multi-file).

Related APIs: xcdat.open_dataset() & xcdat.open_mfdataset()

[2]:
# data used in this tutorial
dataset_links = [
    "https://esgf-data2.llnl.gov/thredds/dodsC/user_pub_work/E3SM/1_0/amip_1850_aeroF/1deg_atm_60-30km_ocean/atmos/180x360/time-series/mon/ens2/v3/TS_187001_189412.nc",
    "https://esgf-data2.llnl.gov/thredds/dodsC/user_pub_work/E3SM/1_0/amip_1850_aeroF/1deg_atm_60-30km_ocean/atmos/180x360/time-series/mon/ens2/v3/TS_189501_191912.nc",
]
[3]:
# NOTE: Opening a multi-file dataset will result in data variables to be dask
# arrays.
ds = xc.open_mfdataset(dataset_links)
# print dataset
ds
[3]:
<xarray.Dataset> Size: 156MB
Dimensions:    (lat: 180, lon: 360, nbnd: 2, time: 600)
Coordinates:
  * lat        (lat) float64 1kB -89.5 -88.5 -87.5 -86.5 ... 86.5 87.5 88.5 89.5
  * lon        (lon) float64 3kB 0.5 1.5 2.5 3.5 4.5 ... 356.5 357.5 358.5 359.5
  * time       (time) object 5kB 1870-02-01 00:00:00 ... 1920-01-01 00:00:00
Dimensions without coordinates: nbnd
Data variables:
    lat_bnds   (lat, nbnd) float64 3kB dask.array<chunksize=(180, 2), meta=np.ndarray>
    lon_bnds   (lon, nbnd) float64 6kB dask.array<chunksize=(360, 2), meta=np.ndarray>
    gw         (lat) float64 1kB dask.array<chunksize=(180,), meta=np.ndarray>
    time_bnds  (time, nbnd) object 10kB dask.array<chunksize=(300, 2), meta=np.ndarray>
    area       (lat, lon) float64 518kB dask.array<chunksize=(180, 360), meta=np.ndarray>
    TS         (time, lat, lon) float32 156MB dask.array<chunksize=(300, 180, 360), meta=np.ndarray>
Attributes: (12/21)
    ne:                              30
    np:                              4
    Conventions:                     CF-1.0
    source:                          CAM
    case:                            20180622.DECKv1b_A2_1850aeroF.ne30_oEC.e...
    title:                           UNSET
    ...                              ...
    remap_script:                    ncremap
    remap_hostname:                  acme1
    remap_version:                   4.9.6
    map_file:                        /export/zender1/data/maps/map_ne30np4_to...
    input_file:                      /p/user_pub/e3sm/baldwin32/workshop/amip...
    DODS_EXTRA.Unlimited_Dimension:  time

Reorient the longitude axis#

Longitude can be represented from 0 to 360 E or as 180 W to 180 E. xcdat allows you to convert between these axes systems.

  • Related API: xcdat.swap_lon_axis()

  • Alternative solution: xcdat.open_mfdataset(dataset_links, lon_orient=(-180, 180))

[4]:
ds.lon
[4]:
<xarray.DataArray 'lon' (lon: 360)> Size: 3kB
array([  0.5,   1.5,   2.5, ..., 357.5, 358.5, 359.5])
Coordinates:
  * lon      (lon) float64 3kB 0.5 1.5 2.5 3.5 4.5 ... 356.5 357.5 358.5 359.5
Attributes:
    long_name:      Longitude of Grid Cell Centers
    standard_name:  longitude
    units:          degrees_east
    axis:           X
    valid_min:      0.0
    valid_max:      360.0
    bounds:         lon_bnds
[5]:
ds2 = xc.swap_lon_axis(ds, to=(-180, 180))
[6]:
ds2.lon
[6]:
<xarray.DataArray 'lon' (lon: 360)> Size: 3kB
array([-179.5, -178.5, -177.5, ...,  177.5,  178.5,  179.5])
Coordinates:
  * lon      (lon) float64 3kB -179.5 -178.5 -177.5 -176.5 ... 177.5 178.5 179.5
Attributes:
    long_name:      Longitude of Grid Cell Centers
    standard_name:  longitude
    units:          degrees_east
    axis:           X
    valid_min:      0.0
    valid_max:      360.0
    bounds:         lon_bnds

Center the time coordinates#

A given point of time often represents some time period (e.g., a monthly average). In this situation, data providers sometimes record the time as the beginning, middle, or end of the period. center_times() places the time coordinate in the center of the time interval (using time bounds to determine the center of the period).

  • Related API: xcdat.center_times()

  • Alternative solution: xcdat.open_mfdataset(dataset_links, center_times=True)

The time bounds used for centering time coordinates:

[7]:
# We access the values with .values because it is a dask array.
ds.time_bnds.values
[7]:
array([[cftime.DatetimeNoLeap(1870, 1, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 2, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 2, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 3, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 3, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 4, 1, 0, 0, 0, 0, has_year_zero=True)],
       ...,
       [cftime.DatetimeNoLeap(1919, 10, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 11, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 11, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 12, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 12, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1920, 1, 1, 0, 0, 0, 0, has_year_zero=True)]],
      dtype=object)

Before centering time coordinates:

[8]:
ds.time
[8]:
<xarray.DataArray 'time' (time: 600)> Size: 5kB
array([cftime.DatetimeNoLeap(1870, 2, 1, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1870, 3, 1, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1870, 4, 1, 0, 0, 0, 0, has_year_zero=True), ...,
       cftime.DatetimeNoLeap(1919, 11, 1, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1919, 12, 1, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1920, 1, 1, 0, 0, 0, 0, has_year_zero=True)],
      dtype=object)
Coordinates:
  * time     (time) object 5kB 1870-02-01 00:00:00 ... 1920-01-01 00:00:00
Attributes:
    long_name:     time
    bounds:        time_bnds
    cell_methods:  time: mean

Now center the time coordinates:

[9]:
ds3 = xc.center_times(ds)

After centering time coordinates:

[10]:
ds3.time
[10]:
<xarray.DataArray 'time' (time: 600)> Size: 5kB
array([cftime.DatetimeNoLeap(1870, 1, 16, 12, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1870, 2, 15, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1870, 3, 16, 12, 0, 0, 0, has_year_zero=True),
       ...,
       cftime.DatetimeNoLeap(1919, 10, 16, 12, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1919, 11, 16, 0, 0, 0, 0, has_year_zero=True),
       cftime.DatetimeNoLeap(1919, 12, 16, 12, 0, 0, 0, has_year_zero=True)],
      dtype=object)
Coordinates:
  * time     (time) object 5kB 1870-01-16 12:00:00 ... 1919-12-16 12:00:00
Attributes:
    long_name:     time
    bounds:        time_bnds
    cell_methods:  time: mean

Add bounds#

Bounds are critical to many xcdat operations. For example, they are used in determining the weights in spatial or temporal averages and in regridding operations. add_bounds() will attempt to produce bounds if they do not exist in the original dataset.

  • Related API: xarray.Dataset.bounds.add_bounds()

  • Alternative solution: xcdat.open_mfdataset(dataset_links, add_bounds=["X", "Y", "T"])

    • (Assuming the file doesn’t already have bounds for your desired axis/axes)

[11]:
# We are dropping the existing bounds to demonstrate adding bounds.
# we are starting with the dataset with centered time points
ds4 = ds3.drop_vars("time_bnds")
[12]:
try:
    ds4.bounds.get_bounds("T")
except KeyError as e:
    print(e)
"No bounds data variables were found for the 'T' 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()`."

There are two options for adding time bounds. The midpoint method places bounds at the midpoints between time bounds and the frequency method creates bounds based on the time stamp of each time point and the frequency of the data. This is the midpoint method:

[13]:
# midpoint method
ds4 = ds4.bounds.add_time_bounds(method="midpoint")
# print results
ds4.bounds.get_bounds("T")
[13]:
<xarray.DataArray 'time_bnds' (time: 600, bnds: 2)> Size: 10kB
array([[cftime.DatetimeNoLeap(1870, 1, 1, 18, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 1, 31, 6, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 1, 31, 6, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 3, 1, 18, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 3, 1, 18, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 3, 31, 18, 0, 0, 0, has_year_zero=True)],
       ...,
       [cftime.DatetimeNoLeap(1919, 10, 1, 6, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 10, 31, 18, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 10, 31, 18, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 12, 1, 6, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 12, 1, 6, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 12, 31, 18, 0, 0, 0, has_year_zero=True)]],
      dtype=object)
Coordinates:
  * time     (time) object 5kB 1870-01-16 12:00:00 ... 1919-12-16 12:00:00
Dimensions without coordinates: bnds
Attributes:
    xcdat_bounds:  True

Notice that the midpoint method does not place the bounds between the last moment of month n and the first moment of month n+1. The frequency method was meant to try to infer the correct bounds by taking into account the time stamps and the frequency of the data. The frequency method (below) is what is used when add_bounds=["T"] is specified in open_dataset or open_mfdataset.

[14]:
# drop time bounds again
ds5 = ds4.drop_vars("time_bnds")
# timestamp / frequency method
ds5 = ds5.bounds.add_time_bounds(method="freq")
# print results
ds5.bounds.get_bounds("T")
[14]:
<xarray.DataArray 'time_bnds' (time: 600, bnds: 2)> Size: 10kB
array([[cftime.DatetimeNoLeap(1870, 1, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 2, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 2, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 3, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1870, 3, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1870, 4, 1, 0, 0, 0, 0, has_year_zero=True)],
       ...,
       [cftime.DatetimeNoLeap(1919, 10, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 11, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 11, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1919, 12, 1, 0, 0, 0, 0, has_year_zero=True)],
       [cftime.DatetimeNoLeap(1919, 12, 1, 0, 0, 0, 0, has_year_zero=True),
        cftime.DatetimeNoLeap(1920, 1, 1, 0, 0, 0, 0, has_year_zero=True)]],
      dtype=object)
Coordinates:
  * time     (time) object 5kB 1870-01-16 12:00:00 ... 1919-12-16 12:00:00
Dimensions without coordinates: bnds
Attributes:
    xcdat_bounds:  True

Note that ds.bounds.add_time_bounds(method="midpoint") is the same as ds.bounds.add_bounds("T"). The latter method can be used to add bounds to other axes (e.g., latitude) as show below.

[15]:
ds6 = ds.drop_vars("lat_bnds")
ds6 = ds6.bounds.add_bounds("Y")
ds6.lat_bnds
[15]:
<xarray.DataArray 'lat_bnds' (lat: 180, bnds: 2)> Size: 3kB
array([[-90., -89.],
       [-89., -88.],
       [-88., -87.],
       [-87., -86.],
       [-86., -85.],
       [-85., -84.],
       [-84., -83.],
       [-83., -82.],
       [-82., -81.],
       [-81., -80.],
       [-80., -79.],
       [-79., -78.],
       [-78., -77.],
       [-77., -76.],
       [-76., -75.],
       [-75., -74.],
       [-74., -73.],
       [-73., -72.],
       [-72., -71.],
       [-71., -70.],
...
       [ 70.,  71.],
       [ 71.,  72.],
       [ 72.,  73.],
       [ 73.,  74.],
       [ 74.,  75.],
       [ 75.,  76.],
       [ 76.,  77.],
       [ 77.,  78.],
       [ 78.,  79.],
       [ 79.,  80.],
       [ 80.,  81.],
       [ 81.,  82.],
       [ 82.,  83.],
       [ 83.,  84.],
       [ 84.,  85.],
       [ 85.,  86.],
       [ 86.,  87.],
       [ 87.,  88.],
       [ 88.,  89.],
       [ 89.,  90.]])
Coordinates:
  * lat      (lat) float64 1kB -89.5 -88.5 -87.5 -86.5 ... 86.5 87.5 88.5 89.5
Dimensions without coordinates: bnds
Attributes:
    xcdat_bounds:  True

Add missing bounds for all axes supported by xcdat (X, Y, T, Z)#

[16]:
# We drop the dataset axes bounds to demonstrate generating missing bounds.
ds7 = ds.drop_vars(["time_bnds", "lat_bnds", "lon_bnds"])
[17]:
ds7
[17]:
<xarray.Dataset> Size: 156MB
Dimensions:  (lat: 180, lon: 360, time: 600)
Coordinates:
  * lat      (lat) float64 1kB -89.5 -88.5 -87.5 -86.5 ... 86.5 87.5 88.5 89.5
  * lon      (lon) float64 3kB 0.5 1.5 2.5 3.5 4.5 ... 356.5 357.5 358.5 359.5
  * time     (time) object 5kB 1870-02-01 00:00:00 ... 1920-01-01 00:00:00
Data variables:
    gw       (lat) float64 1kB dask.array<chunksize=(180,), meta=np.ndarray>
    area     (lat, lon) float64 518kB dask.array<chunksize=(180, 360), meta=np.ndarray>
    TS       (time, lat, lon) float32 156MB dask.array<chunksize=(300, 180, 360), meta=np.ndarray>
Attributes: (12/21)
    ne:                              30
    np:                              4
    Conventions:                     CF-1.0
    source:                          CAM
    case:                            20180622.DECKv1b_A2_1850aeroF.ne30_oEC.e...
    title:                           UNSET
    ...                              ...
    remap_script:                    ncremap
    remap_hostname:                  acme1
    remap_version:                   4.9.6
    map_file:                        /export/zender1/data/maps/map_ne30np4_to...
    input_file:                      /p/user_pub/e3sm/baldwin32/workshop/amip...
    DODS_EXTRA.Unlimited_Dimension:  time
[18]:
# add now-missing bounds
ds7 = ds7.bounds.add_missing_bounds(["X", "Y", "T"])
# print dataset
ds7
[18]:
<xarray.Dataset> Size: 156MB
Dimensions:    (lat: 180, lon: 360, time: 600, bnds: 2)
Coordinates:
  * lat        (lat) float64 1kB -89.5 -88.5 -87.5 -86.5 ... 86.5 87.5 88.5 89.5
  * lon        (lon) float64 3kB 0.5 1.5 2.5 3.5 4.5 ... 356.5 357.5 358.5 359.5
  * time       (time) object 5kB 1870-02-01 00:00:00 ... 1920-01-01 00:00:00
Dimensions without coordinates: bnds
Data variables:
    gw         (lat) float64 1kB dask.array<chunksize=(180,), meta=np.ndarray>
    area       (lat, lon) float64 518kB dask.array<chunksize=(180, 360), meta=np.ndarray>
    TS         (time, lat, lon) float32 156MB dask.array<chunksize=(300, 180, 360), meta=np.ndarray>
    lon_bnds   (lon, bnds) float64 6kB 0.0 1.0 1.0 2.0 ... 359.0 359.0 360.0
    lat_bnds   (lat, bnds) float64 3kB -90.0 -89.0 -89.0 ... 89.0 89.0 90.0
    time_bnds  (time, bnds) object 10kB 1870-02-01 00:00:00 ... 1920-02-01 00...
Attributes: (12/21)
    ne:                              30
    np:                              4
    Conventions:                     CF-1.0
    source:                          CAM
    case:                            20180622.DECKv1b_A2_1850aeroF.ne30_oEC.e...
    title:                           UNSET
    ...                              ...
    remap_script:                    ncremap
    remap_hostname:                  acme1
    remap_version:                   4.9.6
    map_file:                        /export/zender1/data/maps/map_ne30np4_to...
    input_file:                      /p/user_pub/e3sm/baldwin32/workshop/amip...
    DODS_EXTRA.Unlimited_Dimension:  time

Note that ds.bounds.add_missing_bounds uses ds.bounds.add_bounds for the latitude and longitude axes and defaults to the frequency method and add_time_bounds for the time axis. If you click on the database symbol for time_bnds above, the bounds are slightly mis-aligned because the time axis was not centered before adding the time axis. In this case, the user should call xcdat.center_times and then ds.bounds.add_missing_bounds (as shown earlier).

Get the dimension coordinates for an axis.#

In xarray, you can get a dimension coordinates by directly referencing its name (e.g., ds.lat). xcdat provides an alternative way to get dimension coordinates agnostically by simply passing the CF axis key to applicable APIs.

Helpful knowledge:

  • This API uses cf_xarray to interpret CF axis names and coordinate names in the xarray object attributes. Refer to Metadata Interpretation for more information.

Xarray documentation on coordinates (source):

  • There are two types of coordinates in xarray:

    • dimension coordinates are one dimensional coordinates with a name equal to their sole dimension (marked by * when printing a dataset or data array). They are used for label based indexing and alignment, like the index found on a pandas DataFrame or Series. Indeed, these “dimension” coordinates use a pandas.Index internally to store their values.

    • non-dimension coordinates are variables that contain coordinate data, but are not a dimension coordinate. They can be multidimensional (see Working with Multidimensional Coordinates), and there is no relationship between the name of a non-dimension coordinate and the name(s) of its dimension(s). Non-dimension coordinates can be useful for indexing or plotting; otherwise, xarray does not make any direct use of the values associated with them. They are not used for alignment or automatic indexing, nor are they required to match when doing arithmetic (see Coordinates).

  • Xarray’s terminology differs from the CF terminology, where the “dimension coordinates” are called “coordinate variables”, and the “non-dimension coordinates” are called “auxiliary coordinate variables” (see GH1295 for more details).

1. axis attr#

[19]:
ds.lat.attrs["axis"]
[19]:
'Y'

2. standard_name attr#

[20]:
ds.lat.attrs["standard_name"]
[20]:
'latitude'
[21]:
"lat" in ds.dims
[21]:
True

Utilities to get the coordinate axis and coordinate axis key#

[22]:
xc.get_dim_coords(ds, axis="Y")
[22]:
<xarray.DataArray 'lat' (lat: 180)> Size: 1kB
array([-89.5, -88.5, -87.5, -86.5, -85.5, -84.5, -83.5, -82.5, -81.5, -80.5,
       -79.5, -78.5, -77.5, -76.5, -75.5, -74.5, -73.5, -72.5, -71.5, -70.5,
       -69.5, -68.5, -67.5, -66.5, -65.5, -64.5, -63.5, -62.5, -61.5, -60.5,
       -59.5, -58.5, -57.5, -56.5, -55.5, -54.5, -53.5, -52.5, -51.5, -50.5,
       -49.5, -48.5, -47.5, -46.5, -45.5, -44.5, -43.5, -42.5, -41.5, -40.5,
       -39.5, -38.5, -37.5, -36.5, -35.5, -34.5, -33.5, -32.5, -31.5, -30.5,
       -29.5, -28.5, -27.5, -26.5, -25.5, -24.5, -23.5, -22.5, -21.5, -20.5,
       -19.5, -18.5, -17.5, -16.5, -15.5, -14.5, -13.5, -12.5, -11.5, -10.5,
        -9.5,  -8.5,  -7.5,  -6.5,  -5.5,  -4.5,  -3.5,  -2.5,  -1.5,  -0.5,
         0.5,   1.5,   2.5,   3.5,   4.5,   5.5,   6.5,   7.5,   8.5,   9.5,
        10.5,  11.5,  12.5,  13.5,  14.5,  15.5,  16.5,  17.5,  18.5,  19.5,
        20.5,  21.5,  22.5,  23.5,  24.5,  25.5,  26.5,  27.5,  28.5,  29.5,
        30.5,  31.5,  32.5,  33.5,  34.5,  35.5,  36.5,  37.5,  38.5,  39.5,
        40.5,  41.5,  42.5,  43.5,  44.5,  45.5,  46.5,  47.5,  48.5,  49.5,
        50.5,  51.5,  52.5,  53.5,  54.5,  55.5,  56.5,  57.5,  58.5,  59.5,
        60.5,  61.5,  62.5,  63.5,  64.5,  65.5,  66.5,  67.5,  68.5,  69.5,
        70.5,  71.5,  72.5,  73.5,  74.5,  75.5,  76.5,  77.5,  78.5,  79.5,
        80.5,  81.5,  82.5,  83.5,  84.5,  85.5,  86.5,  87.5,  88.5,  89.5])
Coordinates:
  * lat      (lat) float64 1kB -89.5 -88.5 -87.5 -86.5 ... 86.5 87.5 88.5 89.5
Attributes:
    long_name:      Latitude of Grid Cell Centers
    standard_name:  latitude
    units:          degrees_north
    axis:           Y
    valid_min:      -90.0
    valid_max:      90.0
    bounds:         lat_bnds
[23]:
xc.get_dim_keys(ds, axis="X")
[23]:
'lon'