cdshealpix.nested.elliptical_cone_search¶
- cdshealpix.nested.elliptical_cone_search(lon, lat, a, b, pa, depth, delta_depth=2, flat=False)¶
- Get the HEALPix cells contained in an elliptical cone at a given depth. - This method is wrapped around the elliptical_cone_coverage_custom method from the cdshealpix Rust crate. - Parameters:
- lonastropy.coordinates.Longitude
- Longitude of the center of the elliptical cone. 
- latastropy.coordinates.Latitude
- Latitude of the center of the elliptical cone. 
- aastropy.coordinates.Angle
- Semi-major axe angle of the elliptical cone. 
- bastropy.coordinates.Angle
- Semi-minor axe angle of the elliptical cone. 
- paastropy.coordinates.Angle
- The position angle (i.e. the angle between the north and the semi-major axis, east-of-north). 
- depthint
- Maximum depth of the HEALPix cells that will be returned. 
- delta_depthint, optional
- To control the approximation, you can choose to perform the computations at a deeper depth using the - depth_deltaparameter. The depth at which the computations will be made will therefore be equal to- depth+- depth_delta.
- flatboolean, optional
- False by default (i.e. returns a consistent MOC). If True, the HEALPix cells returned will all be at depth indicated by - depth.
 
- lon
- Returns:
- ipix, depth, fully_covered(numpy.ndarray,numpy.ndarray,numpy.ndarray)
- A tuple containing 3 numpy arrays of identical size: - ipixstores HEALPix cell indices.
- depthstores HEALPix cell depths.
- fully_coveredstores flags on whether the HEALPix cells are fully covered by the elliptical cone.
 
 
- ipix, depth, fully_covered(
- Raises:
- ValueError
- If one of - lon,- lat,- major_axe,- minor_axeor- pacontains more that one value.
- ValueError
- If the semi-major axis - aexceeds 90deg (i.e. area of one hemisphere)
- ValueError
- If the semi-minor axis - bis greater than the semi-major axis- a
 
 - Examples - >>> from cdshealpix import elliptical_cone_search >>> from astropy.coordinates import Angle, SkyCoord, Longitude, Latitude >>> import astropy.units as u >>> import numpy as np >>> lon = Longitude(0, u.deg) >>> lat = Latitude(0, u.deg) >>> a = Angle(50, unit="deg") >>> b = Angle(10, unit="deg") >>> pa = Angle(45, unit="deg") >>> max_depth = 12 >>> ipix, depth, fully_covered = elliptical_cone_search(lon, lat, a, b, pa, max_depth)