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 lsst.sims.utils import (_hpid2RaDec, _raDec2Hpid, Site, calcLmstLast, 

3 m5_flat_sed, _approx_RaDec2AltAz, _angularSeparation) 

4import lsst.sims.skybrightness_pre as sb 

5import healpy as hp 

6from datetime import datetime 

7from lsst.sims.downtimeModel import ScheduledDowntimeData, UnscheduledDowntimeData 

8import lsst.sims.downtimeModel as downtimeModel 

9from lsst.sims.seeingModel import SeeingData, SeeingModel 

10from lsst.sims.cloudModel import CloudData 

11from lsst.sims.featureScheduler.features import Conditions 

12from lsst.sims.featureScheduler.utils import set_default_nside, approx_altaz2pa, create_season_offset 

13from lsst.ts.observatory.model import ObservatoryModel, Target 

14from astropy.coordinates import EarthLocation 

15from astropy.time import Time 

16from lsst.sims.almanac import Almanac 

17import warnings 

18import matplotlib.pylab as plt 

19from lsst.ts.observatory.model import ObservatoryState 

20from importlib import import_module 

21 

22__all__ = ['Model_observatory'] 

23 

24 

25class ExtendedObservatoryModel(ObservatoryModel): 

26 """Add some functionality to ObservatoryModel 

27 """ 

28 

29 def expose(self, target): 

30 # Break out the exposure command from observe method 

31 visit_time = sum(target.exp_times) + \ 

32 target.num_exp * self.params.shuttertime + \ 

33 max(target.num_exp - 1, 0) * self.params.readouttime 

34 self.update_state(self.current_state.time + visit_time) 

35 

36 def observe_times(self, target): 

37 """observe a target, and return the slewtime and visit time for the action 

38 Note, slew and expose will update the current_state 

39 """ 

40 t1 = self.current_state.time + 0 

41 # Note, this slew assumes there is a readout that needs to be done. 

42 self.slew(target) 

43 t2 = self.current_state.time + 0 

44 self.expose(target) 

45 t3 = self.current_state.time + 0 

46 if not self.current_state.tracking: 

47 ValueError('Telescope model stopped tracking, that seems bad.') 

48 slewtime = t2 - t1 

49 visitime = t3 - t2 

50 return slewtime, visitime 

51 

52 # Adding wrap_padding to make azimuth slews more intelligent 

53 def get_closest_angle_distance(self, target_rad, current_abs_rad, 

54 min_abs_rad=None, max_abs_rad=None, 

55 wrap_padding=0.873): 

56 """Calculate the closest angular distance including handling \ 

57 cable wrap if necessary. 

58 

59 Parameters 

60 ---------- 

61 target_rad : float 

62 The destination angle (radians). 

63 current_abs_rad : float 

64 The current angle (radians). 

65 min_abs_rad : float, optional 

66 The minimum constraint angle (radians). 

67 max_abs_rad : float, optional 

68 The maximum constraint angle (radians). 

69 wrap_padding : float (0.873) 

70 The amount of padding to use to make sure we don't track into limits (radians). 

71 

72 

73 Returns 

74 ------- 

75 tuple(float, float) 

76 (accumulated angle in radians, distance angle in radians) 

77 """ 

78 # if there are wrap limits, normalizes the target angle 

79 TWOPI = 2 * np.pi 

80 if min_abs_rad is not None: 

81 norm_target_rad = divmod(target_rad - min_abs_rad, TWOPI)[1] + min_abs_rad 

82 if max_abs_rad is not None: 

83 # if the target angle is unreachable 

84 # then sets an arbitrary value 

85 if norm_target_rad > max_abs_rad: 

86 norm_target_rad = max(min_abs_rad, norm_target_rad - np.pi) 

87 else: 

88 norm_target_rad = target_rad 

89 

90 # computes the distance clockwise 

91 distance_rad = divmod(norm_target_rad - current_abs_rad, TWOPI)[1] 

92 

93 # take the counter-clockwise distance if shorter 

94 if distance_rad > np.pi: 

