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

Refactor area boundary sides retrieval with _geographic_sides and _projection_sides methods #566

Merged
merged 9 commits into from
Dec 8, 2023
Merged
198 changes: 123 additions & 75 deletions pyresample/geometry.py
Original file line number Diff line number Diff line change
Expand Up @@ -284,7 +284,7 @@

def get_bbox_lonlats(self, vertices_per_side: Optional[int] = None, force_clockwise: bool = True,
frequency: Optional[int] = None) -> tuple:
"""Return the bounding box lons and lats.
"""Return the bounding box lons and lats sides.

Args:
vertices_per_side:
Expand Down Expand Up @@ -321,44 +321,49 @@
warnings.warn("The `frequency` argument is pending deprecation, use `vertices_per_side` instead",
PendingDeprecationWarning, stacklevel=2)
vertices_per_side = vertices_per_side or frequency
lons, lats = self._get_bbox_elements(self.get_lonlats, vertices_per_side)
lon_sides, lat_sides = self._get_geographic_sides(vertices_per_side=vertices_per_side)
if force_clockwise and not self._corner_is_clockwise(
lons[0][-2], lats[0][-2], lons[0][-1], lats[0][-1], lons[1][1], lats[1][1]):
lon_sides[0][-2], lat_sides[0][-2],
lon_sides[0][-1], lat_sides[0][-1],
lon_sides[1][1], lat_sides[1][1]):
# going counter-clockwise
lons, lats = self._reverse_boundaries(lons, lats)
return lons, lats
lon_sides, lat_sides = self._reverse_boundaries(lon_sides, lat_sides)
return lon_sides, lat_sides

def _get_bbox_elements(self, coord_fun, vertices_per_side: Optional[int] = None) -> tuple:
s1_slice, s2_slice, s3_slice, s4_slice = self._get_bbox_slices(vertices_per_side)
s1_dim1, s1_dim2 = coord_fun(data_slice=s1_slice)
s2_dim1, s2_dim2 = coord_fun(data_slice=s2_slice)
s3_dim1, s3_dim2 = coord_fun(data_slice=s3_slice)
s4_dim1, s4_dim2 = coord_fun(data_slice=s4_slice)
dim1, dim2 = zip(*[(s1_dim1.squeeze(), s1_dim2.squeeze()),
(s2_dim1.squeeze(), s2_dim2.squeeze()),
(s3_dim1.squeeze(), s3_dim2.squeeze()),
(s4_dim1.squeeze(), s4_dim2.squeeze())])
if hasattr(dim1[0], 'compute') and da is not None:
dim1, dim2 = da.compute(dim1, dim2)
clean_dim1, clean_dim2 = self._filter_bbox_nans(dim1, dim2)
return clean_dim1, clean_dim2

def _filter_bbox_nans(
def _get_sides(self, coord_fun, vertices_per_side) -> tuple[list[np.ndarray], list[np.ndarray]]:
"""Return the boundary sides."""
top_slice, right_slice, bottom_slice, left_slice = self._get_bbox_slices(vertices_per_side)
top_dim1, top_dim2 = coord_fun(data_slice=top_slice)
right_dim1, right_dim2 = coord_fun(data_slice=right_slice)
bottom_dim1, bottom_dim2 = coord_fun(data_slice=bottom_slice)
left_dim1, left_dim2 = coord_fun(data_slice=left_slice)
sides_dim1, sides_dim2 = zip(*[(top_dim1.squeeze(), top_dim2.squeeze()),
(right_dim1.squeeze(), right_dim2.squeeze()),
(bottom_dim1.squeeze(), bottom_dim2.squeeze()),
(left_dim1.squeeze(), left_dim2.squeeze())])
if hasattr(sides_dim1[0], 'compute') and da is not None:
sides_dim1, sides_dim2 = da.compute(sides_dim1, sides_dim2)
return self._filter_sides_nans(sides_dim1, sides_dim2)

def _filter_sides_nans(
self,
dim1_sides: tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray],
dim2_sides: tuple[np.ndarray, np.ndarray, np.ndarray, np.ndarray],
) -> tuple[list[np.ndarray], list[np.ndarray]]:
"""Remove nan and inf values present in each side."""
new_dim1_sides = []
new_dim2_sides = []

