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 

2from .baseMetric import BaseMetric 

3import lsst.sims.maf.utils as mafUtils 

4import lsst.sims.utils as utils 

5from scipy.optimize import curve_fit 

6from builtins import str 

7 

8__all__ = ['ParallaxMetric', 'ProperMotionMetric', 'RadiusObsMetric', 

9 'ParallaxCoverageMetric', 'ParallaxDcrDegenMetric'] 

10 

11 

12class ParallaxMetric(BaseMetric): 

13 """Calculate the uncertainty in a parallax measurement given a series of observations. 

14 

15 Uses columns ra_pi_amp and dec_pi_amp, calculated by the ParallaxFactorStacker. 

16 

17 Parameters 

18 ---------- 

19 metricName : str, opt 

20 Default 'parallax'. 

21 m5Col : str, opt 

22 The default column name for m5 information in the input data. Default fiveSigmaDepth. 

23 filterCol : str, opt 

24 The column name for the filter information. Default filter. 

25 seeingCol : str, opt 

26 The column name for the seeing information. Since the astrometry errors are based on the physical 

27 size of the PSF, this should be the FWHM of the physical psf. Default seeingFwhmGeom. 

28 rmag : float, opt 

29 The r magnitude of the fiducial star in r band. Other filters are sclaed using sedTemplate keyword. 

30 Default 20.0 

31 SedTemplate : str, opt 

32 The template to use. This can be 'flat' or 'O','B','A','F','G','K','M'. Default flat. 

33 atm_err : float, opt 

34 The expected centroiding error due to the atmosphere, in arcseconds. Default 0.01. 

35 normalize : boolean, opt 

36 Compare the astrometric uncertainty to the uncertainty that would result if half the observations 

37 were taken at the start and half at the end. A perfect survey will have a value close to 1, while 

38 a poorly scheduled survey will be close to 0. Default False. 

39 badval : float, opt 

40 The value to return when the metric value cannot be calculated. Default -666. 

41 """ 

42 def __init__(self, metricName='parallax', m5Col='fiveSigmaDepth', 

43 filterCol='filter', seeingCol='seeingFwhmGeom', rmag=20., 

44 SedTemplate='flat', badval=-666, 

45 atm_err=0.01, normalize=False, **kwargs): 

46 Cols = [m5Col, filterCol, seeingCol, 'ra_pi_amp', 'dec_pi_amp'] 

47 if normalize: 

48 units = 'ratio' 

49 else: 

50 units = 'mas' 

51 super(ParallaxMetric, self).__init__(Cols, metricName=metricName, units=units, 

52 badval=badval, **kwargs) 

53 # set return type 

54 self.m5Col = m5Col 

55 self.seeingCol = seeingCol 

56 self.filterCol = filterCol 

57 filters = ['u', 'g', 'r', 'i', 'z', 'y'] 

58 self.mags = {} 

59 if SedTemplate == 'flat': 

60 for f in filters: 

61 self.mags[f] = rmag 

62 else: 

63 self.mags = utils.stellarMags(SedTemplate, rmag=rmag) 

64 self.atm_err = atm_err 

65 self.normalize = normalize 

66 self.comment = 'Estimated uncertainty in parallax measurement ' \ 

67 '(assuming no proper motion or that proper motion ' 

68 self.comment += 'is well fit). Uses measurements in all bandpasses, ' \ 

69 'and estimates astrometric error based on SNR ' 

70 self.comment += 'in each visit. ' 

71 if SedTemplate == 'flat': 

72 self.comment += 'Assumes a flat SED. ' 

73 if self.normalize: 

74 self.comment += 'This normalized version of the metric displays the ' \ 

75 'estimated uncertainty in the parallax measurement, ' 

76 self.comment += 'divided by the minimum parallax uncertainty possible ' \ 

77 '(if all visits were six ' 

78 self.comment += 'months apart). Values closer to 1 indicate more optimal ' \ 

79 'scheduling for parallax measurement.' 

80 

81 def _final_sigma(self, position_errors, ra_pi_amp, dec_pi_amp): 

82 """Assume parallax in RA and DEC are fit independently, then combined. 

83 All inputs assumed to be arcsec """ 

84 sigma_A = position_errors/ra_pi_amp 