95 distance_rad = distance_rad - TWOPI 

96 

97 # if there are wrap limits 

98 if (min_abs_rad is not None) and (max_abs_rad is not None): 

99 # compute accumulated angle 

100 accum_abs_rad = current_abs_rad + distance_rad 

101 

102 # if limits reached chose the other direction 

103 if accum_abs_rad > max_abs_rad - wrap_padding: 

104 distance_rad = distance_rad - TWOPI 

105 if accum_abs_rad < min_abs_rad + wrap_padding: 

106 distance_rad = distance_rad + TWOPI 

107 

108 # compute final accumulated angle 

109 final_abs_rad = current_abs_rad + distance_rad 

110 

111 return (final_abs_rad, distance_rad) 

112 

113 # Put in wrap padding kwarg so it's not used on camera rotation. 

114 def get_closest_state(self, targetposition, istracking=False): 

115 """Find the closest observatory state for the given target position. 

116 

117 Parameters 

118 ---------- 

119 targetposition : :class:`.ObservatoryPosition` 

120 A target position instance. 

121 istracking : bool, optional 

122 Flag for saying if the observatory is tracking. Default is False. 

123 

124 Returns 

125 ------- 

126 :class:`.ObservatoryState` 

127 The state that is closest to the current observatory state. 

128 

129 Binary schema 

130 ------------- 

131 The binary schema used to determine the state of a proposed target. A 

132 value of 1 indicates that is it failing. A value of 0 indicates that the 

133 state is passing. 

134 ___ ___ ___ ___ ___ ___ 

135 | | | | | | 

136 rot rot az az alt alt 

137 max min max min max min 

138 

139 For example, if a proposed target exceeds the rotators maximum value, 

140 and is below the minimum azimuth we would have a binary value of; 

141 

142 0 1 0 1 0 0 

143 

144 If the target passed, then no limitations would occur; 

145 

146 0 0 0 0 0 0 

147 """ 

148 TWOPI = 2 * np.pi 

149 

150 valid_state = True 

151 fail_record = self.current_state.fail_record 

152 self.current_state.fail_state = 0 

153 

154 if targetposition.alt_rad < self.params.telalt_minpos_rad: 

155 telalt_rad = self.params.telalt_minpos_rad 

156 domalt_rad = self.params.telalt_minpos_rad 

157 valid_state = False 

158 

159 if "telalt_minpos_rad" in fail_record: 

160 fail_record["telalt_minpos_rad"] += 1 

161 else: 

162 fail_record["telalt_minpos_rad"] = 1 

163 

164 self.current_state.fail_state = self.current_state.fail_state | \ 

165 self.current_state.fail_value_table["altEmin"] 

166 

167 elif targetposition.alt_rad > self.params.telalt_maxpos_rad: 

168 telalt_rad = self.params.telalt_maxpos_rad 

169 domalt_rad = self.params.telalt_maxpos_rad 

170 valid_state = False 

171 if "telalt_maxpos_rad" in fail_record: 

172 fail_record["telalt_maxpos_rad"] += 1 

173 else: 

174 fail_record["telalt_maxpos_rad"] = 1 

175 

176 self.current_state.fail_state = self.current_state.fail_state | \ 

177 self.current_state.fail_value_table["altEmax"] 

178 

179 else: 

180 telalt_rad = targetposition.alt_rad 

181 domalt_rad = targetposition.alt_rad 

182 

183 if istracking: 

184 (telaz_rad, delta) = self.get_closest_angle_distance(targetposition.az_rad, 

185 self.current_state.telaz_rad) 

186 if telaz_rad < self.params.telaz_minpos_rad: 

187 telaz_rad = self.params.telaz_minpos_rad 

188 valid_state = False 

189 if "telaz_minpos_rad" in fail_record: 

190 fail_record["telaz_minpos_rad"] += 1 

191 else: 

192 fail_record["telaz_minpos_rad"] = 1 

193 

194 self.current_state.fail_state = self.current_state.fail_state | \ 

195 self.current_state.fail_value_table["azEmin"] 

