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

1from __future__ import absolute_import 

2from builtins import object 

3import numpy as np 

4import healpy as hp 

5from lsst.sims.featureScheduler import utils 

6from lsst.sims.utils import m5_flat_sed 

7 

8 

9__all__ = ['BaseFeature', 'BaseSurveyFeature', 'N_obs_count', 'N_obs_survey', 

10 'Last_observation', 'LastSequence_observation', 'LastFilterChange', 

11 'N_observations', 'Coadded_depth', 'Last_observed', 'N_obs_night', 'Pair_in_night', 

12 'Rotator_angle', 'N_observations_season', 'N_obs_count_season', 'N_observations_current_season', 

13 'Last_N_obs_times'] 

14 

15 

16class BaseFeature(object): 

17 """ 

18 Base class for features. 

19 """ 

20 def __init__(self, **kwargs): 

21 # self.feature should be a float, bool, or healpix size numpy array, or numpy masked array 

22 self.feature = None 

23 

24 def __call__(self): 

25 return self.feature 

26 

27 

28class BaseSurveyFeature(object): 

29 """ 

30 Feature that tracks progreess of the survey. Takes observations and updates self.feature 

31 """ 

32 def add_observation(self, observation, indx=None, **kwargs): 

33 """ 

34 Parameters 

35 ---------- 

36 obsevation : dict-like 

37 Object that contains the information about the observation (ra, dec, filter, mjd, etc) 

38 indx : ints (None) 

39 The healpixel indices that the observation overlaps. 

40 """ 

41 raise NotImplementedError 

42 

43 

44class N_obs_count(BaseSurveyFeature): 

45 """Count the number of observations. Total number, not tracked over sky 

46 

47 Parameters 

48 ---------- 

49 filtername : str (None) 

50 The filter to count (if None, all filters counted) 

51 """ 

52 def __init__(self, filtername=None, tag=None): 

53 self.feature = 0 

54 self.filtername = filtername 

55 self.tag = tag 

56 

57 def add_observation(self, observation, indx=None): 

58 

59 if (self.filtername is None) and (self.tag is None): 

60 # Track all observations 

61 self.feature += 1 

62 elif (self.filtername is not None) and (self.tag is None) and (observation['filter'][0] in self.filtername): 

63 # Track all observations on a specified filter 

64 self.feature += 1 

65 elif (self.filtername is None) and (self.tag is not None) and (observation['tag'][0] in self.tag): 

66 # Track all observations on a specified tag 

67 self.feature += 1 

68 elif ((self.filtername is None) and (self.tag is not None) and 

69 # Track all observations on a specified filter on a specified tag 

70 (observation['filter'][0] in self.filtername) and (observation['tag'][0] in self.tag)): 

71 self.feature += 1 

72 

73 

74class N_obs_count_season(BaseSurveyFeature): 

75 """Count the number of observations. 

76 

77 Parameters 

78 ---------- 

79 filtername : str (None) 

80 The filter to count (if None, all filters counted) 

81 """ 

82 def __init__(self, season, nside=None, filtername=None, tag=None, 

83 season_modulo=2, offset=None, max_season=None, season_length=365.25): 

84 self.feature = 0 

85 self.filtername = filtername 

86 self.tag = tag 

87 self.season = season 

88 self.season_modulo = season_modulo 

89 if offset is None: 

90 self.offset = np.zeros(hp.nside2npix(nside), dtype=int) 

91 else: 

92 self.offset = offset 

93 self.max_season = max_season 

94 self.season_length = season_length 

95 

96 def add_observation(self, observation, indx=None): 

97 

98 season = utils.season_calc(observation['night'], modulo=self.season_modulo, 

99 offset=self.offset[indx], max_season=self.max_season, 

100 season_length=self.season_length) 

101 if self.season in season: 

102 if (self.filtername is None) and (self.tag is None): 

103 # Track all observations 

104 self.feature += 1 

105 elif (self.filtername is not None) and (self.tag is None) and (observation['filter'][0] in self.filtername): 

106 # Track all observations on a specified filter 

107 self.feature += 1 

108 elif (self.filtername is None) and (self.tag is not None) and (observation['tag'][0] in self.tag): 

109 # Track all observations on a specified tag 

110 self.feature += 1 

111 elif ((self.filtername is None) and (self.tag is not None) and 

112 # Track all observations on a specified filter on a specified tag 

113 (observation['filter'][0] in self.filtername) and (observation['tag'][0] in self.tag)): 

