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

1import numpy as np 

2import healpy as hp 

3from scipy import interpolate 

4from .baseMetric import BaseMetric 

5 

6# A collection of metrics which are primarily intended to be used as summary statistics. 

7 

8__all__ = ['fOArea', 'fONv', 'TableFractionMetric', 'IdentityMetric', 

9 'NormalizeMetric', 'ZeropointMetric', 'TotalPowerMetric', 

10 'StaticProbesFoMEmulatorMetricSimple'] 

11 

12 

13class fONv(BaseMetric): 

14 """ 

15 Metrics based on a specified area, but returning NVISITS related to area: 

16 given Asky, what is the minimum and median number of visits obtained over that much area? 

17 (choose the portion of the sky with the highest number of visits first). 

18 

19 Parameters 

20 ---------- 

21 col : str or list of strs, opt 

22 Name of the column in the numpy recarray passed to the summary metric. 

23 Asky : float, opt 

24 Area of the sky to base the evaluation of number of visits over. 

25 Default 18,0000 sq deg. 

26 nside : int, opt 

27 Nside parameter from healpix slicer, used to set the physical relationship between on-sky area 

28 and number of healpixels. Default 128. 

29 Nvisit : int, opt 

30 Number of visits to use as the benchmark value, if choosing to return a normalized Nvisit value. 

31 norm : boolean, opt 

32 Normalize the returned "nvisit" (min / median) values by Nvisit, if true. 

33 Default False. 

34 metricName : str, opt 

35 Name of the summary metric. Default fONv. 

36 """ 

37 def __init__(self, col='metricdata', Asky=18000., nside=128, Nvisit=825, 

38 norm=False, metricName='fONv', **kwargs): 

39 """Asky = square degrees """ 

40 super().__init__(col=col, metricName=metricName, **kwargs) 

41 self.Nvisit = Nvisit 

42 self.nside = nside 

43 # Determine how many healpixels are included in Asky sq deg. 

44 self.Asky = Asky 

45 self.scale = hp.nside2pixarea(self.nside, degrees=True) 

46 self.npix_Asky = np.int(np.ceil(self.Asky / self.scale)) 

47 self.norm = norm 

48 

49 def run(self, dataSlice, slicePoint=None): 

50 result = np.empty(2, dtype=[('name', np.str_, 20), ('value', float)]) 

51 result['name'][0] = "MedianNvis" 

52 result['name'][1] = "MinNvis" 

53 # If there is not even as much data as needed to cover Asky: 

54 if len(dataSlice) < self.npix_Asky: 

55 # Return the same type of metric value, to make it easier downstream. 

56 result['value'][0] = self.badval 

57 result['value'][1] = self.badval 

58 return result 

59 # Otherwise, calculate median and mean Nvis: 

60 name = dataSlice.dtype.names[0] 

61 nvis_sorted = np.sort(dataSlice[name]) 

62 # Find the Asky's worth of healpixels with the largest # of visits. 

63 nvis_Asky = nvis_sorted[-self.npix_Asky:] 

64 result['value'][0] = np.median(nvis_Asky) 

65 result['value'][1] = np.min(nvis_Asky) 

66 if self.norm: 

67 result['value'] /= float(self.Nvisit) 

68 return result 

69 

70 

71class fOArea(BaseMetric): 

72 """ 

73 Metrics based on a specified number of visits, but returning AREA related to Nvisits: 

74 given Nvisit, what amount of sky is covered with at least that many visits? 

75 

76 Parameters 

77 ---------- 

78 col : str or list of strs, opt 

79 Name of the column in the numpy recarray passed to the summary metric. 

80 Nvisit : int, opt 

81 Number of visits to use as the minimum required -- metric calculated area that has this many visits. 

82 Default 825. 

83 Asky : float, opt 

84 Area to use as the benchmark value, if choosing to returned a normalized Area value. 

85 Default 18,0000 sq deg. 

86 nside : int, opt 

87 Nside parameter from healpix slicer, used to set the physical relationship between on-sky area 

88 and number of healpixels. Default 128. 

89 norm : boolean, opt 

90 Normalize the returned "area" (area with minimum Nvisit visits) value by Asky, if true. 

91 Default False. 

92 metricName : str, opt 

93 Name of the summary metric. Default fOArea. 

94 """ 

95 def __init__(self, col='metricdata', Nvisit=825, Asky = 18000.0, nside=128, 

96 norm=False, metricName='fOArea', **kwargs): 

97 """Asky = square degrees """ 

98 super().__init__(col=col, metricName=metricName, **kwargs) 

99 self.Nvisit = Nvisit 

100 self.nside = nside 

101 self.Asky = Asky 

102 self.scale = hp.nside2pixarea(self.nside, degrees=True) 

103 self.norm = norm 

104 

105 def run(self, dataSlice, slicePoint=None): 

106 name = dataSlice.dtype.names[0] 

