Skip to content

Commit a71ce34

Browse files
Make polygonize work for polygons with holes (#1794)
Fixes #1793 # Description In case of a grid such as: ``` 1 1 1 1 2 1 1 1 1 ``` ``polygonize`` would not set the hole in polygon 1 where polygon 2 is located. Instead polygon 1 would overlap 2 entirely, as we forgot to pass the interiors that were computed by shapely. This fixes this and adds a unit test. # Checklist <!--- Before requesting review, please go through this checklist: --> - [x] Links to correct issue - [x] Update changelog, if changes affect users - [x] PR title starts with ``Issue #nr``, e.g. ``Issue #737`` - [x] Unit tests were added - [ ] **If feature added**: Added/extended example - [ ] **If feature added**: Added feature to API documentation - [ ] **If pixi.lock was changed**: Ran `pixi run generate-sbom` and committed changes
1 parent f52cf7d commit a71ce34

3 files changed

Lines changed: 53 additions & 1 deletion

File tree

docs/api/changelog.rst

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -55,6 +55,7 @@ Fixed
5555
layer, when the stage and bottom_elevation were exactly equal to the bottom of
5656
a layer in the model discretization, which would cause these cells to be
5757
dropped when distributing conductances later.
58+
- Fixed :func:`imod.prepare.spatial.polygonize` for polygons with holes.
5859

5960
Changed
6061
~~~~~~~

imod/tests/test_spatial.py

Lines changed: 45 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -7,6 +7,7 @@
77
import shapely.geometry as sg
88
import xarray as xr
99
from pytest_cases import parametrize_with_cases
10+
from shapely.geometry import Polygon
1011

1112
import imod
1213
from imod.testing import assert_frame_equal
@@ -118,6 +119,50 @@ def test_polygonize():
118119
assert sorted(gdf["value"]) == [1.0, 2.0, 3.0]
119120

120121

122+
def fits_exact_hole(parent: Polygon, child: Polygon, tol=0.0) -> bool:
123+
hole_polys = Polygon(parent.interiors[0])
124+
exterior = Polygon(child.exterior)
125+
return exterior.equals(hole_polys)
126+
127+
128+
def test_polygonize_with_holes():
129+
"""
130+
Polygonize should be able to handle polygons with holes. In this test, we
131+
create a 5x5 grid of values as follows:
132+
133+
1 1 1 1 1
134+
1 2 2 2 1
135+
1 2 3 2 1
136+
1 2 2 2 1
137+
1 1 1 1 1
138+
139+
The polygon for the value 1.0 should be a large square that contains the
140+
entire grid, but has a hole in the middle that corresponds to the polygon
141+
for the value 2.0. The polygon for the value 2.0 should have a hole that
142+
corresponds exactly to the polygon for the value 3.0.
143+
"""
144+
# Arrange
145+
nrow, ncol = 5, 5
146+
dx, dy = 1.0, -1.0
147+
xmin, xmax = 0.0, 5.0
148+
ymin, ymax = 0.0, 5.0
149+
coords = imod.util.spatial._xycoords((xmin, xmax, ymin, ymax), (dx, dy))
150+
kwargs = {"name": "test", "coords": coords, "dims": ("y", "x")}
151+
data = np.ones((nrow, ncol), dtype=np.float32)
152+
data[1:-1, 1:-1] = 2.0
153+
data[2, 2] = 3.0
154+
da = xr.DataArray(data, **kwargs)
155+
# Act
156+
gdf = imod.prepare.polygonize(da)
157+
gdf = gdf.sort_values("value").reset_index(drop=True)
158+
# Assert
159+
assert len(gdf) == 3
160+
assert sorted(gdf["value"]) == [1.0, 2.0, 3.0]
161+
p1, p2, p3 = gdf.geometry
162+
ok = fits_exact_hole(p1, p2) and fits_exact_hole(p2, p3)
163+
assert ok
164+
165+
121166
def test_handle_dtype():
122167
assert imod.prepare.spatial._handle_dtype(np.uint8, None) == (1, 0)
123168
assert imod.prepare.spatial._handle_dtype(np.uint16, None) == (2, 0)

imod/util/spatial.py

Lines changed: 7 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -730,7 +730,13 @@ def _polygonize(da: xr.DataArray) -> "gpd.GeoDataFrame":
730730
geometries = []
731731
colvalues = []
732732
for geom, colval in shapes:
733-
geometries.append(shapely.geometry.Polygon(geom["coordinates"][0]))
733+
coordinates = geom["coordinates"]
734+
exterior = coordinates[0]
735+
if len(coordinates) > 1:
736+
interiors = coordinates[1:]
737+
else:
738+
interiors = None
739+
geometries.append(shapely.geometry.Polygon(exterior, interiors))
734740
colvalues.append(colval)
735741

736742
gdf = gpd.GeoDataFrame({"value": colvalues, "geometry": geometries})

0 commit comments

Comments
 (0)