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 builtins import object 

2import numpy as np 

3import glob 

4import os 

5import healpy as hp 

6from lsst.utils import getPackageDir 

7import warnings 

8from lsst.sims.utils import _angularSeparation 

9 

10__all__ = ['SkyModelPre', 'interp_angle'] 

11 

12 

13def shortAngleDist(a0, a1): 

14 """ 

15 from https://gist.github.com/shaunlebron/8832585 

16 """ 

17 max_angle = 2.*np.pi 

18 da = (a1 - a0) % max_angle 

19 return 2.*da % max_angle - da 

20 

21 

22def interp_angle(x_out, xp, anglep, degrees=False): 

23 """ 

24 Interpolate angle values (handle wrap around properly). Does nearest neighbor 

25 interpolation if values out of range. 

26 

27 Parameters 

28 ---------- 

29 x_out : float (or array) 

30 The points to interpolate to. 

31 xp : array 

32 Points to interpolate between (must be sorted) 

33 anglep : array 

34 The angles ascociated with xp 

35 degrees : bool (False) 

36 Set if anglep is degrees (True) or radidian (False) 

37 """ 

38 

39 # Where are the interpolation points 

40 x = np.atleast_1d(x_out) 

41 left = np.searchsorted(xp, x)-1 

42 right = left+1 

43 

44 # If we are out of bounds, just use the edges 

45 right[np.where(right >= xp.size)] -= 1 

46 left[np.where(left < 0)] += 1 

47 baseline = xp[right] - xp[left] 

48 

49 wterm = (x - xp[left])/baseline 

50 wterm[np.where(baseline == 0)] = 0 

51 if degrees: 

52 result = np.radians(anglep[left]) + shortAngleDist(np.radians(anglep[left]), np.radians(anglep[right]))*wterm 

53 result = result % (2.*np.pi) 

54 result = np.degrees(result) 

55 else: 

56 result = anglep[left] + shortAngleDist(anglep[left], anglep[right])*wterm 

57 result = result % (2.*np.pi) 

58 return result 

59 

60 

61class SkyModelPre(object): 

62 """ 

63 Load pre-computed sky brighntess maps for the LSST site and use them to interpolate to 

64 arbitrary dates. 

65 """ 

66 

67 def __init__(self, data_path=None, opsimFields=False, 

68 speedLoad=True, verbose=False): 

69 """ 

70 Parameters 

71 ---------- 

72 data_path : str (None) 

73 path to the numpy save files. Looks in standard plances if set to None. 

74 opsimFields : bool (False) 

75 Mostly depreciated, if True, loads sky brightnesses computed at field centers. 

76 Otherwise uses healpixels. 

77 speedLoad : bool (True) 

78 If True, use the small 3-day file to load found in the usual spot. 

79 """ 

80 

81 self.info = None 

82 self.sb = None 

83 self.header = None 

84 self.filter_names = None 

85 

86 self.opsimFields = opsimFields 

87 self.verbose = verbose 

88 

89 # Look in default location for .npz files to load 

90 if 'SIMS_SKYBRIGHTNESS_DATA' in os.environ: 

91 data_dir = os.environ['SIMS_SKYBRIGHTNESS_DATA'] 

92 else: 

93 data_dir = os.path.join(getPackageDir('sims_skybrightness_pre'), 'data') 

94 

95 if data_path is None: 

96 if opsimFields: 

97 data_path = os.path.join(data_dir, 'opsimFields') 

98 else: 

99 data_path = os.path.join(data_dir, 'healpix') 

100 # Expect filenames of the form mjd1_mjd2.npz, e.g., 59632.155_59633.2.npz 

101 self.files = glob.glob(os.path.join(data_path, '*.npz')) 

102 if len(self.files) == 0: 

103 errmssg = 'Failed to find pre-computed .npz files. ' 

104 errmssg += 'Copy data from NCSA with sims_skybrightness_pre/data/data_down.sh \n' 

105 errmssg += 'or build by running sims_skybrightness_pre/data/generate_sky.py' 

106 warnings.warn(errmssg) 

107 self.filesizes = np.array([os.path.getsize(filename) for filename in self.files]) 

108 mjd_left = [] 

109 mjd_right = [] 

110 # glob does not always order things I guess? 

111 self.files.sort() 

112 for filename in self.files: 

113 temp = os.path.split(filename)[-1].replace('.npz', '').split('_') 

114 mjd_left.append(float(temp[0])) 

115 mjd_right.append(float(temp[1])) 

116 

117 self.mjd_left = np.array(mjd_left) 

118 self.mjd_right = np.array(mjd_right) 

119 

120 # Set that nothing is loaded at this point 

121 self.loaded_range = np.array([-1]) 