196 

197 elif telaz_rad > self.params.telaz_maxpos_rad: 

198 telaz_rad = self.params.telaz_maxpos_rad 

199 valid_state = False 

200 if "telaz_maxpos_rad" in fail_record: 

201 fail_record["telaz_maxpos_rad"] += 1 

202 else: 

203 fail_record["telaz_maxpos_rad"] = 1 

204 

205 self.current_state.fail_state = self.current_state.fail_state | \ 

206 self.current_state.fail_value_table["azEmax"] 

207 

208 else: 

209 (telaz_rad, delta) = self.get_closest_angle_distance(targetposition.az_rad, 

210 self.current_state.telaz_rad, 

211 self.params.telaz_minpos_rad, 

212 self.params.telaz_maxpos_rad) 

213 

214 (domaz_rad, delta) = self.get_closest_angle_distance(targetposition.az_rad, 

215 self.current_state.domaz_rad) 

216 

217 if istracking: 

218 (telrot_rad, delta) = self.get_closest_angle_distance(targetposition.rot_rad, 

219 self.current_state.telrot_rad, 

220 wrap_padding=0.) 

221 if telrot_rad < self.params.telrot_minpos_rad: 

222 telrot_rad = self.params.telrot_minpos_rad 

223 valid_state = False 

224 if "telrot_minpos_rad" in fail_record: 

225 fail_record["telrot_minpos_rad"] += 1 

226 else: 

227 fail_record["telrot_minpos_rad"] = 1 

228 

229 self.current_state.fail_state = self.current_state.fail_state | \ 

230 self.current_state.fail_value_table["rotEmin"] 

231 

232 elif telrot_rad > self.params.telrot_maxpos_rad: 

233 telrot_rad = self.params.telrot_maxpos_rad 

234 valid_state = False 

235 if "telrot_maxpos_rad" in fail_record: 

236 fail_record["telrot_maxpos_rad"] += 1 

237 else: 

238 fail_record["telrot_maxpos_rad"] = 1 

239 

240 self.current_state.fail_state = self.current_state.fail_state | \ 

241 self.current_state.fail_value_table["rotEmax"] 

242 else: 

243 # if the target rotator angle is unreachable 

244 # then sets an arbitrary value (opposite) 

245 norm_rot_rad = divmod(targetposition.rot_rad - self.params.telrot_minpos_rad, TWOPI)[1] \ 

246 + self.params.telrot_minpos_rad 

247 if norm_rot_rad > self.params.telrot_maxpos_rad: 

248 targetposition.rot_rad = norm_rot_rad - np.pi 

249 (telrot_rad, delta) = self.get_closest_angle_distance(targetposition.rot_rad, 

250 self.current_state.telrot_rad, 

251 self.params.telrot_minpos_rad, 

252 self.params.telrot_maxpos_rad, 

253 wrap_padding=0.) 

254 targetposition.ang_rad = divmod(targetposition.pa_rad - telrot_rad, TWOPI)[1] 

255 

256 targetstate = ObservatoryState() 

257 targetstate.set_position(targetposition) 

258 targetstate.telalt_rad = telalt_rad 

259 targetstate.telaz_rad = telaz_rad 

260 targetstate.telrot_rad = telrot_rad 

261 targetstate.domalt_rad = domalt_rad 

262 targetstate.domaz_rad = domaz_rad 

263 if istracking: 

264 targetstate.tracking = valid_state 

265 

266 self.current_state.fail_record = fail_record 

267 

268 return targetstate 

269 

270 

271class Model_observatory(object): 

272 """A class to generate a realistic telemetry stream for the scheduler 

273 """ 

274 

275 def __init__(self, nside=None, mjd_start=59853.5, seed=42, quickTest=True, 

276 alt_min=5., lax_dome=True, cloud_limit=0.3, sim_ToO=None, 

277 seeing_db=None): 

