/
hpx.py
1896 lines (1564 loc) · 60.1 KB
/
hpx.py
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
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# Licensed under a 3-clause BSD style license - see LICENSE.rst
"""Utilities for dealing with HEALPix projections and mappings."""
from collections import OrderedDict
import re
import copy
import numpy as np
from astropy.io import fits
from astropy.coordinates import SkyCoord
from astropy.units import Quantity
from .utils import INVALID_INDEX
from .wcs import WcsGeom
from .geom import MapGeom, MapCoord, pix_tuple_to_idx
from .geom import coordsys_to_frame, skycoord_to_lonlat
from .geom import find_and_read_bands, make_axes
# Not sure if we should expose this in the docs or not:
# HPX_FITS_CONVENTIONS, HpxConv
__all__ = ["HpxGeom"]
# Approximation of the size of HEALPIX pixels (in degrees) for a particular order.
# Used to convert from HEALPIX to WCS-based projections.
HPX_ORDER_TO_PIXSIZE = np.array(
[32.0, 16.0, 8.0, 4.0, 2.0, 1.0, 0.50, 0.25, 0.1, 0.05, 0.025, 0.01, 0.005, 0.002]
)
class HpxConv:
"""Data structure to define how a HEALPIX map is stored to FITS."""
def __init__(self, convname, **kwargs):
self.convname = convname
self.colstring = kwargs.get("colstring", "CHANNEL")
self.firstcol = kwargs.get("firstcol", 1)
self.hduname = kwargs.get("hduname", "SKYMAP")
self.bands_hdu = kwargs.get("bands_hdu", "EBOUNDS")
self.quantity_type = kwargs.get("quantity_type", "integral")
self.coordsys = kwargs.get("coordsys", "COORDSYS")
def colname(self, indx):
return "{}{}".format(self.colstring, indx)
@classmethod
def create(cls, convname="gadf"):
return copy.deepcopy(HPX_FITS_CONVENTIONS[convname])
HPX_FITS_CONVENTIONS = {}
"""Various conventions for storing HEALPIX maps in FITS files"""
HPX_FITS_CONVENTIONS[None] = HpxConv("gadf", bands_hdu="BANDS")
HPX_FITS_CONVENTIONS["gadf"] = HpxConv("gadf", bands_hdu="BANDS")
HPX_FITS_CONVENTIONS["fgst-ccube"] = HpxConv("fgst-ccube")
HPX_FITS_CONVENTIONS["fgst-ltcube"] = HpxConv(
"fgst-ltcube", colstring="COSBINS", hduname="EXPOSURE", bands_hdu="CTHETABOUNDS"
)
HPX_FITS_CONVENTIONS["fgst-bexpcube"] = HpxConv(
"fgst-bexpcube", colstring="ENERGY", hduname="HPXEXPOSURES", bands_hdu="ENERGIES"
)
HPX_FITS_CONVENTIONS["fgst-srcmap"] = HpxConv(
"fgst-srcmap", hduname=None, quantity_type="differential"
)
HPX_FITS_CONVENTIONS["fgst-template"] = HpxConv(
"fgst-template", colstring="ENERGY", bands_hdu="ENERGIES"
)
HPX_FITS_CONVENTIONS["fgst-srcmap-sparse"] = HpxConv(
"fgst-srcmap-sparse", colstring=None, hduname=None, quantity_type="differential"
)
HPX_FITS_CONVENTIONS["galprop"] = HpxConv(
"galprop",
colstring="Bin",
hduname="SKYMAP2",
bands_hdu="ENERGIES",
quantity_type="differential",
coordsys="COORDTYPE",
)
HPX_FITS_CONVENTIONS["galprop2"] = HpxConv(
"galprop",
colstring="Bin",
hduname="SKYMAP2",
bands_hdu="ENERGIES",
quantity_type="differential",
)
def unravel_hpx_index(idx, npix):
"""Convert flattened global map index to an index tuple.
Parameters
----------
idx : `~numpy.ndarray`
Flat index.
npix : `~numpy.ndarray`
Number of pixels in each band.
Returns
-------
idx : tuple of `~numpy.ndarray`
Index array for each dimension of the map.
"""
if npix.size == 1:
return tuple([idx])
dpix = np.zeros(npix.size, dtype="i")
dpix[1:] = np.cumsum(npix.flat[:-1])
bidx = np.searchsorted(np.cumsum(npix.flat), idx + 1)
pix = idx - dpix[bidx]
return tuple([pix] + list(np.unravel_index(bidx, npix.shape)))
def ravel_hpx_index(idx, npix):
"""Convert map index tuple to a flattened index.
Parameters
----------
idx : tuple of `~numpy.ndarray`
Returns
-------
idx : `~numpy.ndarray`
"""
if len(idx) == 1:
return idx[0]
# TODO: raise exception for indices that are out of bounds
idx0 = idx[0]
idx1 = np.ravel_multi_index(idx[1:], npix.shape, mode="clip")
npix = np.concatenate((np.array([0]), npix.flat[:-1]))
return idx0 + np.cumsum(npix)[idx1]
def coords_to_vec(lon, lat):
"""Converts longitude and latitude coordinates to a unit 3-vector.
Returns
-------
array(3,n) with v_x[i],v_y[i],v_z[i] = directional cosines
"""
phi = np.radians(lon)
theta = (np.pi / 2) - np.radians(lat)
sin_t = np.sin(theta)
cos_t = np.cos(theta)
x = sin_t * np.cos(phi)
y = sin_t * np.sin(phi)
z = cos_t
# Stack them into the output array
out = np.vstack((x, y, z)).swapaxes(0, 1)
return out
def get_nside_from_pix_size(pixsz):
"""Get the NSIDE that is closest to the given pixel size.
Parameters
----------
pix : `~numpy.ndarray`
Pixel size in degrees.
Returns
-------
nside : `~numpy.ndarray`
NSIDE parameter.
"""
import healpy as hp
pixsz = np.array(pixsz, ndmin=1)
nside = 2 ** np.linspace(1, 14, 14, dtype=int)
nside_pixsz = np.degrees(hp.nside2resol(nside))
return nside[np.argmin(np.abs(nside_pixsz - pixsz[..., None]), axis=-1)]
def get_pix_size_from_nside(nside):
"""Estimate of the pixel size from the HEALPIX nside coordinate.
This just uses a lookup table to provide a nice round number
for each HEALPIX order.
"""
order = nside_to_order(nside)
if np.any(order < 0) or np.any(order > 13):
raise ValueError("HEALPIX order must be 0 to 13. Got: {!r}".format(order))
return HPX_ORDER_TO_PIXSIZE[order]
def hpx_to_axes(h, npix):
"""Generate a sequence of bin edge vectors corresponding to the axes of a HPX object.
"""
x = h.ebins
z = np.arange(npix[-1] + 1)
return x, z
def hpx_to_coords(h, shape):
"""Generate an N x D list of pixel center coordinates.
``N`` is the number of pixels and ``D`` is the dimensionality of the map.
"""
x, z = hpx_to_axes(h, shape)
x = np.sqrt(x[0:-1] * x[1:])
z = z[:-1] + 0.5
x = np.ravel(np.ones(shape) * x[:, np.newaxis])
z = np.ravel(np.ones(shape) * z[np.newaxis, :])
return np.vstack((x, z))
def make_hpx_to_wcs_mapping_centers(hpx, wcs):
"""Make the mapping data needed to from from HPX pixelization to a WCS-based array.
Parameters
----------
hpx : `~gammapy.maps.HpxGeom`
The HEALPIX geometry.
wcs : `~gammapy.maps.WcsGeom`
The WCS geometry.
Returns
-------
ipixs : array(nx,ny)
HEALPIX pixel indices for each WCS pixel
-1 indicates the wcs pixel does not contain the center of a HEALPIX pixel
mult_val : array
(nx,ny) of 1.
npix : tuple
(nx,ny) with the shape of the WCS grid
"""
npix = (int(wcs.wcs.crpix[0] * 2), int(wcs.wcs.crpix[1] * 2))
mult_val = np.ones(npix).T.flatten()
sky_crds = hpx.get_sky_coords()
pix_crds = wcs.wcs_world2pix(sky_crds, 0).astype(int)
ipixs = -1 * np.ones(npix, int).T.flatten()
pix_index = npix[1] * pix_crds[0:, 0] + pix_crds[0:, 1]
if hpx._ipix is None:
for ipix, pix_crd in enumerate(pix_index):
ipixs[pix_crd] = ipix
else:
for pix_crd, ipix in zip(pix_index, hpx._ipix):
ipixs[pix_crd] = ipix
ipixs = ipixs.reshape(npix).T.flatten()
return ipixs, mult_val, npix
def make_hpx_to_wcs_mapping(hpx, wcs):
"""Make the pixel mapping from HPX- to a WCS-based geometry.
Parameters
----------
hpx : `~gammapy.maps.HpxGeom`
The HEALPIX geometry
wcs : `~gammapy.maps.WcsGeom`
The WCS geometry
Returns
-------
ipix : `~numpy.ndarray`
array(nx,ny) of HEALPIX pixel indices for each wcs pixel
mult_val : `~numpy.ndarray`
array(nx,ny) of 1./number of WCS pixels pointing at each HEALPIX pixel
npix : tuple
tuple(nx,ny) with the shape of the WCS grid
"""
import healpy as hp
npix = wcs.npix
# FIXME: Calculation of WCS pixel centers should be moved into a
# method of WcsGeom
pix_crds = np.dstack(np.meshgrid(np.arange(npix[0]), np.arange(npix[1])))
pix_crds = pix_crds.swapaxes(0, 1).reshape((-1, 2))
sky_crds = wcs.wcs.wcs_pix2world(pix_crds, 0)
sky_crds *= np.radians(1.0)
sky_crds[0:, 1] = (np.pi / 2) - sky_crds[0:, 1]
mask = ~np.any(np.isnan(sky_crds), axis=1)
ipix = -1 * np.ones((len(hpx.nside), int(npix[0] * npix[1])), int)
m = mask[None, :] * np.ones_like(ipix, dtype=bool)
ipix[m] = hp.ang2pix(
hpx.nside[..., None],
sky_crds[:, 1][mask][None, ...],
sky_crds[:, 0][mask][None, ...],
hpx.nest,
).flatten()
# Here we are counting the number of HEALPIX pixels each WCS pixel
# points to and getting a multiplicative factor that tells use how
# to split up the counts in each HEALPIX pixel (by dividing the
# corresponding WCS pixels by the number of associated HEALPIX
# pixels).
mult_val = np.ones_like(ipix, dtype=float)
for i, t in enumerate(ipix):
count = np.unique(t, return_counts=True)
idx = np.searchsorted(count[0], t)
mult_val[i, ...] = 1.0 / count[1][idx]
if hpx.nside.size == 1:
ipix = np.squeeze(ipix, axis=0)
mult_val = np.squeeze(mult_val, axis=0)
return ipix, mult_val, npix
def match_hpx_pix(nside, nest, nside_pix, ipix_ring):
"""TODO
"""
import healpy as hp
ipix_in = np.arange(12 * nside * nside)
vecs = hp.pix2vec(nside, ipix_in, nest)
pix_match = hp.vec2pix(nside_pix, vecs[0], vecs[1], vecs[2]) == ipix_ring
return ipix_in[pix_match]
def parse_hpxregion(region):
"""Parse the ``HPX_REG`` header keyword into a list of tokens.
"""
m = re.match(r"([A-Za-z\_]*?)\((.*?)\)", region)
if m is None:
raise ValueError("Failed to parse hpx region string: {!r}".format(region))
if not m.group(1):
return re.split(",", m.group(2))
else:
return [m.group(1)] + re.split(",", m.group(2))
def get_hpxregion_dir(region, coordsys):
"""Get the reference direction for a HEALPIX region string.
Parameters
----------
region : str
A string describing a HEALPIX region
coordsys : {'CEL', 'GAL'}
Coordinate system
"""
import healpy as hp
frame = coordsys_to_frame(coordsys)
if region is None:
return SkyCoord(0.0, 0.0, frame=frame, unit="deg")
tokens = parse_hpxregion(region)
if tokens[0] in ["DISK", "DISK_INC"]:
lon, lat = float(tokens[1]), float(tokens[2])
return SkyCoord(lon, lat, frame=frame, unit="deg")
elif tokens[0] == "HPX_PIXEL":
nside_pix = int(tokens[2])
ipix_pix = int(tokens[3])
if tokens[1] == "NESTED":
nest_pix = True
elif tokens[1] == "RING":
nest_pix = False
else:
raise ValueError("Invalid ordering scheme: {!r}".format(tokens[1]))
theta, phi = hp.pix2ang(nside_pix, ipix_pix, nest_pix)
lat = np.degrees((np.pi / 2) - theta)
lon = np.degrees(phi)
return SkyCoord(lon, lat, frame=frame, unit="deg")
else:
raise ValueError("Invalid region type: {!r}".format(tokens[0]))
def get_hpxregion_size(region):
"""Get the approximate size of region (in degrees) from a HEALPIX region string.
"""
tokens = parse_hpxregion(region)
if tokens[0] in {"DISK", "DISK_INC"}:
return float(tokens[3])
elif tokens[0] == "HPX_PIXEL":
pix_size = get_pix_size_from_nside(int(tokens[2]))
return 2.0 * pix_size
else:
raise ValueError("Invalid region type: {!r}".format(tokens[0]))
def is_power2(n):
"""Check if an integer is a power of 2."""
return (n > 0) & ((n & (n - 1)) == 0)
def nside_to_order(nside):
"""Compute the ORDER given NSIDE.
Returns -1 for NSIDE values that are not a power of 2.
"""
nside = np.array(nside, ndmin=1)
order = -1 * np.ones_like(nside)
m = is_power2(nside)
order[m] = np.log2(nside[m]).astype(int)
return order
def upix_to_pix(upix):
"""Get the pixel index and nside from a unique pixel number."""
nside = np.power(2, np.floor(np.log2(upix / 4)) / 2).astype(int)
pix = upix - 4 * np.power(nside, 2)
return pix, nside
def pix_to_upix(pix, nside):
"""Compute the unique pixel number from the pixel number and nside."""
return pix + 4 * np.power(nside, 2)
def get_superpixels(idx, nside_subpix, nside_superpix, nest=True):
"""Compute the indices of superpixels that contain a subpixel.
Parameters
----------
idx : `~numpy.ndarray`
Array of HEALPix pixel indices for subpixels of NSIDE
``nside_subpix``.
nside_subpix : int or `~numpy.ndarray`
NSIDE of subpixel.
nside_superpix : int or `~numpy.ndarray`
NSIDE of superpixel.
nest : bool
If True, assume NESTED pixel ordering, otherwise, RING pixel
ordering.
Returns
-------
idx_super : `~numpy.ndarray`
Indices of HEALpix pixels of nside ``nside_superpix`` that
contain pixel indices ``idx`` of nside ``nside_subpix``.
"""
import healpy as hp
idx = np.array(idx)
nside_superpix = np.asarray(nside_superpix)
nside_subpix = np.asarray(nside_subpix)
if not nest:
idx = hp.ring2nest(nside_subpix, idx)
if np.any(~is_power2(nside_superpix)) or np.any(~is_power2(nside_subpix)):
raise ValueError("NSIDE must be a power of 2.")
ratio = np.array((nside_subpix // nside_superpix) ** 2, ndmin=1)
idx //= ratio
if not nest:
m = idx == INVALID_INDEX.int
idx[m] = 0
idx = hp.nest2ring(nside_superpix, idx)
idx[m] = INVALID_INDEX.int
return idx
def get_subpixels(idx, nside_superpix, nside_subpix, nest=True):
"""Compute the indices of subpixels contained within superpixels.
This function returns an output array with one additional
dimension of size N for subpixel indices where N is the maximum
number of subpixels for any pair of ``nside_superpix`` and
``nside_subpix``. If the number of subpixels is less than N the
remaining subpixel indices will be set to -1.
Parameters
----------
idx : `~numpy.ndarray`
Array of HEALPix pixel indices for superpixels of NSIDE
``nside_superpix``.
nside_superpix : int or `~numpy.ndarray`
NSIDE of superpixel.
nside_subpix : int or `~numpy.ndarray`
NSIDE of subpixel.
nest : bool
If True, assume NESTED pixel ordering, otherwise, RING pixel
ordering.
Returns
-------
idx_sub : `~numpy.ndarray`
Indices of HEALpix pixels of nside ``nside_subpix`` contained
within pixel indices ``idx`` of nside ``nside_superpix``.
"""
import healpy as hp
if not nest:
idx = hp.ring2nest(nside_superpix, idx)
idx = np.asarray(idx)
nside_superpix = np.asarray(nside_superpix)
nside_subpix = np.asarray(nside_subpix)
if np.any(~is_power2(nside_superpix)) or np.any(~is_power2(nside_subpix)):
raise ValueError("NSIDE must be a power of 2.")
# number of subpixels in each superpixel
npix = np.array((nside_subpix // nside_superpix) ** 2, ndmin=1)
x = np.arange(np.max(npix), dtype=int)
idx = idx * npix
if not np.all(npix[0] == npix):
x = np.broadcast_to(x, idx.shape + x.shape)
idx = idx[..., None] + x
idx[x >= np.broadcast_to(npix[..., None], x.shape)] = INVALID_INDEX.int
else:
idx = idx[..., None] + x
if not nest:
m = idx == INVALID_INDEX.int
idx[m] = 0
idx = hp.nest2ring(nside_subpix[..., None], idx)
idx[m] = INVALID_INDEX.int
return idx
class HpxGeom(MapGeom):
"""Geometry class for HEALPIX maps.
This class performs mapping between partial-sky indices (pixel
number within a HEALPIX region) and all-sky indices (pixel number
within an all-sky HEALPIX map). Multi-band HEALPIX geometries use
a global indexing scheme that assigns a unique pixel number based
on the all-sky index and band index. In the single-band case the
global index is the same as the HEALPIX index.
By default the constructor will return an all-sky map.
Partial-sky maps can be defined with the ``region`` argument.
Parameters
----------
nside : `~numpy.ndarray`
HEALPIX nside parameter, the total number of pixels is
12*nside*nside. For multi-dimensional maps one can pass
either a single nside value or a vector of nside values
defining the pixel size for each image plane. If nside is not
a scalar then its dimensionality should match that of the
non-spatial axes.
nest : bool
True -> 'NESTED', False -> 'RING' indexing scheme
coordsys : str
Coordinate system, 'CEL' | 'GAL'
region : str or tuple
Spatial geometry for partial-sky maps. If none the map will
encompass the whole sky. String input will be parsed
according to HPX_REG header keyword conventions. Tuple
input can be used to define an explicit list of pixels
encompassed by the geometry.
axes : list
Axes for non-spatial dimensions.
conv : str
Convention for FITS serialization format.
sparse : bool
If True defer allocation of partial- to all-sky index mapping
arrays. This option is only compatible with partial-sky maps
with an analytic geometry (e.g. DISK).
"""
is_hpx = True
def __init__(
self,
nside,
nest=True,
coordsys="CEL",
region=None,
axes=None,
conv="gadf",
sparse=False,
):
# FIXME: Figure out what to do when sparse=True
# FIXME: Require NSIDE to be power of two when nest=True
self._nside = np.array(nside, ndmin=1)
self._axes = make_axes(axes, conv)
if self.nside.size > 1 and self.nside.shape != self.shape_axes:
raise ValueError(
"Wrong dimensionality for nside. nside must "
"be a scalar or have a dimensionality consistent "
"with the axes argument."
)
self._order = nside_to_order(self._nside)
self._nest = nest
self._coordsys = coordsys
self._maxpix = 12 * self._nside * self._nside
self._maxpix = self._maxpix * np.ones(self.shape_axes, dtype=int)
self._sparse = sparse
self._ipix = None
self._rmap = None
self._region = region
self._create_lookup(region)
if self._ipix is not None:
self._rmap = {}
for i, ipix in enumerate(self._ipix.flat):
self._rmap[ipix] = i
self._npix = self._npix * np.ones(self.shape_axes, dtype=int)
self._conv = conv
self._center_skydir = self._get_ref_dir()
lon, lat, frame = skycoord_to_lonlat(self._center_skydir)
self._center_coord = tuple(
[lon, lat]
+ [ax.pix_to_coord((float(ax.nbin) - 1.0) / 2.0) for ax in self.axes]
)
self._center_pix = self.coord_to_pix(self._center_coord)
@property
def data_shape(self):
"""Shape of the Numpy data array matching this geometry."""
npix_shape = [np.max(self.npix)]
ax_shape = [ax.nbin for ax in self.axes]
return tuple(npix_shape + ax_shape)[::-1]
def _create_lookup(self, region):
"""Create local-to-global pixel lookup table."""
if isinstance(region, str):
ipix = [
self.get_index_list(nside, self._nest, region)
for nside in self._nside.flat
]
self._ibnd = np.concatenate(
[i * np.ones_like(p, dtype="int16") for i, p in enumerate(ipix)]
)
self._ipix = [
ravel_hpx_index((p, i * np.ones_like(p)), np.ravel(self._maxpix))
for i, p in enumerate(ipix)
]
self._region = region
self._indxschm = "EXPLICIT"
self._npix = np.array([len(t) for t in self._ipix])
if self.nside.ndim > 1:
self._npix = self._npix.reshape(self.nside.shape)
self._ipix = np.concatenate(self._ipix)
elif isinstance(region, tuple):
region = [np.asarray(t) for t in region]
m = np.any(np.stack([t >= 0 for t in region]), axis=0)
region = [t[m] for t in region]
self._ipix = ravel_hpx_index(region, self._maxpix)
self._ipix = np.unique(self._ipix)
region = unravel_hpx_index(self._ipix, self._maxpix)
self._region = "explicit"
self._indxschm = "EXPLICIT"
if len(region) == 1:
self._npix = np.array([len(region[0])])
else:
self._npix = np.zeros(self.shape_axes, dtype=int)
idx = np.ravel_multi_index(region[1:], self.shape_axes)
cnt = np.unique(idx, return_counts=True)
self._npix.flat[cnt[0]] = cnt[1]
elif region is None:
self._region = None
self._indxschm = "IMPLICIT"
self._npix = self._maxpix
else:
raise ValueError("Invalid region string: {!r}".format(region))
def local_to_global(self, idx_local):
"""Compute a local index (partial-sky) from a global (all-sky)
index.
Returns
-------
idx_global : tuple
A tuple of pixel index vectors with global HEALPIX pixel
indices.
"""
if self._ipix is None:
return idx_local
if self.nside.size > 1:
idx = ravel_hpx_index(idx_local, self._npix)
else:
idx_tmp = tuple(
[idx_local[0]] + [np.zeros(t.shape, dtype=int) for t in idx_local[1:]]
)
idx = ravel_hpx_index(idx_tmp, self._npix)
idx_global = unravel_hpx_index(self._ipix[idx], self._maxpix)
return idx_global[:1] + tuple(idx_local[1:])
def global_to_local(self, idx_global, ravel=False):
"""Compute global (all-sky) index from a local (partial-sky) index.
Parameters
----------
idx_global : tuple
A tuple of pixel indices with global HEALPix pixel indices.
ravel : bool
Return a raveled index.
Returns
-------
idx_local : tuple
A tuple of pixel indices with local HEALPIX pixel indices.
"""
if self.nside.size == 1:
idx = np.array(idx_global[0], ndmin=1)
else:
idx = ravel_hpx_index(idx_global, self._maxpix)
if self._rmap is not None:
retval = np.full(idx.size, -1, "i")
m = np.in1d(idx.flat, self._ipix)
retval[m] = np.searchsorted(self._ipix, idx.flat[m])
retval = retval.reshape(idx.shape)
else:
retval = idx
if self.nside.size == 1:
idx_local = tuple([retval] + list(idx_global[1:]))
else:
idx_local = unravel_hpx_index(retval, self._npix)
m = np.any(np.stack([t == INVALID_INDEX.int for t in idx_local]), axis=0)
for i, t in enumerate(idx_local):
idx_local[i][m] = INVALID_INDEX.int
if not ravel:
return idx_local
else:
return ravel_hpx_index(idx_local, self.npix)
def __getitem__(self, idx_global):
"""Implement global-to-local index lookup.
For all-sky maps it just returns the input array. For
partial-sky maps it returns the local indices corresponding to
the indices in the input array, and -1 for those pixels that
are outside the selected region. For multi-dimensional maps
with a different ``NSIDE`` in each band the global index is an
unrolled index for both HEALPIX pixel number and image slice.
Parameters
----------
idx_global : `~numpy.ndarray`
An array of global (all-sky) pixel indices. If this is a
tuple, list, or array of integers it will be interpreted
as a global (raveled) index. If this argument is a tuple
of lists or arrays it will be interpreted as a list of
unraveled index vectors.
Returns
-------
idx_local : `~numpy.ndarray`
An array of local HEALPIX pixel indices.
"""
# Convert to tuple representation
if (
isinstance(idx_global, int)
or (isinstance(idx_global, tuple) and isinstance(idx_global[0], int))
or isinstance(idx_global, np.ndarray)
):
idx_global = unravel_hpx_index(np.array(idx_global, ndmin=1), self._maxpix)
return self.global_to_local(idx_global, ravel=True)
def coord_to_pix(self, coords):
import healpy as hp
coords = MapCoord.create(coords, coordsys=self.coordsys)
theta, phi = coords.theta, coords.phi
c = self.coord_to_tuple(coords)
if self.axes:
bins = []
idxs = []
for i, ax in enumerate(self.axes):
bins += [ax.coord_to_pix(c[i + 2])]
idxs += [ax.coord_to_idx(c[i + 2])]
# FIXME: Figure out how to handle coordinates out of
# bounds of non-spatial dimensions
if self.nside.size > 1:
nside = self.nside[tuple(idxs)]
else:
nside = self.nside
m = ~np.isfinite(theta)
theta[m] = 0.0
phi[m] = 0.0
pix = hp.ang2pix(nside, theta, phi, nest=self.nest)
pix = tuple([pix] + bins)
if np.any(m):
for p in pix:
p[m] = INVALID_INDEX.int
else:
pix = (hp.ang2pix(self.nside, theta, phi, nest=self.nest),)
return pix
def pix_to_coord(self, pix):
import healpy as hp
if self.axes:
bins = []
vals = []
for i, ax in enumerate(self.axes):
bins += [pix[1 + i]]
vals += [ax.pix_to_coord(pix[1 + i])]
idxs = pix_tuple_to_idx(bins)
if self.nside.size > 1:
nside = self.nside[idxs]
else:
nside = self.nside
ipix = np.round(pix[0]).astype(int)
m = ipix == INVALID_INDEX.int
ipix[m] = 0
theta, phi = hp.pix2ang(nside, ipix, nest=self.nest)
coords = [np.degrees(phi), np.degrees(np.pi / 2.0 - theta)]
coords = tuple(coords + vals)
if np.any(m):
for c in coords:
c[m] = INVALID_INDEX.float
else:
ipix = np.round(pix[0]).astype(int)
theta, phi = hp.pix2ang(self.nside, ipix, nest=self.nest)
coords = (np.degrees(phi), np.degrees(np.pi / 2.0 - theta))
return coords
def pix_to_idx(self, pix, clip=False):
# FIXME: Look for better method to clip HPX indices
# TODO: copy idx to avoid modifying input pix?
# pix_tuple_to_idx seems to always make a copy!?
idx = pix_tuple_to_idx(pix)
idx_local = self.global_to_local(idx)
for i, _ in enumerate(idx):
if clip:
if i > 0:
np.clip(idx[i], 0, self.axes[i - 1].nbin - 1, out=idx[i])
else:
np.clip(idx[i], 0, None, out=idx[i])
else:
if i > 0:
mask = (idx[i] < 0) | (idx[i] >= self.axes[i - 1].nbin)
np.putmask(idx[i], mask, -1)
else:
mask = (idx_local[i] < 0) | (idx[i] < 0)
np.putmask(idx[i], mask, -1)
return tuple(idx)
def to_slice(self, slices, drop_axes=True):
if len(slices) == 0 and self.ndim == 2:
return copy.deepcopy(self)
if len(slices) != self.ndim - 2:
raise ValueError()
nside = np.ones(self.shape_axes, dtype=int) * self.nside
nside = np.squeeze(nside[slices])
axes = [ax.slice(s) for ax, s in zip(self.axes, slices)]
if drop_axes:
axes = [ax for ax in axes if ax.nbin > 1]
slice_dims = [0] + [i + 1 for i, ax in enumerate(axes) if ax.nbin > 1]
else:
slice_dims = np.arange(self.ndim)
if self.region == "explicit":
idx = self.get_idx()
slices = (slice(None),) + slices
idx = [p[slices[::-1]] for p in idx]
idx = [p[p != INVALID_INDEX.int] for p in idx]
if drop_axes:
idx = [idx[i] for i in range(len(idx)) if i in slice_dims]
region = tuple(idx)
else:
region = self.region
return self.__class__(nside, self.nest, self.coordsys, region, axes, self.conv)
@property
def axes(self):
"""List of non-spatial axes."""
return self._axes
@property
def shape_axes(self):
"""Shape of non-spatial axes."""
return tuple([ax.nbin for ax in self._axes])
@property
def ndim(self):
"""Number of dimensions (int)."""
return len(self._axes) + 2
@property
def ordering(self):
"""HEALPix ordering ('NESTED' or 'RING')."""
return "NESTED" if self.nest else "RING"
@property
def nside(self):
"""NSIDE in each band."""
return self._nside
@property
def order(self):
"""ORDER in each band (NSIDE = 2**ORDER). Set to -1 for bands with
NSIDE that is not a power of 2.
"""
return self._order
@property
def nest(self):
"""Is HEALPix order nested? (bool)."""
return self._nest
@property
def npix(self):
"""Number of pixels in each band.
For partial-sky geometries this can
be less than the number of pixels for the band NSIDE.
"""
return self._npix
@property
def conv(self):
"""Name of default FITS convention associated with this geometry."""
return self._conv
@property
def hpx_conv(self):
return HPX_FITS_CONVENTIONS[self.conv]
@property
def coordsys(self):
return self._coordsys
@property
def projection(self):
"""Map projection."""
return "HPX"
@property
def region(self):
"""Region string."""
return self._region
@property
def is_allsky(self):
"""Flag for all-sky maps."""
if self._region is None:
return True
else:
return False
@property
def is_regular(self):
"""Flag identifying whether this geometry is regular in non-spatial dimensions.
False for multi-resolution or irregular geometries.
If true all image planes have the same pixel geometry.
"""
if self.nside.size > 1 or self.region == "explicit":
return False
else:
return True
@property
def center_coord(self):
"""Map coordinates of the center of the geometry (tuple)."""
return self._center_coord
@property
def center_pix(self):
"""Pixel coordinates of the center of the geometry (tuple)."""
return self._center_pix
@property
def center_skydir(self):
"""Sky coordinate of the center of the geometry.
Returns
-------
pix : `~astropy.coordinates.SkyCoord`
"""
return self._center_skydir
@property