122 

123 # Go ahead and load the small one in the repo by default 

124 if speedLoad: 

125 self._load_data(59853., 

126 filename=os.path.join(data_dir, 'healpix/59853_59856.npz'), 

127 npyfile=os.path.join(data_dir, 'healpix/59853_59856.npy')) 

128 

129 def _load_data(self, mjd, filename=None, npyfile=None): 

130 """ 

131 Load up the .npz file to interpolate things. After python 3 upgrade, numpy.savez refused 

132 to write large .npz files, so data is split between .npz and .npy files. 

133 

134 Parameters 

135 ---------- 

136 mjd : float 

137 The Modified Julian Date that we want to load 

138 filename : str (None) 

139 The filename to restore. If None, it checks the filenames on disk to find one that 

140 should have the requested MJD 

141 npyfile : str (None) 

142 If sky brightness data not in npz file, checks the .npy file with same root name. 

143 """ 

144 del self.info 

145 del self.sb 

146 del self.header 

147 del self.filter_names 

148 

149 if filename is None: 

150 # Figure out which file to load. 

151 file_indx = np.where((mjd >= self.mjd_left) & (mjd <= self.mjd_right))[0] 

152 if np.size(file_indx) == 0: 

153 raise ValueError('MJD = %f is out of range for the files found (%f-%f)' % (mjd, 

154 self.mjd_left.min(), 

155 self.mjd_right.max())) 

156 if np.size(file_indx) > 1: 

157 # Select the smaller file 

158 sizes = self.filesizes[file_indx] 

159 file_indx = np.min(file_indx[np.where(sizes == np.min(sizes))]) 

160 else: 

161 file_indx = np.min(file_indx) 

162 

163 filename = self.files[file_indx] 

164 

165 self.loaded_range = np.array([self.mjd_left[file_indx], self.mjd_right[file_indx]]) 

166 else: 

167 self.loaded_range = None 

168 

169 if self.verbose: 

170 print('Loading file %s' % filename) 

171 # Add encoding kwarg to restore Python 2.7 generated files 

172 data = np.load(filename, encoding='bytes', allow_pickle=True) 

173 self.info = data['dict_of_lists'][()] 

174 self.header = data['header'][()] 

175 if 'sky_brightness' in data.keys(): 

176 self.sb = data['sky_brightness'][()] 

177 data.close() 

178 else: 

179 # the sky brightness had to go in it's own npy file 

180 data.close() 

181 if npyfile is None: 

182 npyfile = filename[:-3]+'npy' 

183 self.sb = np.load(npyfile) 

184 if self.verbose: 

185 print('also loading %s' % npyfile) 

186 

187 # Step to make sure keys are strings not bytes 

188 all_dicts = [self.info, self.sb, self.header] 

189 all_dicts = [single_dict for single_dict in all_dicts if hasattr(single_dict, 'keys')] 

190 for selfDict in all_dicts: 

191 for key in list(selfDict.keys()): 

192 if type(key) != str: 

193 selfDict[key.decode("utf-8")] = selfDict.pop(key) 

194 

195 # Ugh, different versions of the save files could have dicts or np.array. 

196 # Let's hope someone fits some Fourier components to the sky brightnesses and gets rid 

197 # of the giant save files for good. 

198 if hasattr(self.sb, 'keys'): 

199 self.filter_names = list(self.sb.keys()) 

200 else: 

201 self.filter_names = self.sb.dtype.names 

202 

203 if self.verbose: 

204 print('%s loaded' % os.path.split(filename)[1]) 

205 

206 if not self.opsimFields: 

207 self.nside = hp.npix2nside(self.sb[self.filter_names[0]][0, :].size) 

208 

209 if self.loaded_range is None: 

210 self.loaded_range = np.array([self.info['mjds'].min(), self.info['mjds'].max()]) 

211 

212 def returnSunMoon(self, mjd): 

213 """ 

214 Parameters 

215 ---------- 

216 mjd : float 

217 Modified Julian Date(s) to interpolate to 

218 

219 Returns 

220 ------- 

221 sunMoon : dict 

222 Dict with keys for the sun and moon RA and Dec and the 

223 mooon-sun separation. All values in radians, except for moonSunSep 

224 that is in degrees for some reason (that reason is probably because I'm sloppy). 

225 """ 

226 

227 #warnings.warn('Method returnSunMoon to be depreciated. Interpolating angles is bad!') 

228 

229 keys = ['sunAlts', 'moonAlts', 'moonRAs', 'moonDecs', 'sunRAs', 

230 'sunDecs', 'moonSunSep'] 

231 

232 degrees = [False, False, False, False, False, False, True] 

233 

