Skip to content

Commit 020d715

Browse files
authored
sample positions uniformly over HEALPix pixel (#344)
1 parent 24bfce1 commit 020d715

File tree

6 files changed

+159
-2
lines changed

6 files changed

+159
-2
lines changed
Lines changed: 8 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,8 @@
1+
# healpix map information
2+
nside: 2
3+
4+
# sample 500 galaxies in healpix pixel 5
5+
positions: !skypy.position.uniform_in_pixel
6+
nside: $nside
7+
ipix: 5
8+
size: 500

docs/position/index.rst

Lines changed: 26 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -33,6 +33,32 @@ use the `~skypy.position.uniform_around` function:
3333
plt.grid()
3434

3535

36+
.. _skypy.position.uniform_in_pixel:
37+
38+
A more complicated uniform distribution is that over a given healpix pixel,
39+
which is sampled by the `~skypy.position.uniform_in_pixel` function:
40+
41+
.. literalinclude:: examples/uniform_in_pixel.yml
42+
:language: yaml
43+
44+
.. plot::
45+
:include-source: false
46+
47+
from skypy.pipeline import Pipeline
48+
49+
pipeline = Pipeline.read('examples/uniform_in_pixel.yml')
50+
pipeline.execute()
51+
52+
coords = pipeline['positions']
53+
ra, dec = coords.ra.wrap_at('180d').radian, coords.dec.radian
54+
55+
import matplotlib.pyplot as plt
56+
plt.figure(figsize=(8,4.2))
57+
plt.subplot(111, projection="aitoff")
58+
plt.plot(ra, dec, '.', alpha=0.2)
59+
plt.grid()
60+
61+
3662
Reference/API
3763
=============
3864

setup.cfg

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -36,10 +36,12 @@ test =
3636
all =
3737
camb
3838
classy @ git+https://github.com/skypyproject/[email protected]#egg=classy
39+
healpy
3940
skypy-data @ https://github.com/skypyproject/skypy-data/archive/master.tar.gz
4041
specutils
4142
docs =
4243
sphinx-astropy
44+
healpy
4345
matplotlib
4446

4547
[options.package_data]

skypy/position/__init__.py

Lines changed: 2 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2,8 +2,9 @@
22

33
__all__ = []
44

5-
from ._uniform import uniform_around
5+
from ._uniform import uniform_around, uniform_in_pixel
66

77
__all__ += [
88
'uniform_around',
9+
'uniform_in_pixel',
910
]

skypy/position/_uniform.py

Lines changed: 73 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -2,9 +2,9 @@
22

33
import numpy as np
44
from astropy import units
5+
from astropy.coordinates import SkyCoord
56

67
TWO_PI = 2*np.pi
7-
PI_HALF = np.pi/2
88

99

1010
@units.quantity_input(area=units.sr)
@@ -45,3 +45,75 @@ def uniform_around(centre, area, size):
4545

4646
# construct random sky coordinates around centre
4747
return centre.directional_offset_by(phi, theta)
48+
49+
50+
def uniform_in_pixel(nside, ipix, size, nest=False):
51+
'''Uniform distribution of points over healpix pixel.
52+
53+
Draws randomly distributed points from the healpix pixel `ipix` for a map
54+
with a given `nside` parameter.
55+
56+
Parameters
57+
----------
58+
nside : int
59+
Healpix map `nside` parameter.
60+
ipix : int
61+
Healpix map pixel index.
62+
size : int
63+
Number of points to draw.
64+
nest : bool, optional
65+
If True assume ``NESTED`` pixel ordering, otherwise ``RING`` pixel
66+
ordering. Default is ``RING`` pixel ordering.
67+
68+
Returns
69+
-------
70+
coords : `~astropy.coordinates.SkyCoord`
71+
Randomly distributed points over the healpix pixel.
72+
73+
Warnings
74+
--------
75+
This function requires the ``healpy`` package.
76+
77+
Examples
78+
--------
79+
See :ref:`User Documentation <skypy.position.uniform_in_pixel>`.
80+
81+
'''
82+
83+
from healpy import pix2ang, max_pixrad, nside2pixarea, ang2pix
84+
85+
# get the centre of the healpix pixel as a SkyCoord
86+
centre_lon, centre_lat = pix2ang(nside, ipix, nest=nest, lonlat=True)
87+
centre = SkyCoord(centre_lon, centre_lat, unit=units.deg)
88+
89+
# get the maximum radius of a healpix pixel in radian
90+
r = max_pixrad(nside)
91+
92+
# use that radius as the aperture of a spherical area in steradian
93+
area = TWO_PI*(1 - np.cos(r))*units.sr
94+
95+
# oversampling factor = 1/(probability of the draw)
96+
over = area.value/nside2pixarea(nside)
97+
98+
# the array of longitudes and latitudes of the sample
99+
lon, lat = np.empty(0), np.empty(0)
100+
101+
# rejection sampling over irregularly shaped healpix pixels
102+
miss = size
103+
while miss > 0:
104+
# get the coordinates in a circular aperture around centre
105+
sample = uniform_around(centre, area, int(np.ceil(miss*over)))
106+
107+
# get longitude and latitude of the sample
108+
sample_lon, sample_lat = sample.ra.deg, sample.dec.deg
109+
110+
# accept those positions that are inside the correct pixel
111+
accept = ipix == ang2pix(nside, sample_lon, sample_lat, nest=nest, lonlat=True)
112+
113+
# store the new positions
114+
lon = np.append(lon, np.extract(accept, sample_lon))
115+
lat = np.append(lat, np.extract(accept, sample_lat))
116+
miss = size - len(lon)
117+
118+
# construct the coordinates
119+
return SkyCoord(lon[:size], lat[:size], unit=units.deg)

skypy/position/tests/test_healpy.py

Lines changed: 48 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,48 @@
1+
import numpy as np
2+
import pytest
3+
4+
healpy = pytest.importorskip('healpy')
5+
6+
7+
def test_uniform_in_pixel_correctness():
8+
from skypy.position import uniform_in_pixel
9+
10+
# 768 healpix pixels
11+
nside = 8
12+
npix = healpy.nside2npix(nside)
13+
14+
ipix = np.random.randint(npix)
15+
16+
size = 1000
17+
pos = uniform_in_pixel(nside, ipix, size=size)
18+
assert len(pos) == size
19+
20+
lon, lat = pos.ra.deg, pos.dec.deg
21+
apix = healpy.ang2pix(nside, lon, lat, lonlat=True)
22+
np.testing.assert_array_equal(apix, ipix)
23+
24+
25+
@pytest.mark.flaky
26+
def test_uniform_in_pixel_distribution():
27+
from scipy.stats import kstest
28+
29+
from skypy.position import uniform_in_pixel
30+
31+
# 48 healpix pixels
32+
nside = 2
33+
npix = healpy.nside2npix(nside)
34+
35+
# sample entire sky with 20 positions per pixel
36+
size = 20
37+
theta, phi = np.empty(npix*size), np.empty(npix*size)
38+
for ipix in range(npix):
39+
pos = uniform_in_pixel(nside, ipix, size=size)
40+
assert len(pos) == size
41+
theta[ipix*size:(ipix+1)*size] = np.pi/2 - pos.dec.rad
42+
phi[ipix*size:(ipix+1)*size] = pos.ra.rad
43+
44+
# test for uniform distribution
45+
D, p = kstest(theta, lambda t: (1 - np.cos(t))/2)
46+
assert p > 0.01, 'D = {}, p = {}'.format(D, p)
47+
D, p = kstest(phi, 'uniform', args=(0, 2*np.pi))
48+
assert p > 0.01, 'D = {}, p = {}'.format(D, p)

0 commit comments

Comments
 (0)