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"""A HealpixSubsetSlicer - define the subset of healpixels to use to calculate metrics.""" 

2 

3from functools import wraps 

4import numpy as np 

5import healpy as hp 

6 

7import lsst.sims.utils as simsUtils 

8 

9from lsst.sims.maf.plots.spatialPlotters import HealpixSkyMap, HealpixHistogram, HealpixPowerSpectrum 

10 

11from .baseSpatialSlicer import BaseSpatialSlicer 

12 

13 

14__all__ = ['HealpixSubsetSlicer'] 

15 

16 

17class HealpixSubsetSlicer(BaseSpatialSlicer): 

18 """ 

19 A spatial slicer that evaluates pointings on a subset of a healpix-based grid. 

20 The advantage of using this healpixSubsetSlicer (rather than just putting the RA/Dec values into 

21 the UserPointsSlicer, which is another valid approach) is that you preserve the full healpix array. 

22 This means you could do things like calculate the power spectrum and plot without remapping into 

23 healpixels first. The downside is that you must first (externally) define the healpixels that you 

24 wish to use - the lsst.sims.featureScheduler.footprints is a useful add-on here. 

25 

26 When plotting with RA/Dec, the default HealpixSkyMap can be used, corresponding to 

27 {'rot': (0, 0, 0), 'flip': 'astro'}. 

28 

29 Parameters 

30 ---------- 

31 nside : int 

32 The nside parameter of the healpix grid. Must be a power of 2. 

33 hpid : np.ndarray 

34 The subset of healpix id's to use to calculate the metric. 

35 Because the hpid should be defined based on a particular nside, these first two 

36 arguments are not optional for this slicer. 

37 lonCol : str, optional 

38 Name of the longitude (RA equivalent) column to use from the input data. 

39 Default fieldRA 

40 latCol : str, optional 

41 Name of the latitude (Dec equivalent) column to use from the input data. 

42 Default fieldDec 

43 latLonDeg : boolean, optional 

44 Flag indicating whether the lat and lon values in the input data are in 

45 degrees (True) or radians (False). 

46 Default True. 

47 verbose : boolean, optional 

48 Flag to indicate whether or not to write additional information to stdout during runtime. 

49 Default True. 

50 badval : float, optional 

51 Bad value flag, relevant for plotting. Default the hp.UNSEEN value (in order to properly flag 

52 bad data points for plotting with the healpix plotting routines). This should not be changed. 

53 useCache : boolean 

54 Flag allowing the user to indicate whether or not to cache (and reuse) metric results 

55 calculated with the same set of simulated data pointings. 

56 This can be safely set to True for slicers not using maps and will result in increased speed. 

57 When calculating metric results using maps, the metadata at each individual ra/dec point may 

58 influence the metric results and so useCache should be set to False. 

59 Default True. 

60 leafsize : int, optional 

61 Leafsize value for kdtree. Default 100. 

62 radius : float, optional 

63 Radius for matching in the kdtree. Equivalent to the radius of the FOV. Degrees. 

64 Default 1.75. 

65 useCamera : boolean, optional 

66 Flag to indicate whether to use the LSST camera footprint or not. 

67 Default False. 

68 rotSkyPosColName : str, optional 

69 Name of the rotSkyPos column in the input data. Only used if useCamera is True. 

70 Describes the orientation of the camera orientation compared to the sky. 

71 Default rotSkyPos. 

72 mjdColName : str, optional 

73 Name of the exposure time column. Only used if useCamera is True. 

74 Default observationStartMJD. 

75 chipNames : array-like, optional 

76 List of chips to accept, if useCamera is True. This lets users turn 'on' only a subset of chips. 

77 Default 'all' - this uses all chips in the camera. 

78 """ 

79 def __init__(self, nside, hpid, lonCol ='fieldRA', 

80 latCol='fieldDec', latLonDeg=True, verbose=True, badval=hp.UNSEEN, 

81 useCache=True, leafsize=100, radius=1.75, 

82 useCamera=False, rotSkyPosColName='rotSkyPos', 

83 mjdColName='observationStartMJD', chipNames='all'): 

84 """Instantiate and set up healpix slicer object.""" 

85 super().__init__(verbose=verbose, 

86 lonCol=lonCol, latCol=latCol, 

87 badval=badval, radius=radius, leafsize=leafsize, 

88 useCamera=useCamera, rotSkyPosColName=rotSkyPosColName, 

89 mjdColName=mjdColName, chipNames=chipNames, latLonDeg=latLonDeg) 

90 # Valid values of nside are powers of 2. 

91 # nside=64 gives about 1 deg resolution 

92 # nside=256 gives about 13' resolution (~1 CCD) 

93 # nside=1024 gives about 3' resolution 

94 # Check validity of nside: 

95 if not(hp.isnsideok(nside)): 

96 raise ValueError('Valid values of nside are powers of 2.') 

97 if len(hpid) > hp.nside2npix(nside): 

98 raise ValueError('Nside (%d) and length of hpid (%d) seem incompatible.' % (nside, 

99 hp.nside2npix(nside))) 

100 self.nside = int(nside) 

101 self.hpid = hpid 

102 self.pixArea = hp.nside2pixarea(self.nside) 

103 self.nslice = hp.nside2npix(self.nside) 

104 self.spatialExtent = [0, self.nslice-1] 

105 self.shape = self.nslice 

106 if self.verbose: 

107 print('HealpixSubsetSlicer using NSIDE=%d, ' % (self.nside) + \ 

108 'approximate resolution %f arcminutes' % (hp.nside2resol(self.nside, arcmin=True))) 

