File: test_pixelfunc.py

package info (click to toggle)
healpy 1.19.0-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 17,464 kB
  • sloc: ansic: 113,657; cpp: 15,827; python: 10,793; sh: 8,443; yacc: 5,410; fortran: 2,613; lex: 553; makefile: 380
file content (289 lines) | stat: -rw-r--r-- 11,665 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
from healpy.pixelfunc import *
from healpy._query_disc import boundaries
from healpy._pixelfunc import pix2ring, isnsideok
import numpy as np
import unittest
from healpy._query_disc import query_strip
from healpy.pixelfunc import ring2nest


class TestPixelFunc(unittest.TestCase):
    def setUp(self):
        # data fixture
        self.theta0 = [1.52911759, 0.78550497, 1.57079633, 0.05103658, 3.09055608]
        self.phi0 = [0.0, 0.78539816, 1.61988371, 0.78539816, 0.78539816]
        self.lon0 = np.degrees(self.phi0)
        self.lat0 = 90.0 - np.degrees(self.theta0)

    def test_nside2npix(self):
        self.assertEqual(nside2npix(512), 3145728)
        self.assertEqual(nside2npix(1024), 12582912)

    def test_nside2resol(self):
        self.assertAlmostEqual(nside2resol(512, arcmin=True), 6.87097282363)
        self.assertAlmostEqual(nside2resol(1024, arcmin=True), 3.43548641181)

    def test_max_pixrad(self):
        self.assertAlmostEqual(max_pixrad(512), 2.0870552355e-03)
        self.assertAlmostEqual(
            max_pixrad(512, degrees=True), np.rad2deg(2.0870552355e-03)
        )

    def test_nside2pixarea(self):
        self.assertAlmostEqual(nside2pixarea(512), 3.9947416351188569e-06)

    def test_ang2pix_ring(self):
        # ensure nside = 1 << 23 is correctly calculated
        # by comparing the original theta phi are restored.
        # NOTE: nside needs to be sufficiently large!
        id = ang2pix(1048576 * 8, self.theta0, self.phi0, nest=False)
        theta1, phi1 = pix2ang(1048576 * 8, id, nest=False)
        np.testing.assert_array_almost_equal(theta1, self.theta0)
        np.testing.assert_array_almost_equal(phi1, self.phi0)

    def test_ang2pix_ring_outofrange(self):
        # Healpy_Base2 works up to nside = 2**29.
        # Check that a ValueError is raised for nside = 2**30.
        self.assertRaises(
            ValueError, ang2pix, 1 << 30, self.theta0, self.phi0, nest=False
        )

    def test_ang2pix_nest(self):
        # ensure nside = 1 << 23 is correctly calculated
        # by comparing the original theta phi are restored.
        # NOTE: nside needs to be sufficiently large!
        # NOTE: with Healpy_Base this will fail because nside
        #       is limited to 1 << 13 with Healpy_Base.
        id = ang2pix(1048576 * 8, self.theta0, self.phi0, nest=True)
        theta1, phi1 = pix2ang(1048576 * 8, id, nest=True)
        np.testing.assert_array_almost_equal(theta1, self.theta0)
        np.testing.assert_array_almost_equal(phi1, self.phi0)

        self.assertTrue(np.allclose(theta1, self.theta0))
        self.assertTrue(np.allclose(phi1, self.phi0))

    def test_ang2pix_nest_outofrange_doesntcrash(self):
        # Healpy_Base2 works up to nside = 2**29.
        # Check that a ValueError is raised for nside = 2**30.
        self.assertRaises(
            ValueError, ang2pix, 1 << 30, self.theta0, self.phi0, nest=False
        )

    def test_ang2pix_negative_theta(self):
        self.assertRaises(ValueError, ang2pix, 32, -1, 0)

    def test_ang2pix_lonlat(self):
        # Need to decrease the precision of the check because deg not radians
        id = ang2pix(1048576 * 8, self.lon0, self.lat0, nest=False, lonlat=True)
        lon1, lat1 = pix2ang(1048576 * 8, id, nest=False, lonlat=True)
        np.testing.assert_array_almost_equal(lon1, self.lon0, decimal=5)
        np.testing.assert_array_almost_equal(lat1, self.lat0, decimal=5)

        # Now test nested
        id = ang2pix(1048576 * 8, self.theta0, self.phi0, nest=True)
        theta1, phi1 = pix2ang(1048576 * 8, id, nest=True)
        np.testing.assert_array_almost_equal(theta1, self.theta0)
        np.testing.assert_array_almost_equal(phi1, self.phi0)

    def test_vec2pix_lonlat(self):
        # Need to decrease the precision of the check because deg not radians
        vec = ang2vec(self.lon0, self.lat0, lonlat=True)
        lon1, lat1 = vec2ang(vec, lonlat=True)
        np.testing.assert_array_almost_equal(lon1, self.lon0, decimal=5)
        np.testing.assert_array_almost_equal(lat1, self.lat0, decimal=5)

    def test_get_interp_val_lonlat(self):
        m = np.arange(12.0)
        val0 = get_interp_val(m, self.theta0, self.phi0)
        val1 = get_interp_val(m, self.lon0, self.lat0, lonlat=True)
        np.testing.assert_array_almost_equal(val0, val1)

    def test_get_interp_weights(self):
        p0, w0 = (np.array([0, 1, 4, 5]), np.array([1.0, 0.0, 0.0, 0.0]))

        # phi not specified, theta assumed to be pixel
        p1, w1 = get_interp_weights(1, 0)
        np.testing.assert_array_almost_equal(p0, p1)
        np.testing.assert_array_almost_equal(w0, w1)

        # If phi is not specified, lonlat should do nothing
        p1, w1 = get_interp_weights(1, 0, lonlat=True)
        np.testing.assert_array_almost_equal(p0, p1)
        np.testing.assert_array_almost_equal(w0, w1)

        p0, w0 = (np.array([1, 2, 3, 0]), np.array([0.25, 0.25, 0.25, 0.25]))

        p1, w1 = get_interp_weights(1, 0, 0)
        np.testing.assert_array_almost_equal(p0, p1)
        np.testing.assert_array_almost_equal(w0, w1)

        p1, w1 = get_interp_weights(1, 0, 90, lonlat=True)
        np.testing.assert_array_almost_equal(p0, p1)
        np.testing.assert_array_almost_equal(w0, w1)

    def test_get_all_neighbours(self):
        ipix0 = np.array([8, 4, 0, -1, 1, 6, 9, -1])
        ipix1 = get_all_neighbours(1, np.pi / 2, np.pi / 2)
        ipix2 = get_all_neighbours(1, 90, 0, lonlat=True)
        np.testing.assert_array_almost_equal(ipix0, ipix1)
        np.testing.assert_array_almost_equal(ipix0, ipix2)

    def test_fit_dipole(self):
        nside = 32
        npix = nside2npix(nside)
        d = [0.3, 0.5, 0.2]
        vec = np.transpose(pix2vec(nside, np.arange(npix)))
        signal = np.dot(vec, d)
        mono, dipole = fit_dipole(signal)
        self.assertAlmostEqual(mono, 0.0)
        self.assertAlmostEqual(d[0], dipole[0])
        self.assertAlmostEqual(d[1], dipole[1])
        self.assertAlmostEqual(d[2], dipole[2])

    def test_boundaries(self):
        """Test whether the boundary shapes look sane"""
        for lgNside in range(1, 5):
            nside = 1 << lgNside
            for pix in range(nside2npix(nside)):
                for res in range(1, 50, 7):
                    num = 4 * res  # Expected number of points
                    for nest in (True, False):
                        points = boundaries(nside, pix, res, nest=nest)
                        self.assertTrue(points.shape == (3, num))
                        dist = np.linalg.norm(
                            points[:, : num - 1] - points[:, 1:]
                        )  # distance between points
                        self.assertTrue((dist != 0).all())
                        dmin = np.min(dist)
                        dmax = np.max(dist)
                        self.assertTrue(dmax / dmin <= 2.0)

    def test_ring(self):
        for lgNside in range(1, 5):
            nside = 1 << lgNside
            numPix = nside2npix(nside)
            numRings = 4 * nside - 1  # Expected number of rings
            for nest in (True, False):
                pix = np.arange(numPix, dtype=np.int64)
                ring = pix2ring(nside, pix, nest=nest)
                self.assertTrue(pix.shape == ring.shape)
                self.assertTrue(len(set(ring)) == numRings)
                if not nest:
                    first = ring[: numPix - 1]
                    second = ring[1:]
                    self.assertTrue(
                        np.logical_or(first == second, first == second - 1).all()
                    )

    def test_accept_ma_allows_only_keywords(self):
        """Test whether the accept_ma wrapper accepts calls using only keywords."""
        ma = np.zeros(12 * 16 ** 2)
        try:
            ud_grade(map_in=ma, nside_out=32)
        except IndexError:
            self.fail("IndexError raised")

    def test_isnsideok(self):
        """Test the isnsideok."""
        self.assertTrue(isnsideok(nside=1, nest=False))
        self.assertTrue(isnsideok(nside=16, nest=True))

        self.assertTrue(not isnsideok(nside=-16, nest=True))
        self.assertTrue(not isnsideok(nside=-16, nest=False))
        self.assertTrue(not isnsideok(nside=13, nest=True))

    def test_latauto_and_latbounce(self):
        """Test latauto and latbounce features in lonlat2thetaphi."""

        # Test array input
        lat = [-100, -90, 0, 90, 100]
        lon = [0, 60, 90, 0, 180]

        # Test with latauto=False (default)
        theta, phi = lonlat2thetaphi(lon, lat)
        self.assertTrue(
            np.allclose(
                theta,
                [
                    np.pi / 2 - np.deg2rad(-100),
                    np.pi,
                    np.pi / 2,
                    0,
                    np.pi / 2 - np.deg2rad(100),
                ],
            )
        )
        self.assertTrue(
            np.allclose(phi, np.deg2rad([lon[0], lon[1], lon[2], lon[3], lon[4]]))
        )

        # Test with latauto=True and latbounce=True (default)
        theta, phi = lonlat2thetaphi(lon, lat, latauto=True)
        self.assertTrue(
            np.allclose(
                theta,
                [
                    np.pi / 2 + np.deg2rad(80),
                    np.pi,
                    np.pi / 2,
                    0,
                    np.pi / 2 - np.deg2rad(80),
                ],
            )
        )
        self.assertTrue(
            np.allclose(phi, np.deg2rad([lon[0], lon[1], lon[2], lon[3], lon[4]]))
        )

        # Test with latauto=True and latbounce=False
        theta, phi = lonlat2thetaphi(lon, lat, latauto=True, latbounce=False)
        self.assertTrue(
            np.allclose(
                theta,
                [
                    np.pi / 2 + np.deg2rad(80),
                    np.pi / 2 + np.pi / 2,
                    np.pi / 2,
                    0,
                    np.pi / 2 - np.deg2rad(80),
                ],
            )
        )
        self.assertTrue(np.allclose(phi, np.deg2rad([180, 60, 90, 0, 180 + 180])))

        # Test scalar input
        lat_scalar = -100
        lon_scalar = 10

        # Test with latauto=False (default)
        theta_scalar, phi_scalar = lonlat2thetaphi(lon_scalar, lat_scalar)
        self.assertTrue(np.allclose(theta_scalar, np.pi / 2 + np.deg2rad(100)))
        self.assertTrue(np.allclose(phi_scalar, np.deg2rad(lon_scalar)))

        # Test with latauto=True and latbounce=True (default)
        theta_scalar, phi_scalar = lonlat2thetaphi(lon_scalar, lat_scalar, latauto=True)
        self.assertTrue(np.allclose(theta_scalar, np.pi / 2 + np.deg2rad(80)))
        self.assertTrue(np.allclose(phi_scalar, np.radians(lon_scalar)))

        # Test with latauto=True and latbounce=False
        theta_scalar, phi_scalar = lonlat2thetaphi(
            lon_scalar, lat_scalar, latauto=True, latbounce=False
        )
        self.assertTrue(np.allclose(theta_scalar, np.pi / 2 + np.deg2rad(80)))
        self.assertTrue(np.allclose(phi_scalar, np.deg2rad(lon_scalar) + np.pi))

    def test_query_strip_nest(self):
        # Test query_strip with nest=True, which was previously crashing
        nside = 2
        theta1 = 1
        theta2 = 2

        # Expected result using the workaround from the issue description
        expected_result = ring2nest(
            nside, query_strip(nside, theta1, theta2, nest=False)
        )

        # Actual result with nest=True
        actual_result = query_strip(nside, theta1, theta2, nest=True)

        np.testing.assert_array_equal(actual_result, expected_result)