Check notice on line 355 in pyresample/geometry.py

View check run for this annotation

codefactor.io / CodeFactor

pyresample/geometry.py#L355

unresolved comment '# FIXME: ~(~np.isfinite(dim1_side) | ~np.isfinite(dim1_side))' (C100)
for dim1_side, dim2_side in zip(dim1_sides, dim2_sides):
# FIXME: ~(~np.isfinite(dim1_side) | ~np.isfinite(dim1_side))
Copy link
Member

Choose a reason for hiding this comment

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

Does something need to happen with this?

Copy link
Contributor Author

@ghiggi ghiggi Dec 8, 2023

Choose a reason for hiding this comment

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

Right now Inf that are present somewhere on the boundary can pass in the sides (and mess up stuffs in downstream tasks). But I can't set this now, because for all out-of-Earth boundary, removing every Inf and not retrieving boundary sides crash pykdtree and gradient.
I will fix this in the out-of-Earth boundary PR

is_valid_mask = ~(np.isnan(dim1_side) | np.isnan(dim2_side))
if not is_valid_mask.any():
raise ValueError("Can't compute swath bounding coordinates. At least one side is completely invalid.")
raise ValueError("Can't compute boundary coordinates. At least one side is completely invalid.")
new_dim1_sides.append(dim1_side[is_valid_mask])
new_dim2_sides.append(dim2_side[is_valid_mask])
return new_dim1_sides, new_dim2_sides

Check notice on line 364 in pyresample/geometry.py

View check run for this annotation

codefactor.io / CodeFactor

pyresample/geometry.py#L364

unresolved comment '# FIXME: This currently replicate values if heigh/width < row_num/col_num !' (C100)
def _get_bbox_slices(self, vertices_per_side):
# FIXME: This currently replicate values if heigh/width < row_num/col_num !
height, width = self.shape
if vertices_per_side is None:
row_num = height
Expand Down Expand Up @@ -441,6 +446,9 @@
warnings.warn("The `frequency` argument is pending deprecation, use `vertices_per_side` instead",
PendingDeprecationWarning, stacklevel=2)
vertices_per_side = vertices_per_side or frequency
# FIXME:
# - Here return SphericalBoundary ensuring correct vertices ordering
# - Deprecate get_bbox_lonlats and usage of force_clockwise
lon_sides, lat_sides = self.get_bbox_lonlats(vertices_per_side=vertices_per_side,
force_clockwise=force_clockwise)
return AreaBoundary.from_lonlat_sides(lon_sides, lat_sides)
Expand Down Expand Up @@ -500,7 +508,7 @@