114 self.feature += 1 

115 

116 

117class N_obs_survey(BaseSurveyFeature): 

118 """Count the number of observations. 

119 

120 Parameters 

121 ---------- 

122 note : str (None) 

123 Only count observations that have str in their note field 

124 """ 

125 def __init__(self, note=None): 

126 self.feature = 0 

127 self.note = note 

128 

129 def add_observation(self, observation, indx=None): 

130 # Track all observations 

131 if self.note is None: 

132 self.feature += 1 

133 else: 

134 if self.note == observation['note']: 

135 self.feature += 1 

136 

137 

138class Last_observation(BaseSurveyFeature): 

139 """Track the last observation. Useful if you want to see when the 

140 last time a survey took an observation. 

141 

142 Parameters 

143 ---------- 

144 survey_name : str (None) 

145 Only records if the survey name matches (or survey_name set to None) 

146 """ 

147 def __init__(self, survey_name=None): 

148 self.survey_name = survey_name 

149 # Start out with an empty observation 

150 self.feature = utils.empty_observation() 

151 

152 def add_observation(self, observation, indx=None): 

153 if self.survey_name is not None: 

154 if self.survey_name in observation['note']: 

155 self.feature = observation 

156 else: 

157 self.feature = observation 

158 

159 

160class LastSequence_observation(BaseSurveyFeature): 

161 """When was the last observation 

162 """ 

163 def __init__(self, sequence_ids=''): 

164 self.sequence_ids = sequence_ids # The ids of all sequence observations... 

165 # Start out with an empty observation 

166 self.feature = utils.empty_observation() 

167 

168 def add_observation(self, observation, indx=None): 

169 if observation['survey_id'] in self.sequence_ids: 

170 self.feature = observation 

171 

172 

173class LastFilterChange(BaseSurveyFeature): 

174 """Record when the filter last changed. 

175 """ 

176 def __init__(self): 

177 self.feature = {'mjd': 0., 

178 'previous_filter': None, 

179 'current_filter': None} 

180 

181 def add_observation(self, observation, indx=None): 

182 if self.feature['current_filter'] is None: 

183 self.feature['mjd'] = observation['mjd'][0] 

184 self.feature['previous_filter'] = None 

185 self.feature['current_filter'] = observation['filter'][0] 

186 elif observation['filter'][0] != self.feature['current_filter']: 

187 self.feature['mjd'] = observation['mjd'][0] 

188 self.feature['previous_filter'] = self.feature['current_filter'] 

189 self.feature['current_filter'] = observation['filter'][0] 

190 

191 

192class N_observations(BaseSurveyFeature): 

193 """ 

194 Track the number of observations that have been made across the sky. 

195 

196 Parameters 

197 ---------- 

198 filtername : str ('r') 

199 String or list that has all the filters that can count. 

200 nside : int (32) 

201 The nside of the healpixel map to use 

202 mask_indx : list of ints (None) 

203 List of healpixel indices to mask and interpolate over 

204 

205 """ 

206 def __init__(self, filtername=None, nside=None, mask_indx=None, survey_name=None): 

207 if nside is None: 

208 nside = utils.set_default_nside() 

209 

210 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

211 self.filtername = filtername 

212 self.mask_indx = mask_indx 

213 self.survey_name = survey_name 

214 

215 def add_observation(self, observation, indx=None): 

216 """ 

217 Parameters 

218 ---------- 

219 indx : ints 

220 The indices of the healpixel map that have been observed by observation 

221 """ 

222 

223 if self.filtername is None or observation['filter'][0] in self.filtername: 

224 if self.survey_name is None or observation['note'] in self.survey_name: 

225 self.feature[indx] += 1 

226 

227 if self.mask_indx is not None: 

228 overlap = np.intersect1d(indx, self.mask_indx) 

229 if overlap.size > 0: 

230 # interpolate over those pixels that are DD fields. 

231 # XXX. Do I need to kdtree this? Maybe make a dict on init 

232 # to lookup the N closest non-masked pixels, then do weighted average. 

233 pass 

234 

235 

236class N_observations_season(BaseSurveyFeature): 