278 """ 

279 Parameters 

280 ---------- 

281 nside : int (None) 

282 The healpix nside resolution 

283 mjd_start : float (59853.5) 

284 The MJD to start the observatory up at 

285 alt_min : float (5.) 

286 The minimum altitude to compute models at (degrees). 

287 lax_dome : bool (True) 

288 Passed to observatory model. If true, allows dome creep. 

289 cloud_limit : float (0.3) 

290 The limit to stop taking observations if the cloud model returns something equal or higher 

291 sim_ToO : sim_targetoO object (None) 

292 If one would like to inject simulated ToOs into the telemetry stream. 

293 seeing_db : filename of the seeing data database (None) 

294 If one would like to use an alternate seeing database 

295 """ 

296 

297 if nside is None: 

298 nside = set_default_nside() 

299 self.nside = nside 

300 

301 self.cloud_limit = cloud_limit 

302 

303 self.alt_min = np.radians(alt_min) 

304 self.lax_dome = lax_dome 

305 

306 self.mjd_start = mjd_start 

307 

308 self.sim_ToO = sim_ToO 

309 

310 # Create an astropy location 

311 self.site = Site('LSST') 

312 self.location = EarthLocation(lat=self.site.latitude, lon=self.site.longitude, 

313 height=self.site.height) 

314 

315 # Load up all the models we need 

316 

317 mjd_start_time = Time(self.mjd_start, format='mjd') 

318 # Downtime 

319 self.down_nights = [] 

320 self.sched_downtime_data = ScheduledDowntimeData(mjd_start_time) 

321 self.unsched_downtime_data = UnscheduledDowntimeData(mjd_start_time) 

322 

323 sched_downtimes = self.sched_downtime_data() 

324 unsched_downtimes = self.unsched_downtime_data() 

325 

326 down_starts = [] 

327 down_ends = [] 

328 for dt in sched_downtimes: 

329 down_starts.append(dt['start'].mjd) 

330 down_ends.append(dt['end'].mjd) 

331 for dt in unsched_downtimes: 

332 down_starts.append(dt['start'].mjd) 

333 down_ends.append(dt['end'].mjd) 

334 

335 self.downtimes = np.array(list(zip(down_starts, down_ends)), dtype=list(zip(['start', 'end'], [float, float]))) 

336 self.downtimes.sort(order='start') 

337 

338 # Make sure there aren't any overlapping downtimes 

339 diff = self.downtimes['start'][1:] - self.downtimes['end'][0:-1] 

340 while np.min(diff) < 0: 

341 # Should be able to do this wihtout a loop, but this works 

342 for i, dt in enumerate(self.downtimes[0:-1]): 

343 if self.downtimes['start'][i+1] < dt['end']: 

344 new_end = np.max([dt['end'], self.downtimes['end'][i+1]]) 

345 self.downtimes[i]['end'] = new_end 

346 self.downtimes[i+1]['end'] = new_end 

347 

348 good = np.where(self.downtimes['end'] - np.roll(self.downtimes['end'], 1) != 0) 

349 self.downtimes = self.downtimes[good] 

350 diff = self.downtimes['start'][1:] - self.downtimes['end'][0:-1] 

351 

352 self.seeing_data = SeeingData(mjd_start_time, seeing_db=seeing_db) 

353 self.seeing_model = SeeingModel() 

354 self.seeing_indx_dict = {} 

355 for i, filtername in enumerate(self.seeing_model.filter_list): 

356 self.seeing_indx_dict[filtername] = i 

357 

358 self.cloud_data = CloudData(mjd_start_time, offset_year=0) 

359 

360 self.sky_model = sb.SkyModelPre(speedLoad=quickTest) 

361 

362 self.observatory = ExtendedObservatoryModel() 

363 self.observatory.configure_from_module() 

364 # Make it so it respects my requested rotator angles 

365 self.observatory.params.rotator_followsky = True 

366 

367 self.filterlist = ['u', 'g', 'r', 'i', 'z', 'y'] 

368 self.seeing_FWHMeff = {} 

369 for key in self.filterlist: 

370 self.seeing_FWHMeff[key] = np.zeros(hp.nside2npix(self.nside), dtype=float) 

371 

