Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1from __future__ import division 

2import numpy as np 

3import healpy as hp 

4 

5__all__ = ['hpid2RaDec', 'raDec2Hpid', 'healbin', '_hpid2RaDec', '_raDec2Hpid', '_healbin', 'moc2array'] 

6 

7 

8def _hpid2RaDec(nside, hpids, **kwargs): 

9 """ 

10 Correct for healpy being silly and running dec from 0-180. 

11 

12 Parameters 

13 ---------- 

14 nside : int 

15 Must be a value of 2^N. 

16 hpids : np.array 

17 Array (or single value) of healpixel IDs. 

18 

19 Returns 

20 ------- 

21 raRet : float (or np.array) 

22 RA positions of the input healpixel IDs. In radians. 

23 decRet : float (or np.array) 

24 Dec positions of the input healpixel IDs. In radians. 

25 """ 

26 

27 lat, lon = hp.pix2ang(nside, hpids, **kwargs) 

28 decRet = np.pi / 2.0 - lat 

29 raRet = lon 

30 

31 return raRet, decRet 

32 

33 

34def hpid2RaDec(nside, hpids, **kwargs): 

35 """ 

36 Correct for healpy being silly and running dec from 0-180. 

37 

38 Parameters 

39 ---------- 

40 nside : int 

41 Must be a value of 2^N. 

42 hpids : np.array 

43 Array (or single value) of healpixel IDs. 

44 

45 Returns 

46 ------- 

47 raRet : float (or np.array) 

48 RA positions of the input healpixel IDs. In degrees. 

49 decRet : float (or np.array) 

50 Dec positions of the input healpixel IDs. In degrees. 

51 """ 

52 ra, dec = _hpid2RaDec(nside, hpids, **kwargs) 

53 return np.degrees(ra), np.degrees(dec) 

54 

55 

56def _raDec2Hpid(nside, ra, dec, **kwargs): 

57 """ 

58 Assign ra,dec points to the correct healpixel. 

59 

60 Parameters 

61 ---------- 

62 nside : int 

63 Must be a value of 2^N. 

64 ra : np.array 

65 RA values to assign to healpixels. Radians. 

66 dec : np.array 

67 Dec values to assign to healpixels. Radians. 

68 

69 Returns 

70 ------- 

71 hpids : np.array 

72 Healpixel IDs for the input positions. 

73 """ 

74 lat = np.pi / 2.0 - dec 

75 hpids = hp.ang2pix(nside, lat, ra, **kwargs) 

76 return hpids 

77 

78 

79def raDec2Hpid(nside, ra, dec, **kwargs): 

80 """ 

81 Assign ra,dec points to the correct healpixel. 

82 

83 Parameters 

84 ---------- 

85 nside : int 

86 Must be a value of 2^N. 

87 ra : np.array 

88 RA values to assign to healpixels. Degrees. 

89 dec : np.array 

90 Dec values to assign to healpixels. Degrees. 

91 

92 Returns 

93 ------- 

94 hpids : np.array 

95 Healpixel IDs for the input positions. 

96 """ 

97 return _raDec2Hpid(nside, np.radians(ra), np.radians(dec), **kwargs) 

98 

99 

100def _healbin(ra, dec, values, nside=128, reduceFunc=np.mean, dtype=float, fillVal=hp.UNSEEN): 

101 """ 

102 Take arrays of ra's, dec's, and value and bin into healpixels. Like numpy.hexbin but for 

103 bins on a sphere. 

104 

105 Parameters 

106 ---------- 

107 ra : np.array 

108 RA positions of the data points. Radians. 

109 dec : np.array 

110 Dec positions of the data points. Radians 

111 values : np.array 

112 The values at each ra,dec position. 

113 nside : int 

114 Healpixel nside resolution. Must be a value of 2^N. 

115 reduceFunc : function (numpy.mean) 

116 A function that will return a single value given a subset of `values`. 

117 dtype : dtype ('float') 

118 Data type of the resulting mask 

119 

120 Returns 

121 ------- 

122 mapVals : np.array 

123 A numpy array that is a valid Healpixel map. 

124 """ 

125 

126 hpids = _raDec2Hpid(nside, ra, dec) 

127 

128 order = np.argsort(hpids) 

129 hpids = hpids[order] 

130 values = values[order] 

131 pixids = np.unique(hpids) 

132 