85 sigma_B = position_errors/dec_pi_amp 

86 sigma_ra = np.sqrt(1./np.sum(1./sigma_A**2)) 

87 sigma_dec = np.sqrt(1./np.sum(1./sigma_B**2)) 

88 # Combine RA and Dec uncertainties, convert to mas 

89 sigma = np.sqrt(1./(1./sigma_ra**2+1./sigma_dec**2))*1e3 

90 return sigma 

91 

92 def run(self, dataslice, slicePoint=None): 

93 filters = np.unique(dataslice[self.filterCol]) 

94 if hasattr(filters[0], 'decode'): 

95 filters = [str(f.decode('utf-8')) for f in filters] 

96 snr = np.zeros(len(dataslice), dtype='float') 

97 # compute SNR for all observations 

98 for filt in filters: 

99 good = np.where(dataslice[self.filterCol] == filt) 

100 snr[good] = mafUtils.m52snr(self.mags[str(filt)], dataslice[self.m5Col][good]) 

101 position_errors = np.sqrt(mafUtils.astrom_precision(dataslice[self.seeingCol], 

102 snr)**2+self.atm_err**2) 

103 sigma = self._final_sigma(position_errors, dataslice['ra_pi_amp'], dataslice['dec_pi_amp']) 

104 if self.normalize: 

105 # Leave the dec parallax as zero since one can't have ra and dec maximized at the same time. 

106 sigma = self._final_sigma(position_errors, 

107 dataslice['ra_pi_amp']*0+1., dataslice['dec_pi_amp']*0)/sigma 

108 return sigma 

109 

110 

111class ProperMotionMetric(BaseMetric): 

112 """Calculate the uncertainty in the returned proper motion. 

113 

114 This metric assumes gaussian errors in the astrometry measurements. 

115 

116 Parameters 

117 ---------- 

118 metricName : str, opt 

119 Default 'properMotion'. 

120 m5Col : str, opt 

121 The default column name for m5 information in the input data. Default fiveSigmaDepth. 

122 mjdCol : str, opt 

123 The column name for the exposure time. Default observationStartMJD. 

124 filterCol : str, opt 

125 The column name for the filter information. Default filter. 

126 seeingCol : str, opt 

127 The column name for the seeing information. Since the astrometry errors are based on the physical 

128 size of the PSF, this should be the FWHM of the physical psf. Default seeingFwhmGeom. 

129 rmag : float, opt 

130 The r magnitude of the fiducial star in r band. Other filters are sclaed using sedTemplate keyword. 

131 Default 20.0 

132 SedTemplate : str, opt 

133 The template to use. This can be 'flat' or 'O','B','A','F','G','K','M'. Default flat. 

134 atm_err : float, opt 

135 The expected centroiding error due to the atmosphere, in arcseconds. Default 0.01. 

136 normalize : boolean, opt 

137 Compare the astrometric uncertainty to the uncertainty that would result if half the observations 

138 were taken at the start and half at the end. A perfect survey will have a value close to 1, while 

139 a poorly scheduled survey will be close to 0. Default False. 

140 baseline : float, opt 

141 The length of the survey used for the normalization, in years. Default 10. 

142 badval : float, opt 

143 The value to return when the metric value cannot be calculated. Default -666. 

144 """ 

145 def __init__(self, metricName='properMotion', 

146 m5Col='fiveSigmaDepth', mjdCol='observationStartMJD', 

147 filterCol='filter', seeingCol='seeingFwhmGeom', rmag=20., 

148 SedTemplate='flat', badval= -666, 

149 atm_err=0.01, normalize=False, 

150 baseline=10., **kwargs): 

151 cols = [m5Col, mjdCol, filterCol, seeingCol] 

152 if normalize: 

153 units = 'ratio' 

154 else: 

155 units = 'mas/yr' 

156 super(ProperMotionMetric, self).__init__(col=cols, metricName=metricName, units=units, 

157 badval=badval, **kwargs) 

158 # set return type 

159 self.mjdCol = mjdCol 

160 self.seeingCol = seeingCol 

161 self.m5Col = m5Col 

162 filters = ['u', 'g', 'r', 'i', 'z', 'y'] 

163 self.mags = {} 

164 if SedTemplate == 'flat': 