107 nvis_sorted = np.sort(dataSlice[name]) 

108 # Identify the healpixels with more than Nvisits. 

109 nvis_min = nvis_sorted[np.where(nvis_sorted >= self.Nvisit)] 

110 if len(nvis_min) == 0: 

111 result = self.badval 

112 else: 

113 result = nvis_min.size * self.scale 

114 if self.norm: 

115 result /= float(self.Asky) 

116 return result 

117 

118 

119class TableFractionMetric(BaseMetric): 

120 """ 

121 Count the completeness (for many fields) and summarize how many fields have given completeness levels 

122 (within a series of bins). Works with completenessMetric only. 

123 

124 This metric is meant to be used as a summary statistic on something like the completeness metric. 

125 The output is DIFFERENT FROM SSTAR and is: 

126 element matching values 

127 0 0 == P 

128 1 0 < P < .1 

129 2 .1 <= P < .2 

130 3 .2 <= P < .3 

131 ... 

132 10 .9 <= P < 1 

133 11 1 == P 

134 12 1 < P 

135 Note the 1st and last elements do NOT obey the numpy histogram conventions. 

136 """ 

137 def __init__(self, col='metricdata', nbins=10, maskVal=0.): 

138 """ 

139 colname = the column name in the metric data (i.e. 'metricdata' usually). 

140 nbins = number of bins between 0 and 1. Should divide evenly into 100. 

141 """ 

142 super(TableFractionMetric, self).__init__(col=col, maskVal=maskVal, metricDtype='float') 

143 self.nbins = nbins 

144 

145 def run(self, dataSlice, slicePoint=None): 

146 # Calculate histogram of completeness values that fall between 0-1. 

147 goodVals = np.where((dataSlice[self.colname] > 0) & (dataSlice[self.colname] < 1) ) 

148 bins = np.arange(self.nbins+1.)/self.nbins 

149 hist, b = np.histogram(dataSlice[self.colname][goodVals], bins=bins) 

150 # Fill in values for exact 0, exact 1 and >1. 

151 zero = np.size(np.where(dataSlice[self.colname] == 0)[0]) 

152 one = np.size(np.where(dataSlice[self.colname] == 1)[0]) 

153 overone = np.size(np.where(dataSlice[self.colname] > 1)[0]) 

154 hist = np.concatenate((np.array([zero]), hist, np.array([one]), np.array([overone]))) 

155 # Create labels for each value 

156 binNames = ['0 == P'] 

157 binNames.append('0 < P < 0.1') 

158 for i in np.arange(1, self.nbins): 

159 binNames.append('%.2g <= P < %.2g'%(b[i], b[i+1]) ) 

160 binNames.append('1 == P') 

161 binNames.append('1 < P') 

162 # Package the names and values up 

163 result = np.empty(hist.size, dtype=[('name', np.str_, 20), ('value', float)]) 

164 result['name'] = binNames 

165 result['value'] = hist 

166 return result 

167 

168 

169class IdentityMetric(BaseMetric): 

170 """ 

171 Return the metric value itself .. this is primarily useful as a summary statistic for UniSlicer metrics. 

172 """ 

173 def run(self, dataSlice, slicePoint=None): 

174 if len(dataSlice[self.colname]) == 1: 

175 result = dataSlice[self.colname][0] 

176 else: 

177 result = dataSlice[self.colname] 

178 return result 

179 

180 

181class NormalizeMetric(BaseMetric): 

182 """ 

183 Return a metric values divided by 'normVal'. Useful for turning summary statistics into fractions. 

184 """ 

185 def __init__(self, col='metricdata', normVal=1, **kwargs): 

186 super(NormalizeMetric, self).__init__(col=col, **kwargs) 

187 self.normVal = float(normVal) 

188 def run(self, dataSlice, slicePoint=None): 

189 result = dataSlice[self.colname]/self.normVal 

190 if len(result) == 1: 

191 return result[0] 

192 else: 

193 return result 

194 

195class ZeropointMetric(BaseMetric): 

196 """ 

197 Return a metric values with the addition of 'zp'. Useful for altering the zeropoint for summary statistics. 

198 """ 

199 def __init__(self, col='metricdata', zp=0, **kwargs): 

200 super(ZeropointMetric, self).__init__(col=col, **kwargs) 

201 self.zp = zp 

202 def run(self, dataSlice, slicePoint=None): 

203 result = dataSlice[self.colname] + self.zp 

204 if len(result) == 1: 

205 return result[0] 

206 else: 

207 return result 

208 

209class TotalPowerMetric(BaseMetric): 

210 """ 

211 Calculate the total power in the angular power spectrum between lmin/lmax. 

212 """ 

213 def __init__(self, col='metricdata', lmin=100., lmax=300., removeDipole=True, maskVal=hp.UNSEEN, **kwargs): 

214 self.lmin = lmin 