133 left = np.searchsorted(hpids, pixids) 

134 right = np.searchsorted(hpids, pixids, side='right') 

135 

136 mapVals = np.zeros(hp.nside2npix(nside), dtype=dtype) + fillVal 

137 

138 # Wow, I thought histogram would be faster than the loop, but this has been faster! 

139 for i, idx in enumerate(pixids): 

140 mapVals[idx] = reduceFunc(values[left[i]:right[i]]) 

141 

142 # Change any NaNs to fill value 

143 mapVals[np.isnan(mapVals)] = fillVal 

144 

145 return mapVals 

146 

147 

148def healbin(ra, dec, values, nside=128, reduceFunc=np.mean, dtype=float): 

149 """ 

150 Take arrays of ra's, dec's, and value and bin into healpixels. Like numpy.hexbin but for 

151 bins on a sphere. 

152 

153 Parameters 

154 ---------- 

155 ra : np.array 

156 RA positions of the data points. Degrees. 

157 dec : np.array 

158 Dec positions of the data points. Degrees. 

159 values : np.array 

160 The values at each ra,dec position. 

161 nside : int 

162 Healpixel nside resolution. Must be a value of 2^N. 

163 reduceFunc : function (numpy.mean) 

164 A function that will return a single value given a subset of `values`. 

165 dtype : dtype ('float') 

166 Data type of the resulting mask 

167 

168 Returns 

169 ------- 

170 mapVals : np.array 

171 A numpy array that is a valid Healpixel map. 

172 """ 

173 return _healbin(np.radians(ra), np.radians(dec), values, nside=nside, 

174 reduceFunc=reduceFunc, dtype=dtype) 

175 

176 

177def moc2array(data, uniq, nside=128, reduceFunc=np.sum, density=True, fillVal=0.): 

178 """Convert a Multi-Order Coverage Map to a single nside HEALPix array. Useful 

179 for converting maps output by LIGO alerts. Expect that future versions of 

180 healpy or astropy will be able to replace this functionality. Note that this is 

181 a convienence function that will probably degrade portions of the MOC that are 

182 sampled at high resolution. 

183 

184 Details of HEALPix Mulit-Order Coverage map: http://ivoa.net/documents/MOC/20190404/PR-MOC-1.1-20190404.pdf 

185 

186 Parameters 

187 ---------- 

188 data : np.array 

189 Data values for the MOC map 

190 uniq : np.array 

191 The UNIQ values for the MOC map 

192 nside : int (128) 

193 The output map nside 

194 reduceFunc : function (np.sum) 

195 The function to use to combine data into single healpixels. 

196 density : bool (True) 

197 If True, multiplies data values by pixel area before applying reduceFunc, and divides 

198 the final array by the output pixel area. Should be True if working on a probability density MOC. 

199 fillVal : float (0.) 

200 Value to fill empty HEALPixels with. Good choices include 0 (default), hp.UNSEEN, and np.nan. 

201 

202 Returns 

203 ------- 

204 np.array : HEALpy array of nside. Units should be the same as the input map as processed by reduceFunc. 

205 """ 

206 

207 # NUNIQ packing, from page 12 of http://ivoa.net/documents/MOC/20190404/PR-MOC-1.1-20190404.pdf 

208 orders = np.floor(np.log2(uniq / 4) / 2).astype(int) 

209 npixs = (uniq - 4 * 4**orders).astype(int) 

210 

211 nsides = 2**orders 

212 names = ['ra', 'dec', 'area'] 

213 types = [float]*len(names) 

214 data_points = np.zeros(data.size, dtype=list(zip(names, types))) 

215 for order in np.unique(orders): 

216 good = np.where(orders == order) 

217 ra, dec = _hpid2RaDec(nsides[good][0], npixs[good], nest=True) 

218 data_points['ra'][good] = ra 

219 data_points['dec'][good] = dec 

220 data_points['area'][good] = hp.nside2pixarea(nsides[good][0]) 

221 

222 if density: 

223 tobin_data = data*data_points['area'] 

224 else: 

225 tobin_data = data 

226 

227 result = _healbin(data_points['ra'], data_points['dec'], tobin_data, nside=nside, 

228 reduceFunc=reduceFunc, fillVal=fillVal) 

229 

230 if density: 

231 good = np.where(result != fillVal) 

232 result[good] = result[good] / hp.nside2pixarea(nside) 

233 

234 return result