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 warnings 

2import numpy as np 

3import palpy 

4from lsst.sims.utils import Site, m5_flat_sed, xyz_from_ra_dec, xyz_angular_radius, \ 

5 _buildTree, _xyz_from_ra_dec 

6from lsst.sims.survey.fields import FieldsDatabase 

7from .baseStacker import BaseStacker 

8 

9__all__ = ['NormAirmassStacker', 'ParallaxFactorStacker', 'HourAngleStacker', 

10 'FilterColorStacker', 'ZenithDistStacker', 'ParallacticAngleStacker', 

11 'DcrStacker', 'FiveSigmaStacker', 'OpSimFieldStacker', 

12 'SaturationStacker'] 

13 

14# Original stackers by Peter Yoachim (yoachim@uw.edu) 

15# Filter color stacker by Lynne Jones (lynnej@uw.edu) 

16 

17 

18class SaturationStacker(BaseStacker): 

19 """Calculate the saturation limit of a point source. Assumes Guassian PSF. 

20  

21 Parameters 

22 ---------- 

23 pixscale : float, opt (0.2) 

24 Arcsec per pixel 

25 gain : float, opt (2.3) 

26 electrons per adu 

27 saturation_e : float, opt (150e3) 

28 The saturation level in electrons 

29 zeropoints : dict-like, opt (None) 

30 The zeropoints for the telescope. Keys should be str with filter names, values in mags. 

31 If None, will use Rubin-like zeropoints. 

32 km : dict-like, opt (None) 

33 Atmospheric extinction values. Keys should be str with filter names. If None, will use Rubin-like zeropoints. 

34 """ 

35 colsAdded = ['saturation_mag'] 

36 

37 def __init__(self, seeingCol='seeingFwhmEff', skybrightnessCol='skyBrightness', 

38 exptimeCol='visitExposureTime', nexpCol='numExposures', 

39 filterCol='filter', airmassCol='airmass', 

40 saturation_e=150e3, zeropoints=None, km=None, pixscale=0.2, gain=1.0): 

41 self.units = ['mag'] 

42 self.colsReq = [seeingCol, skybrightnessCol, exptimeCol, nexpCol, filterCol, airmassCol] 

43 self.seeingCol = seeingCol 

44 self.skybrightnessCol = skybrightnessCol 

45 self.exptimeCol = exptimeCol 

46 self.nexpCol = nexpCol 

47 self.filterCol = filterCol 

48 self.airmassCol = airmassCol 

49 self.saturation_adu = saturation_e/gain 

50 self.pixscale = 0.2 

51 names = ['u', 'g', 'r', 'i', 'z', 'y'] 

52 types = [float]*6 

53 if zeropoints is None: 53 ↛ 59line 53 didn't jump to line 59, because the condition on line 53 was never false

54 # Note these zeropoints are calculating the number of *electrons* per second (thus gain=1) 

55 # https://github.com/lsst-pst/syseng_throughputs/blob/master/notebooks/Syseng%20Throughputs%20Repo%20Demo.ipynb 

56 self.zeropoints = np.array([27.03, 28.38, 28.15, 27.86, 27.46, 26.68]).view(list(zip(names, types))) 

57 self.saturation_adu = saturation_e 

58 else: 

59 self.zeropoints = zeropoints 

60 

61 if km is None: 61 ↛ 65line 61 didn't jump to line 65, because the condition on line 61 was never false

62 # Also from notebook above 

63 self.km = np.array([0.491, 0.213, 0.126, 0.096, 0.069, 0.170]).view(list(zip(names, types))) 

64 else: 

65 self.km = km 

66 

67 def _run(self, simData, cols_present=False): 

68 for filtername in np.unique(simData[self.filterCol]): 

69 in_filt = np.where(simData[self.filterCol] == filtername)[0] 

70 # Calculate the length of the on-sky time per EXPOSURE 

71 exptime = simData[self.exptimeCol][in_filt] / simData[self.nexpCol][in_filt] 

72 # Calculate sky counts per pixel per second from skybrightness + zeropoint (e/1s) 

73 sky_counts = 10.**(0.4*(self.zeropoints[filtername] - simData[self.skybrightnessCol][in_filt])) * self.pixscale**2 

