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

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

from __future__ import division 

import numpy as np 

import healpy as hp 

 

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

 

 

def _hpid2RaDec(nside, hpids): 

""" 

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

 

Parameters 

---------- 

nside : int 

Must be a value of 2^N. 

hpids : np.array 

Array (or single value) of healpixel IDs. 

 

Returns 

------- 

raRet : float (or np.array) 

RA positions of the input healpixel IDs. In radians. 

decRet : float (or np.array) 

Dec positions of the input healpixel IDs. In radians. 

""" 

 

lat, lon = hp.pix2ang(nside, hpids) 

decRet = np.pi / 2.0 - lat 

raRet = lon 

 

return raRet, decRet 

 

 

def hpid2RaDec(nside, hpids): 

""" 

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

 

Parameters 

---------- 

nside : int 

Must be a value of 2^N. 

hpids : np.array 

Array (or single value) of healpixel IDs. 

 

Returns 

------- 

raRet : float (or np.array) 

RA positions of the input healpixel IDs. In degrees. 

decRet : float (or np.array) 

Dec positions of the input healpixel IDs. In degrees. 

""" 

ra, dec = _hpid2RaDec(nside, hpids) 

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

 

 

def _raDec2Hpid(nside, ra, dec): 

""" 

Assign ra,dec points to the correct healpixel. 

 

Parameters 

---------- 

nside : int 

Must be a value of 2^N. 

ra : np.array 

RA values to assign to healpixels. Radians. 

dec : np.array 

Dec values to assign to healpixels. Radians. 

 

Returns 

------- 

hpids : np.array 

Healpixel IDs for the input positions. 

""" 

lat = np.pi / 2.0 - dec 

hpids = hp.ang2pix(nside, lat, ra) 

return hpids 

 

 

def raDec2Hpid(nside, ra, dec): 

""" 

Assign ra,dec points to the correct healpixel. 

 

Parameters 

---------- 

nside : int 

Must be a value of 2^N. 

ra : np.array 

RA values to assign to healpixels. Degrees. 

dec : np.array 

Dec values to assign to healpixels. Degrees. 

 

Returns 

------- 

hpids : np.array 

Healpixel IDs for the input positions. 

""" 

return _raDec2Hpid(nside, np.radians(ra), np.radians(dec)) 

 

 

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

""" 

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

bins on a sphere. 

 

Parameters 

---------- 

ra : np.array 

RA positions of the data points. Radians. 

dec : np.array 

Dec positions of the data points. Radians 

values : np.array 

The values at each ra,dec position. 

nside : int 

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

reduceFunc : function (numpy.mean) 

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

dtype : dtype ('float') 

Data type of the resulting mask 

 

Returns 

------- 

mapVals : np.array 

A numpy array that is a valid Healpixel map. 

""" 

 

hpids = _raDec2Hpid(nside, ra, dec) 

 

order = np.argsort(hpids) 

hpids = hpids[order] 

values = values[order] 

pixids = np.unique(hpids) 

 

left = np.searchsorted(hpids, pixids) 

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

 

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

 

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

for i, idx in enumerate(pixids): 

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

 

# Change any NaNs to healpy mask value 

mapVals[np.isnan(mapVals)] = hp.UNSEEN 

 

return mapVals 

 

 

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

""" 

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

bins on a sphere. 

 

Parameters 

---------- 

ra : np.array 

RA positions of the data points. Degrees. 

dec : np.array 

Dec positions of the data points. Degrees. 

values : np.array 

The values at each ra,dec position. 

nside : int 

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

reduceFunc : function (numpy.mean) 

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

dtype : dtype ('float') 

Data type of the resulting mask 

 

Returns 

------- 

mapVals : np.array 

A numpy array that is a valid Healpixel map. 

""" 

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

reduceFunc=reduceFunc, dtype=dtype)