109 # Set variables so slicer can be re-constructed 

110 self.slicer_init = {'nside': nside, 'hpid': hpid, 'lonCol': lonCol, 'latCol': latCol, 

111 'radius': radius} 

112 if useCache: 

113 # useCache set the size of the cache for the memoize function in sliceMetric. 

114 binRes = hp.nside2resol(nside) # Pixel size in radians 

115 # Set the cache size to be ~2x the circumference 

116 self.cacheSize = int(np.round(4.*np.pi/binRes)) 

117 # Set up slicePoint metadata. 

118 self.slicePoints['nside'] = nside 

119 self.slicePoints['sid'] = np.arange(self.nslice) 

120 self.slicePoints['ra'], self.slicePoints['dec'] = self._pix2radec(self.slicePoints['sid']) 

121 # Set the default plotting functions. 

122 self.plotFuncs = [HealpixSkyMap, HealpixHistogram, HealpixPowerSpectrum] 

123 

124 def __eq__(self, otherSlicer): 

125 """Evaluate if two slicers are equivalent.""" 

126 # If the two slicers are both healpix slicers, check nsides value. 

127 result = False 

128 if isinstance(otherSlicer, HealpixSubsetSlicer): 

129 if otherSlicer.nside == self.nside: 

130 if np.all(otherSlicer.hpid == self.hpid): 

131 if (otherSlicer.lonCol == self.lonCol and otherSlicer.latCol == self.latCol): 

132 if otherSlicer.radius == self.radius: 

133 if otherSlicer.useCamera == self.useCamera: 

134 if otherSlicer.chipsToUse == self.chipsToUse: 

135 if otherSlicer.rotSkyPosColName == self.rotSkyPosColName: 

136 if np.all(otherSlicer.shape == self.shape): 

137 result = True 

138 return result 

139 

140 def _pix2radec(self, islice): 

141 """Given the pixel number / sliceID, return the RA/Dec of the pointing, in radians.""" 

142 # Calculate RA/Dec in RADIANS of pixel in this healpix slicer. 

143 # Note that ipix could be an array, 

144 # in which case RA/Dec values will be an array also. 

145 lat, ra = hp.pix2ang(self.nside, islice) 

146 # Move dec to +/- 90 degrees 

147 dec = np.pi/2.0 - lat 

148 return ra, dec 

149 

150 # This slicer does iterate over all of the slicepoints - mainly so it can return a masked value for 

151 # non-calculated healpixels. 

152 def setupSlicer(self, simData, maps=None): 

153 """Use simData[self.lonCol] and simData[self.latCol] (in radians) to set up KDTree. 

154 

155 Parameters 

156 ----------- 

157 simData : numpy.recarray 

158 The simulated data, including the location of each pointing. 

159 maps : list of lsst.sims.maf.maps objects, optional 

160 List of maps (such as dust extinction) that will run to build up additional metadata at each 

161 slicePoint. This additional metadata is available to metrics via the slicePoint dictionary. 

162 Default None. 

163 """ 

164 if maps is not None: 

165 if self.cacheSize != 0 and len(maps) > 0: 

166 warnings.warn('Warning: Loading maps but cache on.' 

167 'Should probably set useCache=False in slicer.') 

168 self._runMaps(maps) 

169 self._setRad(self.radius) 

170 if self.useCamera: 

171 self._setupLSSTCamera() 

172 self._presliceFootprint(simData) 

173 else: 

174 if self.latLonDeg: 

175 self._buildTree(np.radians(simData[self.lonCol]), 

176 np.radians(simData[self.latCol]), self.leafsize) 

177 else: 

178 self._buildTree(simData[self.lonCol], simData[self.latCol], self.leafsize) 

179 

180 @wraps(self._sliceSimData) 

181 def _sliceSimData(islice): 

182 """Return indexes for relevant opsim data at slicepoint 

183 (slicepoint=lonCol/latCol value .. usually ra/dec).""" 

184 if islice not in self.hpid: 

185 return {'idxs': [], 'slicePoint': {self.slicePoints['sid'][islice], 

186 self.slicePoints['ra'][islice], 

187 self.slicePoints['dec'][islice]}} 

188 # Build dict for slicePoint info 

189 slicePoint = {} 

190 if self.useCamera: 

191 indices = self.sliceLookup[islice] 

192 slicePoint['chipNames'] = self.chipNames[islice] 

193 else: 

194 sx, sy, sz = simsUtils._xyz_from_ra_dec(self.slicePoints['ra'][islice], 

195 self.slicePoints['dec'][islice]) 

196 # Query against tree. 

197 indices = self.opsimtree.query_ball_point((sx, sy, sz), self.rad) 

198 

199 # Loop through all the slicePoint keys. If the first dimension of slicepoint[key] has 

200 # the same shape as the slicer, assume it is information per slicepoint. 

201 # Otherwise, pass the whole slicePoint[key] information. Useful for stellar LF maps 

202 # where we want to pass only the relevant LF and the bins that go with it. 

203 for key in self.slicePoints: 

204 if len(np.shape(self.slicePoints[key])) == 0: 

205 keyShape = 0 

206 else: 

207 keyShape = np.shape(self.slicePoints[key])[0] 

208 if (keyShape == self.nslice): 

209 slicePoint[key] = self.slicePoints[key][islice] 

210 else: 

211 slicePoint[key] = self.slicePoints[key] 

212 return {'idxs': indices, 'slicePoint': slicePoint} 

213 setattr(self, '_sliceSimData', _sliceSimData)