74 # Total sky counts in each exposure  

75 sky_counts = sky_counts * exptime 

76 # The counts available to the source (at peak) in each exposure is the difference between saturation and sky 

77 remaining_counts_peak = (self.saturation_adu - sky_counts) 

78 # Now to figure out how many counts there would be total, if there are that many in the peak 

79 sigma = simData[self.seeingCol][in_filt]/2.354 

80 source_counts = remaining_counts_peak * 2.*np.pi*(sigma/self.pixscale)**2 

81 # source counts = counts per exposure (expTimeCol / nexp) 

82 # Translate to counts per second, to apply zeropoint  

83 count_rate = source_counts / exptime 

84 simData['saturation_mag'][in_filt] = -2.5*np.log10(count_rate) + self.zeropoints[filtername] 

85 # Airmass correction 

86 simData['saturation_mag'][in_filt] -= self.km[filtername]*(simData[self.airmassCol][in_filt] - 1.) 

87 

88 return simData 

89 

90 

91class FiveSigmaStacker(BaseStacker): 

92 """ 

93 Calculate the 5-sigma limiting depth for a point source in the given conditions. 

94 

95 This is generally not needed, unless the m5 parameters have been updated 

96 or m5 was not previously calculated. 

97 """ 

98 colsAdded = ['m5_simsUtils'] 

99 

100 def __init__(self, airmassCol='airmass', seeingCol='seeingFwhmEff', skybrightnessCol='skyBrightness', 

101 filterCol='filter', exptimeCol='visitExposureTime'): 

102 self.units = ['mag'] 

103 self.colsReq = [airmassCol, seeingCol, skybrightnessCol, filterCol, exptimeCol] 

104 self.airmassCol = airmassCol 

105 self.seeingCol = seeingCol 

106 self.skybrightnessCol = skybrightnessCol 

107 self.filterCol = filterCol 

108 self.exptimeCol = exptimeCol 

109 

110 def _run(self, simData, cols_present=False): 

111 if cols_present: 

112 # Column already present in data; assume it needs updating and recalculate. 

113 return simData 

114 filts = np.unique(simData[self.filterCol]) 

115 for filtername in filts: 

116 infilt = np.where(simData[self.filterCol] == filtername) 

117 simData['m5_simsUtils'][infilt] = m5_flat_sed(filtername, 

118 simData[infilt][self.skybrightnessCol], 

119 simData[infilt][self.seeingCol], 

120 simData[infilt][self.exptimeCol], 

121 simData[infilt][self.airmassCol]) 

122 return simData 

123 

124 

125class NormAirmassStacker(BaseStacker): 

126 """Calculate the normalized airmass for each opsim pointing. 

127 """ 

128 colsAdded = ['normairmass'] 

129 

130 def __init__(self, airmassCol='airmass', decCol='fieldDec', 

131 degrees=True, telescope_lat = -30.2446388): 

132 self.units = ['X / Xmin'] 

133 self.colsReq = [airmassCol, decCol] 

134 self.airmassCol = airmassCol 

135 self.decCol = decCol 

136 self.telescope_lat = telescope_lat 

137 self.degrees = degrees 

138 

139 def _run(self, simData, cols_present=False): 

140 """Calculate new column for normalized airmass.""" 

141 # Run method is required to calculate column. 

142 # Driver runs getColInfo to know what columns are needed from db & which are calculated, 

143 # then gets data from db and then calculates additional columns (via run methods here). 

144 if cols_present: 

145 # Column already present in data; assume it is correct and does not need recalculating. 

146 return simData 

147 dec = simData[self.decCol] 

148 if self.degrees: 

149 dec = np.radians(dec) 

150 min_z_possible = np.abs(dec - np.radians(self.telescope_lat)) 

151 min_airmass_possible = 1./np.cos(min_z_possible) 

152 simData['normairmass'] = simData[self.airmassCol] / min_airmass_possible 

153 return simData 

154 

155 

156class ZenithDistStacker(BaseStacker): 

157 """Calculate the zenith distance for each pointing. 

158 If 'degrees' is True, then assumes altCol is in degrees and returns degrees. 

159 If 'degrees' is False, assumes altCol is in radians and returns radians. 

160 """ 