372 self.almanac = Almanac(mjd_start=mjd_start) 

373 

374 # Let's make sure we're at an openable MJD 

375 good_mjd = False 

376 to_set_mjd = mjd_start 

377 while not good_mjd: 

378 good_mjd, to_set_mjd = self.check_mjd(to_set_mjd) 

379 self.mjd = to_set_mjd 

380 

381 sun_moon_info = self.almanac.get_sun_moon_positions(self.mjd) 

382 season_offset = create_season_offset(self.nside, sun_moon_info['sun_RA']) 

383 self.sun_RA_start = sun_moon_info['sun_RA'] + 0 

384 # Conditions object to update and return on request 

385 self.conditions = Conditions(nside=self.nside, mjd_start=mjd_start, 

386 season_offset=season_offset, sun_RA_start=self.sun_RA_start) 

387 

388 

389 self.obsID_counter = 0 

390 

391 def get_info(self): 

392 """ 

393 Returns 

394 ------- 

395 Array with model versions that were instantiated 

396 """ 

397 

398 # The things we want to get info on 

399 models = {'cloud data': self.cloud_data, 'sky model': self.sky_model, 

400 'seeing data': self.seeing_data, 'seeing model': self.seeing_model, 

401 'observatory model': self.observatory, 

402 'sched downtime data': self.sched_downtime_data, 

403 'unched downtime data': self.unsched_downtime_data} 

404 

405 result = [] 

406 for model_name in models: 

407 try: 

408 module_name = models[model_name].__module__ 

409 module = import_module(module_name) 

410 ver = import_module(module.__package__+'.version') 

411 version = ver.__version__ 

412 fingerprint = ver.__fingerprint__ 

413 except: 

414 version = 'NA' 

415 fingerprint = 'NA' 

416 result.append([model_name+' version', version]) 

417 result.append([model_name+' fingerprint', fingerprint]) 

418 result.append([model_name+' module', models[model_name].__module__]) 

419 try: 

420 info = models[model_name].config_info() 

421 for key in info: 

422 result.append([key, str(info[key])]) 

423 except: 

424 result.append([model_name, 'no config_info']) 

425 

426 return result 

427 

428 def return_conditions(self): 

429 """ 

430 

431 Returns 

432 ------- 

433 lsst.sims.featureScheduler.features.conditions object 

434 """ 

435 

436 self.conditions.mjd = self.mjd 

437 

438 self.conditions.night = self.night 

439 # Current time as astropy time 

440 current_time = Time(self.mjd, format='mjd') 

441 

442 # Clouds. XXX--just the raw value 

443 self.conditions.bulk_cloud = self.cloud_data(current_time) 

444 

445 # use conditions object itself to get aprox altitude of each healpx 

446 alts = self.conditions.alt 

447 azs = self.conditions.az 

448 

449 good = np.where(alts > self.alt_min) 

450 

451 # Compute the airmass at each heapix 

452 airmass = np.zeros(alts.size, dtype=float) 

453 airmass.fill(np.nan) 

454 airmass[good] = 1./np.cos(np.pi/2. - alts[good]) 

455 self.conditions.airmass = airmass 

456 

457 # reset the seeing 

458 for key in self.seeing_FWHMeff: 

459 self.seeing_FWHMeff[key].fill(np.nan) 

460 # Use the model to get the seeing at this time and airmasses. 

461 FWHM_500 = self.seeing_data(current_time) 

462 seeing_dict = self.seeing_model(FWHM_500, airmass[good]) 

463 fwhm_eff = seeing_dict['fwhmEff'] 

464 for i, key in enumerate(self.seeing_model.filter_list): 

465 self.seeing_FWHMeff[key][good] = fwhm_eff[i, :] 

466 self.conditions.FWHMeff = self.seeing_FWHMeff 

467 

468 # sky brightness 

469 self.conditions.skybrightness = self.sky_model.returnMags(self.mjd, airmass_mask=False, 

470 planet_mask=False, 

471 moon_mask=False, zenith_mask=False) 

472 