165 for f in filters: 

166 self.mags[f] = rmag 

167 else: 

168 self.mags = utils.stellarMags(SedTemplate, rmag=rmag) 

169 self.atm_err = atm_err 

170 self.normalize = normalize 

171 self.baseline = baseline 

172 self.comment = 'Estimated uncertainty of the proper motion fit ' \ 

173 '(assuming no parallax or that parallax is well fit). ' 

174 self.comment += 'Uses visits in all bands, and generates approximate ' \ 

175 'astrometric errors using the SNR in each visit. ' 

176 if SedTemplate == 'flat': 

177 self.comment += 'Assumes a flat SED. ' 

178 if self.normalize: 

179 self.comment += 'This normalized version of the metric represents ' \ 

180 'the estimated uncertainty in the proper ' 

181 self.comment += 'motion divided by the minimum uncertainty possible ' \ 

182 '(if all visits were ' 

183 self.comment += 'obtained on the first and last days of the survey). ' 

184 self.comment += 'Values closer to 1 indicate more optimal scheduling.' 

185 

186 def run(self, dataslice, slicePoint=None): 

187 filters = np.unique(dataslice['filter']) 

188 filters = [str(f) for f in filters] 

189 precis = np.zeros(dataslice.size, dtype='float') 

190 for f in filters: 

191 observations = np.where(dataslice['filter'] == f) 

192 if np.size(observations[0]) < 2: 

193 precis[observations] = self.badval 

194 else: 

195 snr = mafUtils.m52snr(self.mags[f], 

196 dataslice[self.m5Col][observations]) 

197 precis[observations] = mafUtils.astrom_precision( 

198 dataslice[self.seeingCol][observations], snr) 

199 precis[observations] = np.sqrt(precis[observations]**2 + self.atm_err**2) 

200 good = np.where(precis != self.badval) 

201 result = mafUtils.sigma_slope(dataslice[self.mjdCol][good], precis[good]) 

202 result = result*365.25*1e3 # Convert to mas/yr 

203 if (self.normalize) & (good[0].size > 0): 

204 new_dates = dataslice[self.mjdCol][good]*0 

205 nDates = new_dates.size 