237 """ 

238 Track the number of observations that have been made across sky 

239 

240 Parameters 

241 ---------- 

242 season : int 

243 Only count observations in this season (year). 

244 filtername : str ('r') 

245 String or list that has all the filters that can count. 

246 nside : int (32) 

247 The nside of the healpixel map to use 

248 offset : int (0) 

249 The offset to use when computing the season (days) 

250 modulo : int (None) 

251 How to mod the years when computing season 

252 

253 """ 

254 def __init__(self, season, filtername=None, nside=None, offset=0, modulo=None, 

255 max_season=None, season_length=365.25): 

256 if offset is None: 

257 offset = np.zeros(hp.nside2npix(nside), dtype=int) 

258 if nside is None: 

259 nside = utils.set_default_nside() 

260 

261 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

262 self.filtername = filtername 

263 self.offset = offset 

264 self.modulo = modulo 

265 self.season = season 

266 self.max_season = max_season 

267 self.season_length = season_length 

268 

269 def add_observation(self, observation, indx=None): 

270 """ 

271 Parameters 

272 ---------- 

273 indx : ints 

274 The indices of the healpixel map that have been observed by observation 

275 """ 

276 

277 observation_season = utils.season_calc(observation['night'], offset=self.offset[indx], 

278 modulo=self.modulo, max_season=self.max_season, 

279 season_length=self.season_length) 

280 if self.season in observation_season: 

281 if self.filtername is None or observation['filter'][0] in self.filtername: 

282 self.feature[indx] += 1 

283 

284 

285class Last_N_obs_times(BaseSurveyFeature): 

286 """Record the last three observations for each healpixel 

287 """ 

288 def __init__(self, filtername=None, n_obs=3, nside=None): 

289 self.filtername = filtername 

290 self.n_obs = n_obs 

291 if nside is None: 

292 nside = utils.set_default_nside() 

293 self.feature = np.zeros((n_obs, hp.nside2npix(nside)), dtype=float) 

294 

295 def add_observation(self, observation, indx=None): 

296 

297 if self.filtername is None or observation['filter'][0] in self.filtername: 

298 self.feature[0:-1, indx] = self.feature[1:, indx] 

299 self.feature[-1, indx] = observation['mjd'] 

300 

301 

302class N_observations_current_season(BaseSurveyFeature): 

303 """Track how many observations have been taken in the current season 

304 XXX--experimental 

305 """ 

306 def __init__(self, filtername=None, nside=None, offset=0, season_length=365.25): 

307 self.filtername = filtername 

308 if nside is None: 

309 nside = utils.set_default_nside() 

310 if offset is None: 

311 offset = np.zeros(hp.nside2npix(nside), dtype=int) 

312 self.offset = offset 

313 self.season_length = season_length 

314 self.season_map = utils.season_calc(0., offset=self.offset, season_length=season_length) 

315 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

316 

317 def add_observation(self, observation, indx=None): 

318 current_season = utils.season_calc(observation['night'], offset=self.offset, 

319 season_length=self.season_length) 

320 # If the season has changed anywhere, set that count to zero 

321 new_season = np.where((self.season_map - current_season) != 0) 

322 self.feature[new_season] = 0 

323 self.season_map = current_season 

324 

325 if self.filtername is None or observation['filter'][0] in self.filtername: 

326 self.feature[indx] += 1 

327 

328 

329class Coadded_depth(BaseSurveyFeature): 

330 """ 

331 Track the co-added depth that has been reached accross the sky 

332 """ 

333 def __init__(self, filtername='r', nside=None): 

334 if nside is None: 

335 nside = utils.set_default_nside() 

336 self.filtername = filtername 

337 # Starting at limiting mag of zero should be fine. 

338 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

339 

340 def add_observation(self, observation, indx=None): 

341 

342 if observation['filter'][0] == self.filtername: 

343 m5 = m5_flat_sed(observation['filter'][0], observation['skybrightness'][0], 

344 observation['FWHMeff'][0], observation['exptime'][0], 

345 observation['airmass'][0]) 

346 self.feature[indx] = 1.25 * np.log10(10.**(0.8*self.feature[indx]) + 10.**(0.8*m5)) 

347 

348 

349class Last_observed(BaseSurveyFeature): 

350 """ 

351 Track when a pixel was last observed. Assumes observations are added in chronological 

352 order. 

353 """ 

354 def __init__(self, filtername='r', nside=None): 

355 if nside is None: 

356 nside = utils.set_default_nside() 

357 

358 self.filtername = filtername 