473 self.conditions.mounted_filters = self.observatory.current_state.mountedfilters 

474 self.conditions.current_filter = self.observatory.current_state.filter[0] 

475 

476 # Compute the slewtimes 

477 slewtimes = np.empty(alts.size, dtype=float) 

478 slewtimes.fill(np.nan) 

479 slewtimes[good] = self.observatory.get_approximate_slew_delay(alts[good], azs[good], 

480 self.observatory.current_state.filter, 

481 lax_dome=self.lax_dome) 

482 # Mask out anything the slewtime says is out of bounds 

483 slewtimes[np.where(slewtimes < 0)] = np.nan 

484 self.conditions.slewtime = slewtimes 

485 

486 # Let's get the sun and moon 

487 sun_moon_info = self.almanac.get_sun_moon_positions(self.mjd) 

488 # convert these to scalars 

489 for key in sun_moon_info: 

490 sun_moon_info[key] = sun_moon_info[key].max() 

491 self.conditions.moonPhase = sun_moon_info['moon_phase'] 

492 

493 self.conditions.moonAlt = sun_moon_info['moon_alt'] 

494 self.conditions.moonAz = sun_moon_info['moon_az'] 

495 self.conditions.moonRA = sun_moon_info['moon_RA'] 

496 self.conditions.moonDec = sun_moon_info['moon_dec'] 

497 self.conditions.sunAlt = sun_moon_info['sun_alt'] 

498 self.conditions.sunRA = sun_moon_info['sun_RA'] 

499 self.conditions.sunDec = sun_moon_info['sun_dec'] 

500 

501 self.conditions.lmst, last = calcLmstLast(self.mjd, self.site.longitude_rad) 

502 

503 self.conditions.telRA = self.observatory.current_state.ra_rad 

504 self.conditions.telDec = self.observatory.current_state.dec_rad 

505 self.conditions.telAlt = self.observatory.current_state.alt_rad 

506 self.conditions.telAz = self.observatory.current_state.az_rad 

507 

508 self.conditions.rotTelPos = self.observatory.current_state.rot_rad 

509 

510 # Add in the almanac information 

511 self.conditions.night = self.night 

512 self.conditions.sunset = self.almanac.sunsets['sunset'][self.almanac_indx] 

513 self.conditions.sun_n12_setting = self.almanac.sunsets['sun_n12_setting'][self.almanac_indx] 

514 self.conditions.sun_n18_setting = self.almanac.sunsets['sun_n18_setting'][self.almanac_indx] 

515 self.conditions.sun_n18_rising = self.almanac.sunsets['sun_n18_rising'][self.almanac_indx] 

516 self.conditions.sun_n12_rising = self.almanac.sunsets['sun_n12_rising'][self.almanac_indx] 

517 self.conditions.sunrise = self.almanac.sunsets['sunrise'][self.almanac_indx] 

518 self.conditions.moonrise = self.almanac.sunsets['moonrise'][self.almanac_indx] 

519 self.conditions.moonset = self.almanac.sunsets['moonset'][self.almanac_indx] 

520 

521 # Planet positions from almanac 

522 self.conditions.planet_positions = self.almanac.get_planet_positions(self.mjd) 

523 

524 # See if there are any ToOs to include 

525 if self.sim_ToO is not None: 

526 toos = self.sim_ToO(self.mjd) 

527 if toos is not None: 

528 self.conditions.targets_of_opportunity = toos 

529 

530 return self.conditions 

531 

532 @property 

533 def mjd(self): 

534 return self._mjd 

535 

536 @mjd.setter 

537 def mjd(self, value): 

538 self._mjd = value 

539 self.almanac_indx = self.almanac.mjd_indx(value) 

540 self.night = self.almanac.sunsets['night'][self.almanac_indx] 

541 

542 def observation_add_data(self, observation): 

543 """ 

544 Fill in the metadata for a completed observation 

545 """ 

546 current_time = Time(self.mjd, format='mjd') 

547 

548 observation['clouds'] = self.cloud_data(current_time) 