161 colsAdded = ['zenithDistance'] 

162 

163 def __init__(self, altCol='altitude', degrees=True): 

164 self.altCol = altCol 

165 self.degrees = degrees 

166 if self.degrees: 166 ↛ 169line 166 didn't jump to line 169, because the condition on line 166 was never false

167 self.units = ['degrees'] 

168 else: 

169 self.unit = ['radians'] 

170 self.colsReq = [self.altCol] 

171 

172 def _run(self, simData, cols_present=False): 

173 """Calculate new column for zenith distance.""" 

174 if cols_present: 

175 # Column already present in data; assume it is correct and does not need recalculating. 

176 return simData 

177 if self.degrees: 

178 simData['zenithDistance'] = 90.0 - simData[self.altCol] 

179 else: 

180 simData['zenithDistance'] = np.pi/2.0 - simData[self.altCol] 

181 return simData 

182 

183 

184class ParallaxFactorStacker(BaseStacker): 

185 """Calculate the parallax factors for each opsim pointing. Output parallax factor in arcseconds. 

186 """ 

187 colsAdded = ['ra_pi_amp', 'dec_pi_amp'] 

188 

189 def __init__(self, raCol='fieldRA', decCol='fieldDec', dateCol='observationStartMJD', degrees=True): 

190 self.raCol = raCol 

191 self.decCol = decCol 

192 self.dateCol = dateCol 

193 self.units = ['arcsec', 'arcsec'] 

194 self.colsReq = [raCol, decCol, dateCol] 

195 self.degrees = degrees 

196 

197 def _gnomonic_project_toxy(self, RA1, Dec1, RAcen, Deccen): 

198 """Calculate x/y projection of RA1/Dec1 in system with center at RAcen, Deccenp. 

199 Input radians. 

200 """ 

201 # also used in Global Telescope Network website 

202 cosc = np.sin(Deccen) * np.sin(Dec1) + np.cos(Deccen) * np.cos(Dec1) * np.cos(RA1-RAcen) 

203 x = np.cos(Dec1) * np.sin(RA1-RAcen) / cosc 

204 y = (np.cos(Deccen)*np.sin(Dec1) - np.sin(Deccen)*np.cos(Dec1)*np.cos(RA1-RAcen)) / cosc 

205 return x, y 

206 

207 def _run(self, simData, cols_present=False): 

208 if cols_present: 

209 # Column already present in data; assume it is correct and does not need recalculating. 

210 return simData 

211 ra_pi_amp = np.zeros(np.size(simData), dtype=[('ra_pi_amp', 'float')]) 

212 dec_pi_amp = np.zeros(np.size(simData), dtype=[('dec_pi_amp', 'float')]) 

213 ra_geo1 = np.zeros(np.size(simData), dtype='float') 

214 dec_geo1 = np.zeros(np.size(simData), dtype='float') 

215 ra_geo = np.zeros(np.size(simData), dtype='float') 

216 dec_geo = np.zeros(np.size(simData), dtype='float') 

217 ra = simData[self.raCol] 

218 dec = simData[self.decCol] 

219 if self.degrees: 

220 ra = np.radians(ra) 

221 dec = np.radians(dec) 

222 

223 for i, ack in enumerate(simData): 

224 mtoa_params = palpy.mappa(2000., simData[self.dateCol][i]) 

225 # Object with a 1 arcsec parallax 

226 ra_geo1[i], dec_geo1[i] = palpy.mapqk(ra[i], dec[i], 

227 0., 0., 1., 0., mtoa_params) 

228 # Object with no parallax 

229 ra_geo[i], dec_geo[i] = palpy.mapqk(ra[i], dec[i], 

230 0., 0., 0., 0., mtoa_params) 

231 x_geo1, y_geo1 = self._gnomonic_project_toxy(ra_geo1, dec_geo1, 

232 ra, dec) 

233 x_geo, y_geo = self._gnomonic_project_toxy(ra_geo, dec_geo, ra, dec) 

234 # Return ra_pi_amp and dec_pi_amp in arcseconds. 

235 ra_pi_amp[:] = np.degrees(x_geo1-x_geo)*3600. 

236 dec_pi_amp[:] = np.degrees(y_geo1-y_geo)*3600. 