234 if (mjd < self.loaded_range.min() or (mjd > self.loaded_range.max())): 

235 self._load_data(mjd) 

236 

237 result = {} 

238 for key, degree in zip(keys, degrees): 

239 if key[-1] == 's': 

240 newkey = key[:-1] 

241 else: 

242 newkey = key 

243 if 'RA' in key: 

244 result[newkey] = interp_angle(mjd, self.info['mjds'], self.info[key], degrees=degree) 

245 # Return a scalar if only doing 1 date. 

246 if np.size(result[newkey]) == 1: 

247 result[newkey] = np.max(result[newkey]) 

248 else: 

249 result[newkey] = np.interp(mjd, self.info['mjds'], self.info[key]) 

250 return result 

251 

252 def returnAirmass(self, mjd, maxAM=10., indx=None, badval=hp.UNSEEN): 

253 """ 

254 

255 Parameters 

256 ---------- 

257 mjd : float 

258 Modified Julian Date to interpolate to 

259 indx : List of int(s) (None) 

260 indices to interpolate the sky values at. Returns full sky if None. If the class was 

261 instatiated with opsimFields, indx is the field ID, otherwise it is the healpix ID. 

262 maxAM : float (10) 

263 The maximum airmass to return, everything above this airmass will be set to badval 

264 

265 Returns 

266 ------- 

267 airmass : np.array 

268 Array of airmass values. If the MJD is between sunrise and sunset, all values are masked. 

269 """ 

270 if (mjd < self.loaded_range.min() or (mjd > self.loaded_range.max())): 

271 self._load_data(mjd) 

272 

273 left = np.searchsorted(self.info['mjds'], mjd)-1 

274 right = left+1 

275 

276 # If we are out of bounds 

277 if right >= self.info['mjds'].size: 

278 right -= 1 

279 baseline = 1. 

280 elif left < 0: 

281 left += 1 

282 baseline = 1. 

283 else: 

284 baseline = self.info['mjds'][right] - self.info['mjds'][left] 

285 

286 if indx is None: 

287 result_size = self.sb[self.filter_names[0]][left, :].size 

288 indx = np.arange(result_size) 

289 else: 

290 result_size = len(indx) 

291 # Check if we are between sunrise/set 

292 if baseline > self.header['timestep_max']: 

293 warnings.warn('Requested MJD between sunrise and sunset, returning closest maps') 

294 diff = np.abs(self.info['mjds'][left.max():right.max()+1]-mjd) 

295 closest_indx = np.array([left, right])[np.where(diff == np.min(diff))] 

296 airmass = self.info['airmass'][closest_indx, indx] 

297 mask = np.where((self.info['airmass'][closest_indx, indx].ravel() < 1.) | 

298 (self.info['airmass'][closest_indx, indx].ravel() > maxAM)) 

299 airmass = airmass.ravel() 

300 

301 else: 

302 wterm = (mjd - self.info['mjds'][left])/baseline 

303 w1 = (1. - wterm) 

304 w2 = wterm 

305 airmass = self.info['airmass'][left, indx] * w1 + self.info['airmass'][right, indx] * w2 

306 mask = np.where((self.info['airmass'][left, indx] < 1.) | 

307 (self.info['airmass'][left, indx] > maxAM) | 

308 (self.info['airmass'][right, indx] < 1.) | 

309 (self.info['airmass'][right, indx] > maxAM)) 

310 airmass[mask] = badval 

311 

312 return airmass 

313 

314 def returnMags(self, mjd, indx=None, airmass_mask=True, planet_mask=True, 

315 moon_mask=True, zenith_mask=True, badval=hp.UNSEEN, 

316 filters=['u', 'g', 'r', 'i', 'z', 'y'], extrapolate=False): 

317 """ 

318 Return a full sky map or individual pixels for the input mjd 

319 

320 Parameters 

321 ---------- 

322 mjd : float 

323 Modified Julian Date to interpolate to 

324 indx : List of int(s) (None) 

325 indices to interpolate the sky values at. Returns full sky if None. If the class was 

326 instatiated with opsimFields, indx is the field ID, otherwise it is the healpix ID. 

327 airmass_mask : bool (True) 

328 Set high (>2.5) airmass pixels to badval. 

329 planet_mask : bool (True) 

330 Set sky maps to badval near (2 degrees) bright planets. 

331 moon_mask : bool (True) 

332 Set sky maps near (10 degrees) the moon to badval. 

333 zenith_mask : bool (True) 

334 Set sky maps at high altitude (>86.5) to badval. 

335 badval : float (-1.6375e30) 

336 Mask value. Defaults to the healpy mask value. 

337 filters : list 

338 List of strings for the filters that should be returned. 

339 extrapolate : bool (False) 

340 In indx is set, extrapolate any masked pixels to be the same as the nearest non-masked 

341 value from the full sky map. 

342 

343 Returns 

344 ------- 

345 sbs : dict 

346 A dictionary with filter names as keys and np.arrays as values which 

347 hold the sky brightness maps in mag/sq arcsec. 

348 """ 