549 observation['airmass'] = 1./np.cos(np.pi/2. - observation['alt']) 

550 # Seeing 

551 fwhm_500 = self.seeing_data(current_time) 

552 seeing_dict = self.seeing_model(fwhm_500, observation['airmass']) 

553 observation['FWHMeff'] = seeing_dict['fwhmEff'][self.seeing_indx_dict[observation['filter'][0]]] 

554 observation['FWHM_geometric'] = seeing_dict['fwhmGeom'][self.seeing_indx_dict[observation['filter'][0]]] 

555 observation['FWHM_500'] = fwhm_500 

556 

557 observation['night'] = self.night 

558 observation['mjd'] = self.mjd 

559 

560 hpid = _raDec2Hpid(self.sky_model.nside, observation['RA'], observation['dec']) 

561 observation['skybrightness'] = self.sky_model.returnMags(self.mjd, 

562 indx=[hpid], 

563 extrapolate=True)[observation['filter'][0]] 

564 

565 observation['fivesigmadepth'] = m5_flat_sed(observation['filter'][0], observation['skybrightness'], 

566 observation['FWHMeff'], 

567 observation['exptime']/observation['nexp'], 

568 observation['airmass'], nexp=observation['nexp']) 

569 

570 lmst, last = calcLmstLast(self.mjd, self.site.longitude_rad) 

571 observation['lmst'] = lmst 

572 

573 sun_moon_info = self.almanac.get_sun_moon_positions(self.mjd) 

574 observation['sunAlt'] = sun_moon_info['sun_alt'] 

575 observation['sunAz'] = sun_moon_info['sun_az'] 

576 observation['sunRA'] = sun_moon_info['sun_RA'] 

577 observation['sunDec'] = sun_moon_info['sun_dec'] 

578 observation['moonAlt'] = sun_moon_info['moon_alt'] 

579 observation['moonAz'] = sun_moon_info['moon_az'] 

580 observation['moonRA'] = sun_moon_info['moon_RA'] 

581 observation['moonDec'] = sun_moon_info['moon_dec'] 

582 observation['moonDist'] = _angularSeparation(observation['RA'], observation['dec'], 

583 observation['moonRA'], observation['moonDec']) 

584 observation['solarElong'] = _angularSeparation(observation['RA'], observation['dec'], 

585 observation['sunRA'], observation['sunDec']) 

586 observation['moonPhase'] = sun_moon_info['moon_phase'] 

587 

588 observation['ID'] = self.obsID_counter 

589 self.obsID_counter += 1 

590 

591 return observation 

592 

593 def check_up(self, mjd): 

594 """See if we are in downtime 

595 

596 True if telescope is up 

597 False if in downtime 

598 """ 

599 

600 result = True 

601 indx = np.searchsorted(self.downtimes['start'], mjd, side='right')-1 

602 if mjd < self.downtimes['end'][indx]: 

603 result = False 

604 return result 

605 

606 def check_mjd(self, mjd, cloud_skip=20.): 

607 """See if an mjd is ok to observe 

608 Parameters 

609 ---------- 

610 cloud_skip : float (20) 

611 How much time to skip ahead if it's cloudy (minutes) 

612 

613 

614 Returns 

615 ------- 

616 bool 

617 

618 mdj : float 

619 If True, the input mjd. If false, a good mjd to skip forward to. 

620 """ 

621 passed = True 

622 new_mjd = mjd + 0 

623 

624 # Maybe set this to a while loop to make sure we don't land on another cloudy time? 

625 # or just make this an entire recursive call? 

626 clouds = self.cloud_data(Time(mjd, format='mjd')) 

627 if clouds > self.cloud_limit: 

628 passed = False 

629 while clouds > self.cloud_limit: 

630 new_mjd = new_mjd + cloud_skip/60./24. 

631 clouds = self.cloud_data(Time(new_mjd, format='mjd')) 

632 alm_indx = np.searchsorted(self.almanac.sunsets['sunset'], mjd) - 1 

633 # at the end of the night, advance to the next setting twilight 