237 simData['ra_pi_amp'] = ra_pi_amp 

238 simData['dec_pi_amp'] = dec_pi_amp 

239 return simData 

240 

241 

242class DcrStacker(BaseStacker): 

243 """Calculate the RA,Dec offset expected for an object due to differential chromatic refraction. 

244 

245 For DCR calculation, we also need zenithDistance, HA, and PA -- but these will be explicitly 

246 handled within this stacker so that setup is consistent and they run in order. If those values 

247 have already been calculated elsewhere, they will not be overwritten. 

248 

249 Parameters 

250 ---------- 

251 filterCol : str 

252 The name of the column with filter names. Default 'fitler'. 

253 altCol : str 

254 Name of the column with altitude info. Default 'altitude'. 

255 raCol : str 

256 Name of the column with RA. Default 'fieldRA'. 

257 decCol : str 

258 Name of the column with Dec. Default 'fieldDec'. 

259 lstCol : str 

260 Name of the column with local sidereal time. Default 'observationStartLST'. 

261 site : str or lsst.sims.utils.Site 

262 Name of the observory or a lsst.sims.utils.Site object. Default 'LSST'. 

263 mjdCol : str 

264 Name of column with modified julian date. Default 'observationStartMJD' 

265 dcr_magnitudes : dict 

266 Magitude of the DCR offset for each filter at altitude/zenith distance of 45 degrees. 

267 Defaults u=0.07, g=0.07, r=0.50, i=0.045, z=0.042, y=0.04 (all arcseconds). 

268 

269 Returns 

270 ------- 

271 numpy.array 

272 Returns array with additional columns 'ra_dcr_amp' and 'dec_dcr_amp' with the DCR offsets 

273 for each observation. Also runs ZenithDistStacker and ParallacticAngleStacker. 

274 """ 

275 colsAdded = ['ra_dcr_amp', 'dec_dcr_amp'] # zenithDist, HA, PA 

276 

277 def __init__(self, filterCol='filter', altCol='altitude', degrees=True, 

278 raCol='fieldRA', decCol='fieldDec', lstCol='observationStartLST', 

279 site='LSST', mjdCol='observationStartMJD', 

280 dcr_magnitudes=None): 

281 self.units = ['arcsec', 'arcsec'] 

282 if dcr_magnitudes is None: 282 ↛ 286line 282 didn't jump to line 286, because the condition on line 282 was never false

283 # DCR amplitudes are in arcseconds. 

284 self.dcr_magnitudes = {'u': 0.07, 'g': 0.07, 'r': 0.050, 'i': 0.045, 'z': 0.042, 'y': 0.04} 

285 else: 

286 self.dcr_magnitudes = dcr_magnitudes 

287 self.zdCol = 'zenithDistance' 

288 self.paCol = 'PA' 

289 self.filterCol = filterCol 

290 self.raCol = raCol 

291 self.decCol = decCol 

292 self.degrees = degrees 

293 self.colsReq = [filterCol, raCol, decCol, altCol, lstCol] 

294 # 'zenithDist', 'PA', 'HA' are additional columns required, coming from other stackers which must 

295 # also be configured -- so we handle this explicitly here. 

296 self.zstacker = ZenithDistStacker(altCol=altCol, degrees=self.degrees) 

297 self.pastacker = ParallacticAngleStacker(raCol=raCol, decCol=decCol, mjdCol=mjdCol, 

298 degrees=self.degrees, 

299 lstCol=lstCol, site=site) 

300 # Note that RA/Dec could be coming from a dither stacker! 

301 # But we will assume that coord stackers will be handled separately. 

302 

303 

304 def _run(self, simData, cols_present=False): 

305 if cols_present: 

306 # Column already present in data; assume it is correct and does not need recalculating. 

307 return simData 

308 # Need to make sure the Zenith stacker gets run first 

309 # Call _run method because already added these columns due to 'colsAdded' line. 

310 simData = self.zstacker.run(simData) 

311 simData = self.pastacker.run(simData) 

312 if self.degrees: 

313 zenithTan = np.tan(np.radians(simData[self.zdCol])) 

314 parallacticAngle = np.radians(simData[self.paCol]) 

315 else: 