215 self.lmax = lmax 

216 self.removeDipole = removeDipole 

217 super(TotalPowerMetric, self).__init__(col=col, maskVal=maskVal, **kwargs) 

218 

219 def run(self, dataSlice, slicePoint=None): 

220 # Calculate the power spectrum. 

221 if self.removeDipole: 

222 cl = hp.anafast(hp.remove_dipole(dataSlice[self.colname], verbose=False)) 

223 else: 

224 cl = hp.anafast(dataSlice[self.colname]) 

225 ell = np.arange(np.size(cl)) 

226 condition = np.where((ell <= self.lmax) & (ell >= self.lmin))[0] 

227 totalpower = np.sum(cl[condition]*(2*ell[condition]+1)) 

228 return totalpower 

229 

230 

231class StaticProbesFoMEmulatorMetricSimple(BaseMetric): 

232 """This calculates the Figure of Merit for the combined 

233 static probes (3x2pt, i.e., Weak Lensing, LSS, Clustering). 

234 This FoM is purely statistical and does not factor in systematics. 

235  

236 This version of the emulator was used to generate the results in 

237 https://ui.adsabs.harvard.edu/abs/2018arXiv181200515L/abstract 

238  

239 A newer version is being created. This version has been renamed  

240 Simple in anticipation of the newer, more sophisticated metric 

241 replacing it. 

242 """ 

243 def __init__(self, nside=128, year=10, col=None, **kwargs): 

244 

245 """ 

246 Args: 

247 nside (int): healpix resolution 

248 year (int): year of the FoM emulated values,  

249 can be one of [1, 3, 6, 10] 

250 col (str): column name of metric data. 

251 """ 

252 self.nside = nside 

253 super().__init__(metricName='StaticProbesFoMEmulatorMetric', col=col, **kwargs) 

254 if col is None: 

255 self.col = 'metricdata' 

256 self.year = year 

257 

258 def run(self, dataSlice, slicePoint=None): 

259 """ 

260 Args: 

261 dataSlice (ndarray): Values passed to metric by the slicer,  

262 which the metric will use to calculate metric values  

263 at each slicePoint. 

264 slicePoint (Dict): Dictionary of slicePoint metadata passed 

265 to each metric. 

266 Returns: 

267 float: Interpolated static-probe statistical Figure-of-Merit. 

268 Raises: 

269 ValueError: If year is not one of the 4 for which a FoM is calculated 

270 """ 

271 # Chop off any outliers 

272 good_pix = np.where(dataSlice[self.col] > 0)[0] 

273 

274 # Calculate area and med depth from 

275 area = hp.nside2pixarea(self.nside, degrees=True) * np.size(good_pix) 

276 median_depth = np.median(dataSlice[self.col][good_pix]) 

277 

278 # FoM is calculated at the following values 

279 if self.year == 1: 

280 areas = [7500, 13000, 16000] 

281 depths = [24.9, 25.2, 25.5] 

282 fom_arr = [ 

283 [1.212257e+02, 1.462689e+02, 1.744913e+02], 

284 [1.930906e+02, 2.365094e+02, 2.849131e+02], 

285 [2.316956e+02, 2.851547e+02, 3.445717e+02] 

286 ] 

287 elif self.year == 3: 

288 areas = [10000, 15000, 20000] 

289 depths = [25.5, 25.8, 26.1] 

290 fom_arr = [ 

291 [1.710645e+02, 2.246047e+02, 2.431472e+02], 

292 [2.445209e+02, 3.250737e+02, 3.516395e+02], 

293 [3.173144e+02, 4.249317e+02, 4.595133e+02] 

294 ] 

295 

296 elif self.year == 6: 

297 areas = [10000, 15000, 2000] 

298 depths = [25.9, 26.1, 26.3] 

299 fom_arr = [ 

300 [2.346060e+02, 2.414678e+02, 2.852043e+02], 

301 [3.402318e+02, 3.493120e+02, 4.148814e+02], 

302 [4.452766e+02, 4.565497e+02, 5.436992e+02] 

303 ] 

304 

305 elif self.year == 10: 

306 areas = [10000, 15000, 20000] 

307 depths = [26.3, 26.5, 26.7] 

308 fom_arr = [ 

309 [2.887266e+02, 2.953230e+02, 3.361616e+02], 

310 [4.200093e+02, 4.292111e+02, 4.905306e+02], 

311 [5.504419e+02, 5.624697e+02, 6.441837e+02] 

312 ] 

313 else: 

314 raise ValueError('FoMEmulator is not defined for this year') 

315 

316 # Interpolate FoM to the actual values for this sim 

317 areas = [[i]*3 for i in areas] 

318 depths = [depths]*3 

319 f = interpolate.interp2d(areas, depths, fom_arr, bounds_error=False) 

320 fom = f(area, median_depth)[0] 

321 return fom 

322