359 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

360 

361 def add_observation(self, observation, indx=None): 

362 if self.filtername is None: 

363 self.feature[indx] = observation['mjd'] 

364 elif observation['filter'][0] in self.filtername: 

365 self.feature[indx] = observation['mjd'] 

366 

367 

368class N_obs_night(BaseSurveyFeature): 

369 """ 

370 Track how many times something has been observed in a night 

371 (Note, even if there are two, it might not be a good pair.) 

372 

373 Parameters 

374 ---------- 

375 filtername : string ('r') 

376 Filter to track. 

377 nside : int (32) 

378 Scale of the healpix map 

379 

380 """ 

381 def __init__(self, filtername='r', nside=None): 

382 if nside is None: 

383 nside = utils.set_default_nside() 

384 

385 self.filtername = filtername 

386 self.feature = np.zeros(hp.nside2npix(nside), dtype=int) 

387 self.night = None 

388 

389 def add_observation(self, observation, indx=None): 

390 if observation['night'][0] != self.night: 

391 self.feature *= 0 

392 self.night = observation['night'][0] 

393 if observation['filter'][0] in self.filtername: 

394 self.feature[indx] += 1 

395 

396 

397class Pair_in_night(BaseSurveyFeature): 

398 """ 

399 Track how many pairs have been observed within a night 

400 

401 Parameters 

402 ---------- 

403 gap_min : float (25.) 

404 The minimum time gap to consider a successful pair in minutes 

405 gap_max : float (45.) 

406 The maximum time gap to consider a successful pair (minutes) 

407 """ 

408 def __init__(self, filtername='r', nside=None, gap_min=25., gap_max=45.): 

409 if nside is None: 

410 nside = utils.set_default_nside() 

411 

412 self.filtername = filtername 

413 self.feature = np.zeros(hp.nside2npix(nside), dtype=float) 

414 self.indx = np.arange(self.feature.size) 

415 self.last_observed = Last_observed(filtername=filtername) 

416 self.gap_min = gap_min / (24.*60) # Days 

417 self.gap_max = gap_max / (24.*60) # Days 

418 self.night = 0 

419 # Need to keep a full record of times and healpixels observed in a night. 

420 self.mjd_log = [] 

421 self.hpid_log = [] 

422 

423 def add_observation(self, observation, indx=None): 

424 if observation['filter'][0] in self.filtername: 

425 if indx is None: 

426 indx = self.indx 

427 # Clear values if on a new night 

428 if self.night != observation['night']: 

429 self.feature *= 0. 

430 self.night = observation['night'] 

431 self.mjd_log = [] 

432 self.hpid_log = [] 

433 

434 # record the mjds and healpixels that were observed 

435 self.mjd_log.extend([np.max(observation['mjd'])]*np.size(indx)) 

436 self.hpid_log.extend(list(indx)) 

437 

438 # Look for the mjds that could possibly pair with observation 

439 tmin = observation['mjd'] - self.gap_max 

440 tmax = observation['mjd'] - self.gap_min 

441 mjd_log = np.array(self.mjd_log) 

442 left = np.searchsorted(mjd_log, tmin) 

443 right = np.searchsorted(mjd_log, tmax, side='right') 

444 # Now check if any of the healpixels taken in the time gap 

445 # match the healpixels of the observation. 

446 matches = np.in1d(indx, self.hpid_log[int(left):int(right)]) 

447 # XXX--should think if this is the correct (fastest) order to check things in. 

448 self.feature[indx[matches]] += 1 

449 

450 

451class Rotator_angle(BaseSurveyFeature): 

452 """ 

453 Track what rotation angles things are observed with. 

454 XXX-under construction 

455 """ 

456 def __init__(self, filtername='r', binsize=10., nside=None): 

457 """ 

458 

459 """ 

460 if nside is None: 

461 nside = utils.set_default_nside() 

462 

463 self.filtername = filtername 

464 # Actually keep a histogram at each healpixel 

465 self.feature = np.zeros((hp.nside2npix(nside), 360./binsize), dtype=float) 

466 self.bins = np.arange(0, 360+binsize, binsize) 

467 

468 def add_observation(self, observation, indx=None): 

469 if observation['filter'][0] == self.filtername: 

470 # I think this is how to broadcast things properly. 

471 self.feature[indx, :] += np.histogram(observation.rotSkyPos, bins=self.bins)[0]