316 zenithTan = np.tan(simData[self.zdCol]) 

317 parallacticAngle = simData[self.paCol] 

318 dcr_in_ra = zenithTan * np.sin(parallacticAngle) 

319 dcr_in_dec = zenithTan * np.cos(parallacticAngle) 

320 for filtername in np.unique(simData[self.filterCol]): 

321 fmatch = np.where(simData[self.filterCol] == filtername) 

322 dcr_in_ra[fmatch] = self.dcr_magnitudes[filtername] * dcr_in_ra[fmatch] 

323 dcr_in_dec[fmatch] = self.dcr_magnitudes[filtername] * dcr_in_dec[fmatch] 

324 simData['ra_dcr_amp'] = dcr_in_ra 

325 simData['dec_dcr_amp'] = dcr_in_dec 

326 return simData 

327 

328 

329class HourAngleStacker(BaseStacker): 

330 """Add the Hour Angle for each observation. 

331 Always in HOURS. 

332 """ 

333 colsAdded = ['HA'] 

334 

335 def __init__(self, lstCol='observationStartLST', raCol='fieldRA', degrees=True): 

336 self.units = ['Hours'] 

337 self.colsReq = [lstCol, raCol] 

338 self.lstCol = lstCol 

339 self.raCol = raCol 

340 self.degrees = degrees 

341 

342 def _run(self, simData, cols_present=False): 

343 """HA = LST - RA """ 

344 if cols_present: 

345 # Column already present in data; assume it is correct and does not need recalculating. 

346 return simData 

347 if len(simData) == 0: 

348 return simData 

349 if self.degrees: 

350 ra = np.radians(simData[self.raCol]) 

351 lst = np.radians(simData[self.lstCol]) 

352 else: 

353 ra = simData[self.raCol] 

354 lst = simData[self.lstCol] 

355 # Check that LST is reasonable 

356 if (np.min(lst) < 0) | (np.max(lst) > 2.*np.pi): 

357 warnings.warn('LST values are not between 0 and 2 pi') 

358 # Check that RA is reasonable 

359 if (np.min(ra) < 0) | (np.max(ra) > 2.*np.pi): 

360 warnings.warn('RA values are not between 0 and 2 pi') 

361 ha = lst - ra 

362 # Wrap the results so HA between -pi and pi 

363 ha = np.where(ha < -np.pi, ha + 2. * np.pi, ha) 

364 ha = np.where(ha > np.pi, ha - 2. * np.pi, ha) 

365 # Convert radians to hours 

366 simData['HA'] = ha*12/np.pi 

367 return simData 

368 

369 

370class ParallacticAngleStacker(BaseStacker): 

371 """Add the parallactic angle to each visit. 

372 If 'degrees' is True, this will be in degrees (as are all other angles). If False, then in radians. 

373 """ 

374 colsAdded = ['PA'] 

375 

376 def __init__(self, raCol='fieldRA', decCol='fieldDec', degrees=True, mjdCol='observationStartMJD', 

377 lstCol='observationStartLST', site='LSST'): 

378 

379 self.lstCol = lstCol 

380 self.raCol = raCol 

381 self.decCol = decCol 

382 self.degrees = degrees 

383 self.mjdCol = mjdCol 

384 self.site = Site(name=site) 

385 self.units = ['radians'] 

386 self.colsReq = [self.raCol, self.decCol, self.mjdCol, self.lstCol] 

387 self.haStacker = HourAngleStacker(lstCol=lstCol, raCol=raCol, degrees=self.degrees) 

388 

389 def _run(self, simData, cols_present=False): 

390 # Equation from: 

391 # http://www.gb.nrao.edu/~rcreager/GBTMetrology/140ft/l0058/gbtmemo52/memo52.html 

392 # or 

393 # http://www.gb.nrao.edu/GBT/DA/gbtidl/release2pt9/contrib/contrib/parangle.pro 

394 if cols_present: 

395 # Column already present in data; assume it is correct and does not need recalculating. 

396 return simData 

397 # Using the run method (not _run) means that if HA is present, it will not be recalculated. 

398 simData = self.haStacker.run(simData) 

399 if self.degrees: 

400 dec = np.radians(simData[self.decCol]) 

