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 'SeasonStacker', 'DcrStacker', 'FiveSigmaStacker', 'OpSimFieldStacker'] 

12 

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

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

15# Season stacker by Phil Marshall (dr.phil.marshall@gmail.com), 

16# modified by Humna Awan (humna.awan@rutgers.edu) 

17 

18 

19class FiveSigmaStacker(BaseStacker): 

20 """ 

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

22 

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

24 or m5 was not previously calculated. 

25 """ 

26 colsAdded = ['m5_simsUtils'] 

27 

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

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

30 self.units = ['mag'] 

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

32 self.airmassCol = airmassCol 

33 self.seeingCol = seeingCol 

34 self.skybrightnessCol = skybrightnessCol 

35 self.filterCol = filterCol 

36 self.exptimeCol = exptimeCol 

37 

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

39 if cols_present: 

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

41 return simData 

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

43 for filtername in filts: 

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

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

46 simData[infilt][self.skybrightnessCol], 

47 simData[infilt][self.seeingCol], 

48 simData[infilt][self.exptimeCol], 

49 simData[infilt][self.airmassCol]) 

50 return simData 

51 

52 

53class NormAirmassStacker(BaseStacker): 

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

55 """ 

56 colsAdded = ['normairmass'] 

57 

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

59 degrees=True, telescope_lat = -30.2446388): 

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

61 self.colsReq = [airmassCol, decCol] 

62 self.airmassCol = airmassCol 

63 self.decCol = decCol 

64 self.telescope_lat = telescope_lat 

65 self.degrees = degrees 

66 

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

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

69 # Run method is required to calculate column. 

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

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

72 if cols_present: 

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

74 return simData 

75 dec = simData[self.decCol] 

76 if self.degrees: 

77 dec = np.radians(dec) 

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

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

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

81 return simData 

82 

83 

84class ZenithDistStacker(BaseStacker): 

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

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

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

88 """ 

89 colsAdded = ['zenithDistance'] 

90 

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

92 self.altCol = altCol 

93 self.degrees = degrees 

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

95 self.units = ['degrees'] 

96 else: 

97 self.unit = ['radians'] 

98 self.colsReq = [self.altCol] 

99 

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

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

102 if cols_present: 

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

104 return simData 

105 if self.degrees: 

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

107 else: 

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

109 return simData 

110 

111 

112class ParallaxFactorStacker(BaseStacker): 

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

114 """ 

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

116 

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

118 self.raCol = raCol 

119 self.decCol = decCol 

120 self.dateCol = dateCol 

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

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

123 self.degrees = degrees 

124 

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

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

127 Input radians. 

128 """ 

129 # also used in Global Telescope Network website 

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

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

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

133 return x, y 

134 

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

136 if cols_present: 

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

138 return simData 

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

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

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

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

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

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

145 ra = simData[self.raCol] 

146 dec = simData[self.decCol] 

147 if self.degrees: 

148 ra = np.radians(ra) 

149 dec = np.radians(dec) 

150 

151 for i, ack in enumerate(simData): 

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

153 # Object with a 1 arcsec parallax 

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

155 0., 0., 1., 0., mtoa_params) 

156 # Object with no parallax 

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

158 0., 0., 0., 0., mtoa_params) 

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

160 ra, dec) 

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

162 # Return ra_pi_amp and dec_pi_amp in arcseconds. 

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

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

165 simData['ra_pi_amp'] = ra_pi_amp 

166 simData['dec_pi_amp'] = dec_pi_amp 

167 return simData 

168 

169 

170class DcrStacker(BaseStacker): 

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

172 

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

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

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

176 

177 Parameters 

178 ---------- 

179 filterCol : str 

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

181 altCol : str 

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

183 raCol : str 

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

185 decCol : str 

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

187 lstCol : str 

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

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

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

191 mjdCol : str 

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

193 dcr_magnitudes : dict 

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

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

196 

197 Returns 

198 ------- 

199 numpy.array 

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

201 for each observation. Also runs ZenithDistStacker and ParallacticAngleStacker. 

202 """ 

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

204 

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

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

207 site='LSST', mjdCol='observationStartMJD', 

208 dcr_magnitudes=None): 

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

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

211 # DCR amplitudes are in arcseconds. 

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

213 else: 

214 self.dcr_magnitudes = dcr_magnitudes 

215 self.zdCol = 'zenithDistance' 

216 self.paCol = 'PA' 

217 self.filterCol = filterCol 

218 self.raCol = raCol 

219 self.decCol = decCol 

220 self.degrees = degrees 

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

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

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

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

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

226 degrees=self.degrees, 

227 lstCol=lstCol, site=site) 

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

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

230 

231 

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

233 if cols_present: 

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

235 return simData 

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

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

238 simData = self.zstacker.run(simData) 

239 simData = self.pastacker.run(simData) 

240 if self.degrees: 

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

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

243 else: 

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

245 parallacticAngle = simData[self.paCol] 

246 dcr_in_ra = zenithTan * np.sin(parallacticAngle) 

247 dcr_in_dec = zenithTan * np.cos(parallacticAngle) 

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

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

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

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

252 simData['ra_dcr_amp'] = dcr_in_ra 

253 simData['dec_dcr_amp'] = dcr_in_dec 

254 return simData 

255 

256 

257class HourAngleStacker(BaseStacker): 

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

259 Always in HOURS. 

260 """ 

