From 2851540da112ffdef710c9880a78af0017fa2f89 Mon Sep 17 00:00:00 2001 From: MARCHAND MANON Date: Tue, 10 Sep 2024 12:17:12 +0200 Subject: [PATCH] feat: add support of Regions in from_astropy_regions --- CHANGELOG.md | 1 + python/mocpy/moc/moc.py | 18 +++++++++++++++++- python/mocpy/tests/test_moc.py | 14 +++++++++----- 3 files changed, 27 insertions(+), 6 deletions(-) diff --git a/CHANGELOG.md b/CHANGELOG.md index dfbf7aef..f0d6cf59 100644 --- a/CHANGELOG.md +++ b/CHANGELOG.md @@ -18,6 +18,7 @@ and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0 ### Added +* Add support of `regions.Regions` [#163] * Creation of a MOC from a zone (defined by min/max ra and dec)`MOC.from_zone` * Creation of a single MOC from a lot of small cones is faster with the new option in `MOC.from_cones`: the keyword 'union_strategy' can now take the value 'small_cones'. diff --git a/python/mocpy/moc/moc.py b/python/mocpy/moc/moc.py index 963f8266..7898a90c 100644 --- a/python/mocpy/moc/moc.py +++ b/python/mocpy/moc/moc.py @@ -1844,7 +1844,7 @@ def from_astropy_regions(cls, region, max_depth: int): The supported sky regions are `~regions.CircleSkyRegion`, `~regions.CircleAnnulusSkyRegion`, `~regions.EllipseSkyRegion`, `~regions.RectangleSkyRegion`, `~regions.PolygonSkyRegion`, - `~regions.PointSkyRegion`. + `~regions.PointSkyRegion`, `~regions.Regions`. max_depth : int The maximum HEALPix cell resolution of the MOC. Should be comprised between 0 and 29. @@ -1853,6 +1853,14 @@ def from_astropy_regions(cls, region, max_depth: int): ------- `~mocpy.moc.MOC` + Notes + ----- + - For the `~regions.Regions`, the returned MOC will be the union of all the regions. + - For the `~regions.PolygonSkyRegion` and the `~regions.RectangleSkyRegion`, the MOC + will consider the sides to follow great circles on the sky sphere while in + astropy-regions the sides follow straight lines in the projected space (depending on + a given WCS, see issue https://github.com/astropy/regions/issues/564). + Examples -------- >>> from astropy.coordinates import SkyCoord @@ -1869,6 +1877,7 @@ def from_astropy_regions(cls, region, max_depth: int): regions.RectangleSkyRegion, regions.PolygonSkyRegion, regions.PointSkyRegion, + regions.Regions, ) if isinstance(region, regions.CircleSkyRegion): center = region.center.icrs @@ -1927,6 +1936,13 @@ def from_astropy_regions(cls, region, max_depth: int): return cls.from_polygon_skycoord(region.vertices, max_depth=max_depth) if isinstance(region, regions.PointSkyRegion): return cls.from_skycoords(region.center, max_norder=max_depth) + if isinstance(region, regions.Regions): + mocs = [ + cls.from_astropy_regions(reg, max_depth=max_depth) for reg in region + ] + if len(mocs) == 1: + return mocs[0] + return mocs[0].union(*mocs[1:]) # fastest multi-union raise ValueError( "'from_astropy_regions' does not support this region type." diff --git a/python/mocpy/tests/test_moc.py b/python/mocpy/tests/test_moc.py index f6538274..f26d7ca3 100644 --- a/python/mocpy/tests/test_moc.py +++ b/python/mocpy/tests/test_moc.py @@ -689,14 +689,18 @@ def test_from_astropy_regions(): # polygons vertices = SkyCoord([1, 2, 2], [1, 1, 2], unit="deg", frame="fk5") polygon = regions.PolygonSkyRegion(vertices) - moc = MOC.from_astropy_regions(polygon, max_depth=10) - assert all(moc.contains_skycoords(vertices)) + moc_polygon = MOC.from_astropy_regions(polygon, max_depth=10) + assert all(moc_polygon.contains_skycoords(vertices)) # points point = SkyCoord("+23h20m48.3s +61d12m06s") region_point = regions.PointSkyRegion(point) - moc = MOC.from_astropy_regions(region_point, max_depth=10) - assert moc.max_order == 10 - assert moc.contains_skycoords(point) + moc_point = MOC.from_astropy_regions(region_point, max_depth=10) + assert moc_point.max_order == 10 + assert moc_point.contains_skycoords(point) + # multi regions + multi_region = regions.Regions([region_point, polygon]) + moc = MOC.from_astropy_regions(multi_region, max_depth=10) + assert moc == moc_polygon + moc_point # TODO: IMPROVE THE ALGO