401 else: 

402 dec = simData[self.decCol] 

403 simData['PA'] = np.arctan2(np.sin(simData['HA']*np.pi/12.), (np.cos(dec) * 

404 np.tan(self.site.latitude_rad) - np.sin(dec) * 

405 np.cos(simData['HA']*np.pi/12.))) 

406 if self.degrees: 

407 simData['PA'] = np.degrees(simData['PA']) 

408 return simData 

409 

410 

411class FilterColorStacker(BaseStacker): 

412 """Translate filters ('u', 'g', 'r' ..) into RGB tuples. 

413 

414 This is useful for making movies if you want to make the pointing have a related color-tuple for a plot. 

415 """ 

416 colsAdded = ['rRGB', 'gRGB', 'bRGB'] 

417 

418 def __init__(self, filterCol='filter'): 

419 self.filter_rgb_map = {'u': (0, 0, 1), # dark blue 

420 'g': (0, 1, 1), # cyan 

421 'r': (0, 1, 0), # green 

422 'i': (1, 0.5, 0.3), # orange 

423 'z': (1, 0, 0), # red 

424 'y': (1, 0, 1)} # magenta 

425 self.filterCol = filterCol 

426 # self.units used for plot labels 

427 self.units = ['rChan', 'gChan', 'bChan'] 

428 # Values required for framework operation: this specifies the data columns required from the database. 

429 self.colsReq = [self.filterCol] 

430 

431 def _run(self, simData, cols_present=False): 

432 # Translate filter names into numbers. 

433 if cols_present: 

434 # Column already present in data; assume it is correct and does not need recalculating. 

435 return simData 

436 filtersUsed = np.unique(simData[self.filterCol]) 

437 for f in filtersUsed: 

438 if f not in self.filter_rgb_map: 

439 raise IndexError('Filter %s not in filter_rgb_map' % (f)) 

440 match = np.where(simData[self.filterCol] == f)[0] 

441 simData['rRGB'][match] = self.filter_rgb_map[f][0] 

442 simData['gRGB'][match] = self.filter_rgb_map[f][1] 

443 simData['bRGB'][match] = self.filter_rgb_map[f][2] 

444 return simData 

445 

446 

447class OpSimFieldStacker(BaseStacker): 

448 """Add the fieldId of the closest OpSim field for each RA/Dec pointing. 

449 

450 Parameters 

451 ---------- 

452 raCol : str, opt 

453 Name of the RA column. Default fieldRA. 

454 decCol : str, opt 

455 Name of the Dec column. Default fieldDec. 

456 

457 """ 

458 colsAdded = ['opsimFieldId'] 

459 

460 def __init__(self, raCol='fieldRA', decCol='fieldDec', degrees=True): 

461 self.colsReq = [raCol, decCol] 

462 self.units = ['#'] 

463 self.raCol = raCol 

464 self.decCol = decCol 

465 self.degrees = degrees 

466 fields_db = FieldsDatabase() 

467 # Returned RA/Dec coordinates in degrees 

468 fieldid, ra, dec = fields_db.get_id_ra_dec_arrays("select * from Field;") 

469 asort = np.argsort(fieldid) 

470 self.tree = _buildTree(np.radians(ra[asort]), 

471 np.radians(dec[asort])) 

472 

473 def _run(self, simData, cols_present=False): 

474 if cols_present: 

475 # Column already present in data; assume it is correct and does not need recalculating. 

476 return simData 

477 

478 if self.degrees: 

479 coord_x, coord_y, coord_z = xyz_from_ra_dec(simData[self.raCol], 

480 simData[self.decCol]) 

481 field_ids = self.tree.query_ball_point(list(zip(coord_x, coord_y, coord_z)), 

482 xyz_angular_radius()) 

483 

484 else: 

485 # use _xyz private method (sending radians) 

486 coord_x, coord_y, coord_z = _xyz_from_ra_dec(simData[self.raCol], 

487 simData[self.decCol]) 

488 field_ids = self.tree.query_ball_point(list(zip(coord_x, coord_y, coord_z)), 

489 xyz_angular_radius()) 

490 

491 simData['opsimFieldId'] = np.array([ids[0] for ids in field_ids]) + 1 

492 return simData