206 new_dates[nDates//2:] = self.baseline*365.25 

207 result = (mafUtils.sigma_slope(new_dates, precis[good])*365.25*1e3)/result 

208 # Observations that are very close together can still fail 

209 if np.isnan(result): 

210 result = self.badval 

211 return result 

212 

213 

214class ParallaxCoverageMetric(BaseMetric): 

215 """ 

216 Check how well the parallax factor is distributed. Subtracts the weighted mean position of the 

217 parallax offsets, then computes the weighted mean radius of the points. 

218 If points are well distributed, the mean radius will be near 1. If phase coverage is bad, 

219 radius will be close to zero. 

220 

221 For points on the Ecliptic, uniform sampling should result in a metric value of ~0.5. 

222 At the poles, uniform sampling would result in a metric value of ~1. 

223 Conceptually, it is helpful to remember that the parallax motion of a star at the pole is 

224 a (nearly circular) ellipse while the motion of a star on the ecliptic is a straight line. Thus, any 

225 pair of observations separated by 6 months will give the full parallax range for a star on the pole 

226 but only observations on very specific dates will give the full range for a star on the ecliptic. 

227 

228 Optionally also demand that there are observations above the snrLimit kwarg spanning thetaRange radians. 

229 

230 Parameters 

231 ---------- 

232 m5Col: str, opt 

233 Column name for individual visit m5. Default fiveSigmaDepth. 

234 mjdCol: str, opt 

235 Column name for exposure time dates. Default observationStartMJD. 

236 filterCol: str, opt 

237 Column name for filter. Default filter. 

238 seeingCol: str, opt 

239 Column name for seeing (assumed FWHM). Default seeingFwhmGeom. 

240 rmag: float, opt 

241 Magnitude of fiducial star in r filter. Other filters are scaled using sedTemplate keyword. 

242 Default 20.0 

243 sedTemplate: str, opt 

244 Template to use (can be 'flat' or 'O','B','A','F','G','K','M'). Default 'flat'. 

245 atm_err: float, opt 

246 Centroiding error due to atmosphere in arcsec. Default 0.01 (arcseconds). 

247 thetaRange: float, opt 

248 Range of parallax offset angles to demand (in radians). Default=0 (means no range requirement). 

249 snrLimit: float, opt 

250 Only include points above the snrLimit when computing thetaRange. Default 5. 

251 

252 Returns 

253 -------- 

254 metricValue: float 

255 Returns a weighted mean of the length of the parallax factor vectors. 

256 Values near 1 imply that the points are well distributed. 

257 Values near 0 imply that the parallax phase coverage is bad. 

258 Near the ecliptic, uniform sampling results in metric values of about 0.5. 

259 

260 Notes 

261 ----- 

262 Uses the ParallaxFactor stacker to calculate ra_pi_amp and dec_pi_amp. 

263 """ 

264 def __init__(self, metricName='ParallaxCoverageMetric', m5Col='fiveSigmaDepth', 

265 mjdCol='observationStartMJD', filterCol='filter', seeingCol='seeingFwhmGeom', 

266 rmag=20., SedTemplate='flat', 

267 atm_err=0.01, thetaRange=0., snrLimit=5, **kwargs): 

268 cols = ['ra_pi_amp', 'dec_pi_amp', m5Col, mjdCol, filterCol, seeingCol] 

269 units = 'ratio' 

270 super(ParallaxCoverageMetric, self).__init__(cols, 

271 metricName=metricName, units=units, 

272 **kwargs) 

273 self.m5Col = m5Col 

274 self.seeingCol = seeingCol 

275 self.filterCol = filterCol 

276 self.mjdCol = mjdCol 

277 

278 # Demand the range of theta values 

279 self.thetaRange = thetaRange 

280 self.snrLimit = snrLimit 

281 

282 filters = ['u', 'g', 'r', 'i', 'z', 'y'] 

283 self.mags = {} 

284 if SedTemplate == 'flat': 

285 for f in filters: 

286 self.mags[f] = rmag 

287 else: 

288 self.mags = utils.stellarMags(SedTemplate, rmag=rmag) 

289 self.atm_err = atm_err 

290 caption = "Parallax factor coverage for an r=%.2f star (0 is bad, 0.5-1 is good). " % (rmag) 

291 caption += "One expects the parallax factor coverage to vary because stars on the ecliptic " 

292 caption += "can be observed when they have no parallax offset while stars at the pole are always " 

293 caption += "offset by the full parallax offset.""" 

294 self.comment = caption 

295 

296 def _thetaCheck(self, ra_pi_amp, dec_pi_amp, snr): 

297 good = np.where(snr >= self.snrLimit) 

298 theta = np.arctan2(dec_pi_amp[good], ra_pi_amp[good]) 

299 # Make values between 0 and 2pi 

300 theta = theta-np.min(theta) 

301 result = 0. 

302 if np.max(theta) >= self.thetaRange: 

303 # Check that things are in differnet quadrants 

304 theta = (theta+np.pi) % 2.*np.pi 

305 theta = theta-np.min(theta) 

306 if np.max(theta) >= self.thetaRange: 

307 result = 1 

308 return result 

309 

310 def _computeWeights(self, dataSlice, snr): 

311 # Compute centroid uncertainty in each visit 

312 position_errors = np.sqrt(mafUtils.astrom_precision(dataSlice[self.seeingCol], 

313 snr)**2+self.atm_err**2) 

314 weights = 1./position_errors**2 

315 return weights 

316 

317 def _weightedR(self, dec_pi_amp, ra_pi_amp, weights): 

318 ycoord = dec_pi_amp-np.average(dec_pi_amp, weights=weights) 

319 xcoord = ra_pi_amp-np.average(ra_pi_amp, weights=weights) 

320 radius = np.sqrt(xcoord**2+ycoord**2) 

321 aveRad = np.average(radius, weights=weights) 

322 return aveRad 

323 

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

325 if np.size(dataSlice) < 2: 

326 return self.badval 

327 

328 filters = np.unique(dataSlice[self.filterCol]) 

329 filters = [str(f) for f in filters] 

330 snr = np.zeros(len(dataSlice), dtype='float') 

331 # compute SNR for all observations 

332 for filt in filters: 

333 inFilt = np.where(dataSlice[self.filterCol] == filt) 

334 snr[inFilt] = mafUtils.m52snr(self.mags[str(filt)], dataSlice[self.m5Col][inFilt]) 

335 

336 weights = self._computeWeights(dataSlice, snr) 

337 aveR = self._weightedR(dataSlice['ra_pi_amp'], dataSlice['dec_pi_amp'], weights) 

338 if self.thetaRange > 0: 

339 thetaCheck = self._thetaCheck(dataSlice['ra_pi_amp'], dataSlice['dec_pi_amp'], snr) 

340 else: 

341 thetaCheck = 1. 

342 result = aveR*thetaCheck 

343 return result 

344 

345 

346class ParallaxDcrDegenMetric(BaseMetric): 

347 """Use the full parallax and DCR displacement vectors to find if they are degenerate. 

348 

349 Parameters 

350 ---------- 

351 metricName : str, opt 

352 Default 'ParallaxDcrDegenMetric'. 

353 seeingCol : str, opt 

354 Default 'FWHMgeom' 

355 m5Col : str, opt 

356 Default 'fiveSigmaDepth' 

357 filterCol : str 

358 Default 'filter' 

359 atm_err : float 

360 Minimum error in photometry centroids introduced by the atmosphere (arcseconds). Default 0.01. 

361 rmag : float 

362 r-band magnitude of the fiducual star that is being used (mag). 

363 SedTemplate : str 

364 The SED template to use for fiducia star colors, passed to lsst.sims.utils.stellarMags. 

365 Default 'flat' 

366 tol : float 

367 Tolerance for how well curve_fit needs to work before believing the covariance result. 

368 Default 0.05. 

369 

370 Returns 

371 ------- 

372 metricValue : float 

373 Returns the correlation coefficient between the best-fit parallax amplitude and DCR amplitude. 

374 The RA and Dec offsets are fit simultaneously. Values close to zero are good, values close to +/- 1 

375 are bad. Experience with fitting Monte Carlo simulations suggests the astrometric fits start 

376 becoming poor around a correlation of 0.7. 

377 """ 

378 def __init__(self, metricName='ParallaxDcrDegenMetric', seeingCol='seeingFwhmGeom', 

379 m5Col='fiveSigmaDepth', atm_err=0.01, rmag=20., SedTemplate='flat', 

380 filterCol='filter', tol=0.05, **kwargs): 

381 self.m5Col = m5Col 

382 self.seeingCol = seeingCol 

383 self.filterCol = filterCol 

384 self.tol = tol 

385 units = 'Correlation' 

386 # just put all the columns that all the stackers will need here? 

387 cols = ['ra_pi_amp', 'dec_pi_amp', 'ra_dcr_amp', 'dec_dcr_amp', 

388 seeingCol, m5Col] 

389 super(ParallaxDcrDegenMetric, self).__init__(cols, metricName=metricName, units=units, 

390 **kwargs) 

391 self.filters = ['u', 'g', 'r', 'i', 'z', 'y'] 

392 self.mags = {} 

393 if SedTemplate == 'flat': 

394 for f in self.filters: 

395 self.mags[f] = rmag 

396 else: 

397 self.mags = utils.stellarMags(SedTemplate, rmag=rmag) 

398 self.atm_err = atm_err 

399 

400 def _positions(self, x, a, b): 

401 """ 

402 Function to find parallax and dcr amplitudes 

403 

404 x should be a vector with [[parallax_x1, parallax_x2..., parallax_y1, parallax_y2...], 

405 [dcr_x1, dcr_x2..., dcr_y1, dcr_y2...]] 

406 """ 

407 result = a*x[0, :] + b*x[1, :] 

408 return result 

409 

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

411 # The idea here is that we calculate position errors (in RA and Dec) for all observations. 

412 # Then we generate arrays of the parallax offsets (delta RA parallax = ra_pi_amp, etc) 

413 # and the DCR offsets (delta RA DCR = ra_dcr_amp, etc), and just add them together into one 

414 # RA (and Dec) offset. Then, we try to fit for how we combined these offsets, but while 

415 # considering the astrometric noise. If we can figure out that we just added them together 

416 # (i.e. the curve_fit result is [a=1, b=1] for the function _positions above) 

417 # then we should be able to disentangle the parallax and DCR offsets when fitting 'for real'. 

418 # compute SNR for all observations 

419 snr = np.zeros(len(dataSlice), dtype='float') 

420 for filt in self.filters: 

421 inFilt = np.where(dataSlice[self.filterCol] == filt) 

422 snr[inFilt] = mafUtils.m52snr(self.mags[filt], dataSlice[self.m5Col][inFilt]) 

423 # Compute the centroiding uncertainties 

424 # Note that these centroiding uncertainties depend on the physical size of the PSF, thus 

425 # we are using seeingFwhmGeom for these metrics, not seeingFwhmEff. 

426 position_errors = np.sqrt(mafUtils.astrom_precision(dataSlice[self.seeingCol], snr)**2 + 

427 self.atm_err**2) 

428 # Construct the vectors of RA/Dec offsets. xdata is the "input data". ydata is the "output". 

429 xdata = np.empty((2, dataSlice.size * 2), dtype=float) 

430 xdata[0, :] = np.concatenate((dataSlice['ra_pi_amp'], dataSlice['dec_pi_amp'])) 

431 xdata[1, :] = np.concatenate((dataSlice['ra_dcr_amp'], dataSlice['dec_dcr_amp'])) 

432 ydata = np.sum(xdata, axis=0) 

433 # Use curve_fit to compute covariance between parallax and dcr amplitudes 

434 # Set the initial guess slightly off from the correct [1,1] to make sure it iterates. 

435 popt, pcov = curve_fit(self._positions, xdata, ydata, p0=[1.1, 0.9], 

436 sigma=np.concatenate((position_errors, position_errors)), 

437 absolute_sigma=True) 

438 # Catch if the fit failed to converge on the correct solution. 

439 if np.max(np.abs(popt - np.array([1., 1.]))) > self.tol: 

440 return self.badval 

441 # Covariance between best fit parallax amplitude and DCR amplitude. 

442 cov = pcov[1, 0] 

443 # Convert covarience between parallax and DCR amplitudes to normalized correlation 

444 perr = np.sqrt(np.diag(pcov)) 

445 correlation = cov/(perr[0]*perr[1]) 

446 result = correlation 

447 # This can throw infs. 

448 if np.isinf(result): 

449 result = self.badval 

450 return result 

451 

452 

453def calcDist_cosines(RA1, Dec1, RA2, Dec2): 

454 # Taken from simSelfCalib.py 

455 """Calculates distance on a sphere using spherical law of cosines. 

456 

457 Give this function RA/Dec values in radians. Returns angular distance(s), in radians. 

458 Note that since this is all numpy, you could input arrays of RA/Decs.""" 

459 # This formula can have rounding errors for case where distances are small. 

460 # Oh, the joys of wikipedia - http://en.wikipedia.org/wiki/Great-circle_distance 

461 # For the purposes of these calculations, this is probably accurate enough. 

462 D = np.sin(Dec2)*np.sin(Dec1) + np.cos(Dec1)*np.cos(Dec2)*np.cos(RA2-RA1) 

463 D = np.arccos(D) 

464 return D 

465 

466 

467class RadiusObsMetric(BaseMetric): 

468 """find the radius in the focal plane. returns things in degrees.""" 

469 

470 def __init__(self, metricName='radiusObs', raCol='fieldRA', decCol='fieldDec', 

471 units='radians', **kwargs): 

472 self.raCol = raCol 

473 self.decCol = decCol 

474 super(RadiusObsMetric, self).__init__(col=[self.raCol, self.decCol], 

475 metricName=metricName, units=units, **kwargs) 

476 

477 def run(self, dataSlice, slicePoint): 

478 ra = slicePoint['ra'] 

479 dec = slicePoint['dec'] 

480 distances = calcDist_cosines(ra, dec, np.radians(dataSlice[self.raCol]), 

481 np.radians(dataSlice[self.decCol])) 

482 distances = np.degrees(distances) 

483 return distances 

484 

485 def reduceMean(self, distances): 

486 return np.mean(distances) 

487 

488 def reduceRMS(self, distances): 

489 return np.std(distances) 

490 

491 def reduceFullRange(self, distances): 

492 return np.max(distances)-np.min(distances)