Skip to content

Commit

Permalink
feat: add from_box to create rectangular MOCs
Browse files Browse the repository at this point in the history
  • Loading branch information
ManonMarchand committed May 29, 2024
1 parent a32c9a1 commit d1be81f
Show file tree
Hide file tree
Showing 4 changed files with 101 additions and 0 deletions.
1 change: 1 addition & 0 deletions CHANGELOG.md
Original file line number Diff line number Diff line change
Expand Up @@ -28,6 +28,7 @@ a column name to sum. It returns the sum of the column in the intersection betwe
MOC and the Multi-order-map. `MOC.probability_in_multiordermap` has a similar
behavior but also converts a probability-density into a probability.
* `STMOC.new_empty()` allows to create a new empty Space-Time MOC.
* `MOC.from_box` to create rectangular MOCs

## [0.13.1]

Expand Down
52 changes: 52 additions & 0 deletions python/mocpy/moc/moc.py
Original file line number Diff line number Diff line change
Expand Up @@ -1031,6 +1031,58 @@ def from_cone(cls, lon, lat, radius, max_depth, delta_depth=2):
np.uint8(delta_depth),
)
return cls(index)

@classmethod
@validate_lonlat
def from_box(cls, lon, lat, a, b, angle, max_depth):
"""
Create a MOC from a cone.
The cone is centered around the (`lon`, `lat`) position with a radius expressed by
`radius`.
Parameters
----------
lon : `~astropy.coordinates.Longitude` or its supertype `~astropy.units.Quantity`
The longitude of the center of the cone.
lat : `~astropy.coordinates.Latitude` or its supertype `~astropy.units.Quantity`
The latitude of the center of the cone.
a : `~astropy.coordinates.Angle`
The semi-major axis of the box, i.e. half of the box's width.
b : `~astropy.coordinates.Angle`
The semi-minor axis of the box, i.e. half of the box's height.
angle : `astropy.coordinates.Angle`
The tilt angle between the north and the semi-major axis, east of north.
max_depth : int
Maximum HEALPix cell resolution.
Returns
-------
result : `~mocpy.moc.MOC`
The resulting MOC
Examples
--------
>>> from mocpy import MOC
>>> from astropy.coordinates import Angle, Longitude, Latitude
>>> moc = MOC.from_box(
... lon=Longitude("0d"),
... lat=Latitude("0d"),
... a=Angle("10d"),
... b=Angle("2d"),
... angle=Angle("30d"),
... max_depth=10
... )
"""
index = mocpy.from_box(
lon[0],
lat[0],
np.float64(a.to_value(u.deg)),
np.float64(b.to_value(u.deg)),
np.float64(angle.to_value(u.deg)),
np.uint8(max_depth)
)
return cls(index)

@classmethod
@validate_lonlat
Expand Down
12 changes: 12 additions & 0 deletions python/mocpy/tests/test_moc.py
Original file line number Diff line number Diff line change
@@ -1,4 +1,5 @@
import copy
import math
import re

import astropy.units as u
Expand Down Expand Up @@ -525,6 +526,17 @@ def test_from_ring():
max_depth=10,
)

def test_from_box():
a = Angle("10d")
b = Angle("2d")
moc = MOC.from_box(
lon=Longitude("0d"), lat=Latitude("0d"),
a=a, b=b, angle=Angle("30d"),max_depth=10)
area = moc.sky_fraction * 4 * math.pi * u.steradian
# the moc covers a slightly bigger area than the region defined by the
# parameters
assert area.to(u.deg**2).value > 80
assert area.to(u.deg**2).value < 90

# TODO: IMPROVE THE ALGO
"""
Expand Down
36 changes: 36 additions & 0 deletions src/lib.rs
Original file line number Diff line number Diff line change
Expand Up @@ -158,6 +158,42 @@ fn mocpy(_py: Python, m: &PyModule) -> PyResult<()> {
.map_err(PyValueError::new_err)
}

/// # Input
/// - `lon` the longitude of the center of the box, in radians
/// - `lat` the latitude of the center of the box, in radians
/// - `a` the semi-major axis of the box (half the box width), in radians
/// - `b` the semi-minor axis of the box (half the box height), in radians
/// - `pa` the position angle (i.e. the angle between the north and the semi-major axis, east-of-north), in radians
/// - `depth`: the MOC depth
/// - `selection`: select BMOC cells to keep in the MOC
///
/// # Panics
/// - if `a` not in `]0, pi/2]`
/// - if `b` not in `]0, a]`
/// - if `pa` not in `[0, pi[`
///
#[pyfn(m)]
pub fn from_box(
lon: f64,
lat: f64,
a: f64,
b: f64,
angle: f64,
depth: u8,
) -> PyResult<usize> {
U64MocStore::get_global_store()
.from_box(
lon,
lat,
a,
b,
angle,
depth,
CellSelection::All
)
.map_err(PyValueError::new_err)
}

/// Create a spatial coverage from a given ring.
///
/// # Arguments
Expand Down

0 comments on commit d1be81f

Please sign in to comment.