261 colsAdded = ['HA'] 

262 

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

264 self.units = ['Hours'] 

265 self.colsReq = [lstCol, raCol] 

266 self.lstCol = lstCol 

267 self.raCol = raCol 

268 self.degrees = degrees 

269 

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

271 """HA = LST - RA """ 

272 if cols_present: 

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

274 return simData 

275 if len(simData) == 0: 

276 return simData 

277 if self.degrees: 

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

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

280 else: 

281 ra = simData[self.raCol] 

282 lst = simData[self.lstCol] 

283 # Check that LST is reasonable 

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

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

286 # Check that RA is reasonable 

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

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

289 ha = lst - ra 

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

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

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

293 # Convert radians to hours 

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

295 return simData 

296 

297 

298class ParallacticAngleStacker(BaseStacker): 

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

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

301 """ 

302 colsAdded = ['PA'] 

303 

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

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

306 

307 self.lstCol = lstCol 

308 self.raCol = raCol 

309 self.decCol = decCol 

310 self.degrees = degrees 

311 self.mjdCol = mjdCol 

312 self.site = Site(name=site) 

313 self.units = ['radians'] 

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

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

316 

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

318 # Equation from: 

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

320 # or 

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

322 if cols_present: 

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

324 return simData 

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

326 simData = self.haStacker.run(simData) 

327 if self.degrees: 

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

329 else: 

330 dec = simData[self.decCol] 

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

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

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

334 if self.degrees: 

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

336 return simData 

337 

338 

339class FilterColorStacker(BaseStacker): 

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

341 """ 

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

343 

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

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

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

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

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

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

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

351 self.filterCol = filterCol 

352 # self.units used for plot labels 

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

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

355 self.colsReq = [self.filterCol] 

356 

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

358 # Translate filter names into numbers. 

359 if cols_present: 

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

361 return simData 

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

363 for f in filtersUsed: 

364 if f not in self.filter_rgb_map: 

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

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

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

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

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

370 return simData 

371 

372 

373class SeasonStacker(BaseStacker): 

374 """Add an integer label to show which season a given visit is in. 

375 

376 The season only depends on the RA of the object: we compute the MJD 

377 when each object is on the meridian at midnight, and subtract 6 

378 months to get the start date of each season. 

379 The season index range is 0-10. 

380 Must wrap 0th and 10th to get a total of 10 seasons. 

381 """ 

382 colsAdded = ['year', 'season'] 

383 

384 def __init__(self, mjdCol='observationStartMJD', RACol='fieldRA', degrees=True): 

385 # Names of columns we need from database. 

386 self.colsReq = [mjdCol, RACol] 

387 # List of units for our new columns. 

388 self.units = ['yr', ''] 

389 # And save the column names. 

390 self.mjdCol = mjdCol 

391 self.RACol = RACol 

392 self.degrees = degrees 

393 

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

395 if cols_present: 

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

397 return simData 

398 # Define year number: (note that opsim defines "years" in flat 365 days). 

399 year = np.floor((simData[self.mjdCol] - simData[self.mjdCol][0]) / 365) 

400 # Convert RA to Hours 

401 if self.degrees: 

402 objRA = simData[self.RACol]/15.0 

403 else: 

404 objRA = np.degrees(simData[self.RACol])/15.0 

405 # objRA=0 on autumnal equinox. 

406 # autumnal equinox 2014 happened on Sept 23 --> Equinox MJD 

407 Equinox = 2456923.5 - 2400000.5 

408 # Use 365.25 for the length of a year here, because we're dealing with real seasons. 

409 daysSinceEquinox = 0.5*objRA*(365.25/12.0) # 0.5 to go from RA to month; 365.25/12.0 months to days 

410 firstSeasonBegan = Equinox + daysSinceEquinox - 0.5*365.25 # in MJD 

411 # Now we can compute the number of years since the first season 

412 # began, and so assign a global integer season number: 

413 globalSeason = np.floor((simData[self.mjdCol] - firstSeasonBegan)/365.25) 

414 # Subtract off season number of first observation: 

415 season = globalSeason - np.min(globalSeason) 

416 simData['year'] = year 

417 simData['season'] = season 

418 return simData 

419 

420 

421class OpSimFieldStacker(BaseStacker): 

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

423 

424 Parameters 

425 ---------- 

426 raCol : str, opt 

427 Name of the RA column. Default fieldRA. 

428 decCol : str, opt 

429 Name of the Dec column. Default fieldDec. 

430 

431 """ 

432 colsAdded = ['opsimFieldId'] 

433 

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

435 self.colsReq = [raCol, decCol] 

436 self.units = ['#'] 

437 self.raCol = raCol 

438 self.decCol = decCol 

439 self.degrees = degrees 

440 fields_db = FieldsDatabase() 

441 # Returned RA/Dec coordinates in degrees 

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

443 asort = np.argsort(fieldid) 

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

445 np.radians(dec[asort])) 

446 

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

448 if cols_present: 

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

450 return simData 

451 

452 if self.degrees: 

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

454 simData[self.decCol]) 

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

456 xyz_angular_radius()) 

457 

458 else: 

459 # use _xyz private method (sending radians) 

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

461 simData[self.decCol]) 

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

463 xyz_angular_radius()) 

464 

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

466 return simData