634 if mjd > self.almanac.sunsets['sun_n12_rising'][alm_indx]: 

635 passed = False 

636 new_mjd = self.almanac.sunsets['sun_n12_setting'][alm_indx+1] 

637 if mjd < self.almanac.sunsets['sun_n12_setting'][alm_indx]: 

638 passed = False 

639 new_mjd = self.almanac.sunsets['sun_n12_setting'][alm_indx+1] 

640 # We're in a down night, advance to next night 

641 if not self.check_up(mjd): 

642 passed = False 

643 new_mjd = self.almanac.sunsets['sun_n12_setting'][alm_indx+1] 

644 # recursive call to make sure we skip far enough ahead 

645 if not passed: 

646 while not passed: 

647 passed, new_mjd = self.check_mjd(new_mjd) 

648 return False, new_mjd 

649 else: 

650 return True, mjd 

651 

652 def _update_rotSkyPos(self, observation): 

653 """If we have an undefined rotSkyPos, try to fill it out. 

654 """ 

655 alt, az = _approx_RaDec2AltAz(observation['RA'], observation['dec'], self.site.latitude_rad, 

656 self.site.longitude_rad, self.mjd) 

657 obs_pa = approx_altaz2pa(alt, az, self.site.latitude_rad) 

658 observation['rotSkyPos'] = (obs_pa + observation['rotTelPos']) % (2*np.pi) 

659 observation['rotTelPos'] = 0. 

660 

661 return observation 

662 

663 def observe(self, observation): 

664 """Try to make an observation 

665 

666 Returns 

667 ------- 

668 observation : observation object 

669 None if there was no observation taken. Completed observation with meta data filled in. 

670 new_night : bool 

671 Have we started a new night. 

672 """ 

673 

674 start_night = self.night.copy() 

675 

676 # Make sure the kinematic model is set to the correct mjd 

677 t = Time(self.mjd, format='mjd') 

678 self.observatory.update_state(t.unix) 

679 

680 if np.isnan(observation['rotSkyPos']): 

681 observation = self._update_rotSkyPos(observation) 

682 

683 target = Target(band_filter=observation['filter'], ra_rad=observation['RA'], 

684 dec_rad=observation['dec'], ang_rad=observation['rotSkyPos'], 

685 num_exp=observation['nexp'], exp_times=[observation['exptime']]) 

686 start_ra = self.observatory.current_state.ra_rad 

687 start_dec = self.observatory.current_state.dec_rad 

688 slewtime, visittime = self.observatory.observe_times(target) 

689 

690 # Check if the mjd after slewtime and visitime is fine: 

691 observation_worked, new_mjd = self.check_mjd(self.mjd + (slewtime + visittime)/24./3600.) 

692 if observation_worked: 

693 observation['visittime'] = visittime 

694 observation['slewtime'] = slewtime 

695 observation['slewdist'] = _angularSeparation(start_ra, start_dec, 

696 self.observatory.current_state.ra_rad, 

697 self.observatory.current_state.dec_rad) 

698 self.mjd = self.mjd + slewtime/24./3600. 

699 # Reach into the observatory model to pull out the relevant data it has calculated 

700 # Note, this might be after the observation has been completed. 

701 observation['alt'] = self.observatory.current_state.alt_rad 

702 observation['az'] = self.observatory.current_state.az_rad 

703 observation['pa'] = self.observatory.current_state.pa_rad 

704 observation['rotTelPos'] = self.observatory.current_state.rot_rad 

705 observation['rotSkyPos'] = self.observatory.current_state.ang_rad 

706 

707 # Metadata on observation is after slew and settle, so at start of exposure. 

708 result = self.observation_add_data(observation) 

709 self.mjd = self.mjd + visittime/24./3600. 

710 new_night = False 

711 else: 

712 result = None 

713 self.observatory.park() 

714 # Skip to next legitimate mjd 

715 self.mjd = new_mjd 

716 now_night = self.night 

717 if now_night == start_night: 

718 new_night = False 

719 else: 

720 new_night = True 

721 

722 return result, new_night