Skip to content
Open
3 changes: 3 additions & 0 deletions cf_xarray/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -9,4 +9,7 @@
from .options import set_options # noqa
from .utils import _get_version

from . import geometry # noqa

#
__version__ = _get_version()
63 changes: 63 additions & 0 deletions cf_xarray/geometry.py
Original file line number Diff line number Diff line change
Expand Up @@ -438,3 +438,66 @@ def cf_to_lines(ds: xr.Dataset):
geoms = np.where(np.diff(offset2) == 1, lines[offset2[:-1]], multilines)

return xr.DataArray(geoms, dims=part_node_count.dims, coords=part_node_count.coords)


def grid_to_polygons(ds: xr.Dataset) -> xr.DataArray:
"""
Converts a regular 2D lat/lon grid to a 2D array of shapely polygons.

Modified from https://notebooksharing.space/view/c6c1f3a7d0c260724115eaa2bf78f3738b275f7f633c1558639e7bbd75b31456.

Parameters
----------
ds : xr.Dataset
Dataset with "latitude" and "longitude" variables as well as their bounds variables.
1D "latitude" and "longitude" variables are supported. This function will automatically
broadcast them against each other.

Returns
-------
DataArray
DataArray with shapely polygon per grid cell.
"""
import shapely

grid = ds.cf[["latitude", "longitude"]].load()
bounds = grid.cf.bounds
dims = grid.cf.dims

if "latitude" in dims or "longitude" in dims:
# for 1D lat, lon, this allows them to be
# broadcast against each other
grid = grid.reset_coords()

assert "latitude" in bounds
assert "longitude" in bounds
(lon_bounds,) = bounds["longitude"]
(lat_bounds,) = bounds["latitude"]

with xr.set_options(keep_attrs=True):
(points,) = xr.broadcast(grid)

bounds_dim = grid.cf.get_bounds_dim_name("latitude")
points = points.transpose(..., bounds_dim)
lonbnd = points[lon_bounds].data
latbnd = points[lat_bounds].data

if points.sizes[bounds_dim] == 2:
lonbnd = lonbnd[..., [0, 0, 1, 1]]
latbnd = latbnd[..., [0, 1, 1, 0]]

elif points.sizes[bounds_dim] != 4:
raise ValueError(
f"The size of the detected bounds or vertex dimension {bounds_dim} is not 2 or 4."
)

# geopandas needs this
mask = lonbnd[..., 0] >= 180
lonbnd[mask, :] = lonbnd[mask, :] - 360
Comment on lines +1023 to +1025
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I think the reason for this is that most geographic data is given in WGS84 which is -180...180 (I may be wrong here). So if you want to compare it to regional polygons it is probably a good idea to wrap the data. Still, it may be confusing to users.

Maybe this could be optional. However, I am not sure what a good name is (In regionmask I call it wrap_lon with the options 180 (= what you do here), 360 and False. That works but there may be better names).


polyarray = shapely.polygons(shapely.linearrings(lonbnd, latbnd))

# 'geometry' is a blessed name in geopandas.
boxes = points[lon_bounds][..., 0].copy(data=polyarray).rename("geometry")

return boxes
3 changes: 1 addition & 2 deletions pyproject.toml
Original file line number Diff line number Diff line change
Expand Up @@ -140,8 +140,7 @@ module=[
"pint",
"matplotlib",
"pytest",
"shapely",
"shapely.geometry",
"shapely.*",
"xarray.core.pycompat",
]
ignore_missing_imports = true
Expand Down