@property
def corners(self):
"""Return the corners of the current area."""
"""Return the corners centroids of the current area."""
from pyresample.spherical_geometry import Coordinate
return [Coordinate(*self.get_lonlat(0, 0)),
Coordinate(*self.get_lonlat(0, -1)),
Expand Down Expand Up @@ -606,6 +614,32 @@
"""Compute the slice to read based on an `area_to_cover`."""
raise NotImplementedError

@property
def is_geostationary(self):
"""Whether this geometry is in a geostationary satellite projection or not."""
return False

def _get_geographic_sides(self, vertices_per_side: Optional[int] = None) -> tuple:
"""Return the geographic boundary sides of the current area.

Args:
vertices_per_side:
The number of points to provide for each side.
By default (None) the full width and height will be provided.
If the area object is an AreaDefinition with any corner out of the Earth disk
(i.e. full disc geostationary area, Robinson projection, polar projections, ...)
by default only 50 points are selected.
"""
# FIXME: Add logic for out-of-earth disk projections
if self.is_geostationary:
return self._get_geostationary_boundary_sides(vertices_per_side=vertices_per_side, coordinates="geographic")
sides_lons, sides_lats = self._get_sides(coord_fun=self.get_lonlats, vertices_per_side=vertices_per_side)
return sides_lons, sides_lats

def _get_geostationary_boundary_sides(self, vertices_per_side, coordinates):
class_name = self.__class__.__name__
raise NotImplementedError(f"'_get_geostationary_boundary_sides' is not implemented for {class_name}")

Check notice on line 641 in pyresample/geometry.py

View check run for this annotation

codefactor.io / CodeFactor

pyresample/geometry.py#L641

unresolved comment '# FIXME: Add logic for out-of-earth disk projections' (C100)

Check warning on line 641 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L640-L641

Added lines #L640 - L641 were not covered by tests


class CoordinateDefinition(BaseDefinition):
"""Base class for geometry definitions defined by lons and lats only."""
Expand Down Expand Up @@ -1563,8 +1597,16 @@
return False
return 'geostationary' in coord_operation.method_name.lower()

def _get_geo_boundary_sides(self, vertices_per_side=None):
"""Retrieve the boundary sides list for geostationary projections."""
def _get_geostationary_boundary_sides(self, vertices_per_side=None, coordinates="geographic"):
Copy link
Member

Choose a reason for hiding this comment

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

I'm not a huge fan of this keyword argument, but this is also a private method at the moment so maybe it's OK for now. I can see how the beginning and end of the method are the same so splitting it into two separate methods based on coordinate type would be annoying. @mraspaud other ideas?

Copy link
Contributor Author

Choose a reason for hiding this comment

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

I would suggest to merge meanwhile, and we can split this in two functions later on ...
With this PR merged, I can open other 3 parallel PRs.

"""Retrieve the boundary sides list for geostationary projections with out-of-Earth disk coordinates.

The boundary sides right (1) and side left (3) are set to length 2.
"""
# FIXME:
# - If vertices_per_side is too small, there is the risk to loose boundary side points
# at the intersection corners between the CRS bounds polygon and the area
# extent polygon (which could exclude relevant regions of the geos area).
# - After fixing this, evaluate nb_points required for FULL DISC and CONUS area !
# Define default frequency
if vertices_per_side is None:
vertices_per_side = 50
Expand All @@ -1574,53 +1616,34 @@
# Ensure an even number of vertices for side creation
if (vertices_per_side % 2) != 0:
vertices_per_side = vertices_per_side + 1
lons, lats = get_geostationary_bounding_box_in_lonlats(self, nb_points=vertices_per_side)
# Retrieve dummy sides for GEO (side1 and side3 always of length 2)
side02_step = int(vertices_per_side / 2) - 1
lon_sides = [lons[slice(0, side02_step + 1)],
lons[slice(side02_step, side02_step + 1 + 1)],
lons[slice(side02_step + 1, side02_step * 2 + 1 + 1)],
np.append(lons[side02_step * 2 + 1], lons[0])
]
lat_sides = [lats[slice(0, side02_step + 1)],
lats[slice(side02_step, side02_step + 1 + 1)],
lats[slice(side02_step + 1, side02_step * 2 + 1 + 1)],
np.append(lats[side02_step * 2 + 1], lats[0])
]
return lon_sides, lat_sides

def boundary(self, *, vertices_per_side=None, force_clockwise=False, frequency=None):
"""Retrieve the AreaBoundary object.

Parameters
----------
vertices_per_side:
The number of points to provide for each side. By default (None)
the full width and height will be provided, except for geostationary
projection where by default only 50 points are selected.
frequency:
Deprecated, use vertices_per_side
force_clockwise:
Perform minimal checks and reordering of coordinates to ensure
that the returned coordinates follow a clockwise direction.
This is important for compatibility with
:class:`pyresample.spherical.SphPolygon` where operations depend
on knowing the inside versus the outside of a polygon. These
operations assume that coordinates are clockwise.
Default is False.
"""
from pyresample.boundary import AreaBoundary
if frequency is not None:
warnings.warn("The `frequency` argument is pending deprecation, use `vertices_per_side` instead",
PendingDeprecationWarning, stacklevel=2)
vertices_per_side = vertices_per_side or frequency
if self.is_geostationary:
lon_sides, lat_sides = self._get_geo_boundary_sides(vertices_per_side=vertices_per_side)
# Retrieve coordinates (x,y) or (lon, lat)
if coordinates == "geographic":
x, y = get_geostationary_bounding_box_in_lonlats(self, nb_points=vertices_per_side)
else:
lon_sides, lat_sides = self.get_bbox_lonlats(vertices_per_side=vertices_per_side,
force_clockwise=force_clockwise)
boundary = AreaBoundary.from_lonlat_sides(lon_sides, lat_sides)
return boundary
x, y = get_geostationary_bounding_box_in_proj_coords(self, nb_points=vertices_per_side)

Check warning on line 1623 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L1623

Added line #L1623 was not covered by tests
# Ensure that a portion of the area is within the Earth disk.
if x.shape[0] < 4:
raise ValueError("The geostationary projection area is entirely out of the Earth disk.")

Check warning on line 1626 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L1626

Added line #L1626 was not covered by tests
# Retrieve dummy sides for GEO
# FIXME:
# - _get_geostationary_bounding_box_* does not guarantee to return nb_points and even points!
# - if odd nb_points, above can go out of index
# --> sides_x = self._get_dummy_sides(x, vertices_per_side=vertices_per_side)
# --> sides_y = self._get_dummy_sides(y, vertices_per_side=vertices_per_side)
side02_step = int(vertices_per_side / 2) - 1
sides_x = [
x[slice(0, side02_step + 1)],
x[slice(side02_step, side02_step + 1 + 1)],
x[slice(side02_step + 1, side02_step * 2 + 1 + 1)],
np.append(x[side02_step * 2 + 1], x[0])
]
sides_y = [
y[slice(0, side02_step + 1)],
y[slice(side02_step, side02_step + 1 + 1)],
y[slice(side02_step + 1, side02_step * 2 + 1 + 1)],
np.append(y[side02_step * 2 + 1], y[0])
]
return sides_x, sides_y

def get_edge_bbox_in_projection_coordinates(self, vertices_per_side: Optional[int] = None,
frequency: Optional[int] = None):
Expand All @@ -1629,7 +1652,7 @@
warnings.warn("The `frequency` argument is pending deprecation, use `vertices_per_side` instead",
PendingDeprecationWarning, stacklevel=2)
vertices_per_side = vertices_per_side or frequency
x, y = self._get_bbox_elements(self.get_proj_coords, vertices_per_side)
x, y = self._get_sides(self.get_proj_coords, vertices_per_side=vertices_per_side)
return np.hstack(x), np.hstack(y)

@property
Expand Down Expand Up @@ -2682,6 +2705,32 @@
# return np.max(np.concatenate(vert_dist, hor_dist)) # alternative to histogram
return res

def _get_projection_sides(self, vertices_per_side: Optional[int] = None) -> tuple:
"""Return the projection boundary sides of the current area.

Args:
vertices_per_side:
The number of points to provide for each side.
By default (None) the full width and height will be provided.
If the area object is an AreaDefinition with any corner out of the Earth disk
(i.e. full disc geostationary area, Robinson projection, polar projections, ...)
by default only 50 points are selected.

Returns:
The output structure is a tuple of two lists of four elements each.
The first list contains the projection x coordinates.
The second list contains the projection y coordinates.
Each list element is a numpy array representing a specific side of the geometry.
The order of the sides are [top", "right", "bottom", "left"]
"""
# FIXME: Add logic for out-of-earth-disk
if self.is_geostationary:
return self._get_geostationary_boundary_sides(vertices_per_side=vertices_per_side,

Check warning on line 2728 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L2727-L2728

Added lines #L2727 - L2728 were not covered by tests
coordinates="projection")
sides_x, sides_y = self._get_sides(coord_fun=self.get_proj_coords,

Check warning on line 2730 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L2730

Added line #L2730 was not covered by tests
vertices_per_side=vertices_per_side)
return sides_x, sides_y

Check warning on line 2732 in pyresample/geometry.py

View check run for this annotation

Codecov / codecov/patch

pyresample/geometry.py#L2732

Added line #L2732 was not covered by tests


def get_geostationary_angle_extent(geos_area):
"""Get the max earth (vs space) viewing angles in x and y."""
Expand Down Expand Up @@ -2719,19 +2768,19 @@
try:
x, y = intersection.boundary.xy
except NotImplementedError:
return [], []
return np.array([]), np.array([])
return np.asanyarray(x[:-1]), np.asanyarray(y[:-1])


def get_full_geostationary_bounding_box_in_proj_coords(geos_area, nb_points=50):
"""Get the bbox in geos projection coordinates of the full disk in `geos_area` projection.
"""Get the valid boundary geos projection coordinates of the full disk.

Args:
geos_area: Geostationary area definition to get the bounding box for.
nb_points: Number of points on the polygon
"""
x_max_angle, y_max_angle = get_geostationary_angle_extent(geos_area)
h = get_geostationary_height(geos_area.crs)

Check notice on line 2783 in pyresample/geometry.py

View check run for this annotation

codefactor.io / CodeFactor

pyresample/geometry.py#L2783

unresolved comment '# FIXME: Add logic for out-of-earth-disk' (C100)

# generate points around the north hemisphere in satellite projection
# make it a bit smaller so that we stay inside the valid area
Expand All @@ -2740,7 +2789,6 @@
y = -np.sin(points_around) * (y_max_angle - 0.0001)
x *= h
y *= h

return x, y


Expand Down
35 changes: 18 additions & 17 deletions pyresample/gradient/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -133,32 +133,33 @@
src_prj=src_crs, dst_prj=dst_crs)
self.prj = pyproj.Proj(self.target_geo_def.crs)

def _get_prj_poly(self, geo_def):
# - None if out of Earth Disk
# - False is SwathDefinition
if isinstance(geo_def, SwathDefinition):
return False
try:
poly = get_polygon(self.prj, geo_def)
except (NotImplementedError, ValueError): # out-of-earth disk area or any valid projected boundary coordinates
poly = None

Check warning on line 144 in pyresample/gradient/__init__.py

View check run for this annotation

Codecov / codecov/patch

pyresample/gradient/__init__.py#L143-L144

Added lines #L143 - L144 were not covered by tests
return poly

def _get_src_poly(self, src_y_start, src_y_end, src_x_start, src_x_end):
"""Get bounding polygon for source chunk."""
geo_def = self.source_geo_def[src_y_start:src_y_end,
src_x_start:src_x_end]
try:
src_poly = get_polygon(self.prj, geo_def)
except AttributeError:
# Can't create polygons for SwathDefinition
src_poly = False

return src_poly
return self._get_prj_poly(geo_def)

def _get_dst_poly(self, idx, dst_x_start, dst_x_end,
def _get_dst_poly(self, idx,
dst_x_start, dst_x_end,
dst_y_start, dst_y_end):
"""Get target chunk polygon."""
dst_poly = self.dst_polys.get(idx, None)
if dst_poly is None:
geo_def = self.target_geo_def[dst_y_start:dst_y_end,
dst_x_start:dst_x_end]
try:
dst_poly = get_polygon(self.prj, geo_def)
except AttributeError:
# Can't create polygons for SwathDefinition
dst_poly = False
dst_poly = self._get_prj_poly(geo_def)
self.dst_polys[idx] = dst_poly

return dst_poly

def get_chunk_mappings(self):
Expand Down Expand Up @@ -294,19 +295,20 @@
res = res.squeeze()

res = xr.DataArray(res, dims=data_dims, coords=coords)

return res


def check_overlap(src_poly, dst_poly):
"""Check if the two polygons overlap."""
# swath definition case
if dst_poly is False or src_poly is False:
covers = True
# area / area case
elif dst_poly is not None and src_poly is not None:
covers = src_poly.intersects(dst_poly)
# out of earth disk case
else:
covers = False

return covers


Expand All @@ -328,7 +330,6 @@
src_gradient_yl, src_gradient_yp,
dst_x, dst_y,
method=method)

return image


Expand Down
Loading
Loading