349 if (mjd < self.loaded_range.min() or (mjd > self.loaded_range.max())): 

350 self._load_data(mjd) 

351 

352 mask_rules = {'airmass': airmass_mask, 'planet': planet_mask, 

353 'moon': moon_mask, 'zenith': zenith_mask} 

354 

355 left = np.searchsorted(self.info['mjds'], mjd)-1 

356 right = left+1 

357 

358 # Do full sky by default 

359 if indx is None: 

360 indx = np.arange(self.sb['r'].shape[1]) 

361 full_sky = True 

362 else: 

363 full_sky = False 

364 

365 # If we are out of bounds 

366 if right >= self.info['mjds'].size: 

367 right -= 1 

368 baseline = 1. 

369 elif left < 0: 

370 left += 1 

371 baseline = 1. 

372 else: 

373 baseline = self.info['mjds'][right] - self.info['mjds'][left] 

374 

375 # Check if we are between sunrise/set 

376 if baseline > self.header['timestep_max']: 

377 warnings.warn('Requested MJD between sunrise and sunset, returning closest maps') 

378 diff = np.abs(self.info['mjds'][left.max():right.max()+1]-mjd) 

379 closest_indx = np.array([left, right])[np.where(diff == np.min(diff))].min() 

380 sbs = {} 

381 for filter_name in filters: 

382 sbs[filter_name] = self.sb[filter_name][closest_indx, indx] 

383 for mask_name in mask_rules: 

384 if mask_rules[mask_name]: 

385 toMask = np.where(self.info[mask_name+'_masks'][closest_indx, indx]) 

386 sbs[filter_name][toMask] = badval 

387 sbs[filter_name][np.isinf(sbs[filter_name])] = badval 

388 sbs[filter_name][np.where(sbs[filter_name] == hp.UNSEEN)] = badval 

389 else: 

390 wterm = (mjd - self.info['mjds'][left])/baseline 

391 w1 = (1. - wterm) 

392 w2 = wterm 

393 sbs = {} 

394 for filter_name in filters: 

395 sbs[filter_name] = self.sb[filter_name][left, indx] * w1 + \ 

396 self.sb[filter_name][right, indx] * w2 

397 for mask_name in mask_rules: 

398 if mask_rules[mask_name]: 

399 toMask = np.where(self.info[mask_name+'_masks'][left, indx] | 

400 self.info[mask_name+'_masks'][right, indx] | 

401 np.isinf(sbs[filter_name])) 

402 sbs[filter_name][toMask] = badval 

403 sbs[filter_name][np.where(sbs[filter_name] == hp.UNSEEN)] = badval 

404 sbs[filter_name][np.where(sbs[filter_name] == hp.UNSEEN)] = badval 

405 

406 # If requested a certain pixel(s), and want to extrapolate. 

407 if (not full_sky) & extrapolate: 

408 masked_pix = False 

409 for filter_name in filters: 

410 if (badval in sbs[filter_name]) | (True in np.isnan(sbs[filter_name])): 

411 masked_pix = True 

412 if masked_pix: 

413 # We have pixels that are masked that we want reasonable values for 

414 full_sky_sb = self.returnMags(mjd, airmass_mask=False, planet_mask=False, moon_mask=False, 

415 zenith_mask=False, filters=filters) 

416 good = np.where((full_sky_sb[filters[0]] != badval) & ~np.isnan(full_sky_sb[filters[0]]))[0] 

417 ra_full = np.radians(self.header['ra'][good]) 

418 dec_full = np.radians(self.header['dec'][good]) 

419 for filtername in filters: 

420 full_sky_sb[filtername] = full_sky_sb[filtername][good] 

421 # Going to assume the masked pixels are the same in all filters 

422 masked_indx = np.where((sbs[filters[0]].ravel() == badval) | 

423 np.isnan(sbs[filters[0]].ravel()))[0] 

424 for i, mi in enumerate(masked_indx): 

425 # Note, this is going to be really slow for many pixels, should use a kdtree 

426 dist = _angularSeparation(np.radians(self.header['ra'][indx[i]]), 

427 np.radians(self.header['dec'][indx[i]]), 

428 ra_full, dec_full) 

429 closest = np.where(dist == dist.min())[0] 

430 for filtername in filters: 

431 sbs[filtername].ravel()[mi] = np.min(full_sky_sb[filtername][closest]) 

432 

433 return sbs