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 print_function 

2from builtins import zip 

3from builtins import range 

4from builtins import object 

5import numpy as np 

6import os 

7import glob 

8import healpy as hp 

9from lsst.sims.photUtils import Sed, Bandpass 

10from lsst.sims.skybrightness.twilightFunc import twilightFunc 

11from scipy.interpolate import InterpolatedUnivariateSpline, interp1d 

12from lsst.utils import getPackageDir 

13 

14# Make backwards compatible with healpy 

15if hasattr(hp, 'get_interp_weights'): 15 ↛ 17line 15 didn't jump to line 17, because the condition on line 15 was never false

16 get_neighbours = hp.get_interp_weights 

17elif hasattr(hp, 'get_neighbours'): 

18 get_neighbours = hp.get_neighbours 

19else: 

20 print("Could not find appropriate healpy function for get_interp_weight or get_neighbours") 

21 

22 

23__all__ = ['id2intid', 'intid2id', 'loadSpecFiles', 'BaseSingleInterp', 'ScatteredStar', 'LowerAtm', 

24 'UpperAtm', 'MergedSpec', 'Airglow', 'TwilightInterp', 'MoonInterp', 

25 'ZodiacalInterp'] 

26 

27 

28def id2intid(ids): 

29 """ 

30 take an array of ids, and convert them to an integer id. 

31 Handy if you want to put things into a sparse array. 

32 """ 

33 uids = np.unique(ids) 

34 order = np.argsort(ids) 

35 oids = ids[order] 

36 uintids = np.arange(np.size(uids), dtype=int) 

37 left = np.searchsorted(oids, uids) 

38 right = np.searchsorted(oids, uids, side='right') 

39 intids = np.empty(ids.size, dtype=int) 

40 for i in range(np.size(left)): 

41 intids[left[i]:right[i]] = uintids[i] 

42 result = intids*0 

43 result[order] = intids 

44 return result, uids, uintids 

45 

46 

47def intid2id(intids, uintids, uids, dtype=int): 

48 """ 

49 convert an int back to an id 

50 """ 

51 ids = np.zeros(np.size(intids)) 

52 

53 order = np.argsort(intids) 

54 ointids = intids[order] 

55 left = np.searchsorted(ointids, uintids, side='left') 

56 right = np.searchsorted(ointids, uintids, side='right') 

57 for i, (le, ri) in enumerate(zip(left, right)): 

58 ids[le:ri] = uids[i] 

59 result = np.zeros(np.size(intids), dtype=dtype) 

60 result[order] = ids 

61 

62 return result 

63 

64 

65def loadSpecFiles(filenames, mags=False): 

66 """ 

67 Load up the ESO spectra. 

68 

69 The ESO npz files contain the following arrays: 

70 filterWave: The central wavelengths of the pre-computed magnitudes 

71 wave: wavelengths for the spectra 

72 spec: array of spectra and magnitudes along with the relevant variable inputs. For example, 

73 airglow has dtype = [('airmass', '<f8'), ('solarFlux', '<f8'), ('spectra', '<f8', (17001,)), 

74 ('mags', '<f8', (6,)] 

75 For each unique airmass and solarFlux value, there is a 17001 elements spectra and 6 magnitudes. 

76 """ 

77 

78 if len(filenames) == 1: 

79 temp = np.load(filenames[0]) 

80 wave = temp['wave'].copy() 

81 filterWave = temp['filterWave'].copy() 

82 if mags: 

83 # don't copy the spectra to save memory space 

84 dt = np.dtype([(key, temp['spec'].dtype[i]) for 

85 i, key in enumerate(temp['spec'].dtype.names) if key != 'spectra']) 

86 spec = np.zeros(temp['spec'].size, dtype=dt) 

87 for key in temp['spec'].dtype.names: 

88 if key != 'spectra': 

89 spec[key] = temp['spec'][key].copy() 

90 else: 

91 spec = temp['spec'].copy() 

92 else: 

93 temp = np.load(filenames[0]) 

94 wave = temp['wave'].copy() 

95 filterWave = temp['filterWave'].copy() 

96 if mags: 

97 # don't copy the spectra to save memory space 

98 dt = np.dtype([(key, temp['spec'].dtype[i]) for 

99 i, key in enumerate(temp['spec'].dtype.names) if key != 'spectra']) 

100 spec = np.zeros(temp['spec'].size, dtype=dt) 

101 for key in temp['spec'].dtype.names: 

102 if key != 'spectra': 

103 spec[key] = temp['spec'][key].copy() 

104 else: 

105 spec = temp['spec'].copy() 

106 for filename in filenames[1:]: 

107 temp = np.load(filename) 

108 if mags: 

109 # don't copy the spectra to save memory space 

110 dt = np.dtype([(key, temp['spec'].dtype[i]) for 

111 i, key in enumerate(temp['spec'].dtype.names) if key != 'spectra']) 

112 tempspec = np.zeros(temp['spec'].size, dtype=dt) 

113 for key in temp['spec'].dtype.names: 

114 if key != 'spectra': 

115 tempspec[key] = temp['spec'][key].copy() 

116 else: 

117 tempspec = temp['spec'] 

118 spec = np.append(spec, tempspec) 

119 return spec, wave, filterWave 

120 

121 

122class BaseSingleInterp(object): 

123 """ 

124 Base class for sky components that only need to be interpolated on airmass 

125 """ 

126 

127 def __init__(self, compName=None, sortedOrder=['airmass', 'nightTimes'], mags=False): 

128 """ 

129 mags: Rather than the full spectrum, return the LSST ugrizy magnitudes. 

130 """ 

131 

132 self.mags = mags 

133 

134 dataDir = os.path.join(getPackageDir('sims_skybrightness_data'), 'ESO_Spectra/'+compName) 

135 

136 filenames = sorted(glob.glob(dataDir+'/*.npz')) 

137 self.spec, self.wave, self.filterWave = loadSpecFiles(filenames, mags=self.mags) 

138 

139 # Take the log of the spectra in case we want to interp in log space. 

140 if not mags: 

141 self.logSpec = np.zeros(self.spec['spectra'].shape, dtype=float) 

142 good = np.where(self.spec['spectra'] != 0) 

143 self.logSpec[good] = np.log10(self.spec['spectra'][good]) 

144 self.specSize = self.spec['spectra'][0].size 

145 else: 

146 self.specSize = 0 

147 

148 # What order are the dimesions sorted by (from how the .npz was packaged) 

149 self.sortedOrder = sortedOrder 

150 self.dimDict = {} 

151 self.dimSizes = {} 

152 for dt in self.sortedOrder: 

153 self.dimDict[dt] = np.unique(self.spec[dt]) 

154 self.dimSizes[dt] = np.size(np.unique(self.spec[dt])) 

155 

156 # Set up and save the dict to order the filters once. 

157 self.filterNameDict = {'u': 0, 'g': 1, 'r': 2, 'i': 3, 'z': 4, 'y': 5} 

158 

159 def __call__(self, intepPoints, filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

160 if self.mags: 

161 return self.interpMag(intepPoints, filterNames=filterNames) 

162 else: 

163 return self.interpSpec(intepPoints) 

164 

165 def indxAndWeights(self, points, grid): 

166 """ 

167 for given 1-D points, find the grid points on either side and return the weights 

168 assume grid is sorted 

169 """ 

170 

171 order = np.argsort(points) 

172 

173 indxL = np.empty(points.size, dtype=int) 

174 indxR = np.empty(points.size, dtype=int) 

175 

176 indxR[order] = np.searchsorted(grid, points[order]) 

177 indxL = indxR-1 

178 

179 # If points off the grid were requested, just use the edge grid point 

180 offGrid = np.where(indxR == grid.size) 

181 indxR[offGrid] = grid.size-1 

182 fullRange = grid[indxR]-grid[indxL] 

183 

184 wL = np.zeros(fullRange.size, dtype=float) 

185 wR = np.ones(fullRange.size, dtype=float) 

186 

187 good = np.where(fullRange != 0) 

188 wL[good] = (grid[indxR][good] - points[good])/fullRange[good] 

189 wR[good] = (points[good] - grid[indxL[good]])/fullRange[good] 

190 

191 return indxR, indxL, wR, wL 

192 

193 def _weighting(self, interpPoints, values): 

194 """ 

195 given a list/array of airmass values, return a dict with the interpolated 

196 spectrum at each airmass and the wavelength array. 

197 

198 Input interpPoints should be sorted 

199 """ 

200 results = np.zeros((interpPoints.size, np.size(values[0])), dtype=float) 

201 

202 inRange = np.where((interpPoints['airmass'] <= np.max(self.dimDict['airmass'])) & 

203 (interpPoints['airmass'] >= np.min(self.dimDict['airmass']))) 

204 indxR, indxL, wR, wL = self.indxAndWeights(interpPoints['airmass'][inRange], 

205 self.dimDict['airmass']) 

206 

207 nextra = 1 

208 

209 # XXX--should I use the log spectra? Make a check and switch back and forth? 

210 results[inRange] = wR[:, np.newaxis]*values[indxR*nextra] + \ 

211 wL[:, np.newaxis]*values[indxL*nextra] 

212 

213 return results 

214 

215 def interpSpec(self, interpPoints): 

216 result = self._weighting(interpPoints, self.logSpec) 

217 mask = np.where(result == 0.) 

218 result = 10.**result 

219 result[mask] = 0. 

220 return {'spec': result, 'wave': self.wave} 

221 

222 def interpMag(self, interpPoints, filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

223 filterindx = [self.filterNameDict[key] for key in filterNames] 

224 result = self._weighting(interpPoints, self.spec['mags'][:, filterindx]) 

225 mask = np.where(result == 0.) 

226 result = 10.**(-0.4*(result-np.log10(3631.))) 

227 result[mask] = 0. 

228 return {'spec': result, 'wave': self.filterWave} 

229 

230 

231class ScatteredStar(BaseSingleInterp): 

232 """ 

233 Interpolate the spectra caused by scattered starlight. 

234 """ 

235 

236 def __init__(self, compName='ScatteredStarLight', mags=False): 

237 super(ScatteredStar, self).__init__(compName=compName, mags=mags) 

238 

239 

240class LowerAtm(BaseSingleInterp): 

241 """ 

242 Interpolate the spectra caused by the lower atmosphere. 

243 """ 

244 

245 def __init__(self, compName='LowerAtm', mags=False): 

246 super(LowerAtm, self).__init__(compName=compName, mags=mags) 

247 

248 

249class UpperAtm(BaseSingleInterp): 

250 """ 

251 Interpolate the spectra caused by the upper atmosphere. 

252 """ 

253 

254 def __init__(self, compName='UpperAtm', mags=False): 

255 super(UpperAtm, self).__init__(compName=compName, mags=mags) 

256 

257 

258class MergedSpec(BaseSingleInterp): 

259 """ 

260 Interpolate the spectra caused by the sum of the scattered starlight, airglow, upper and lower atmosphere. 

261 """ 

262 

263 def __init__(self, compName='MergedSpec', mags=False): 

264 super(MergedSpec, self).__init__(compName=compName, mags=mags) 

265 

266 

267class Airglow(BaseSingleInterp): 

268 """ 

269 Interpolate the spectra caused by airglow. 

270 """ 

271 

272 def __init__(self, compName='Airglow', sortedOrder=['airmass', 'solarFlux'], mags=False): 

273 super(Airglow, self).__init__(compName=compName, mags=mags, sortedOrder=sortedOrder) 

274 self.nSolarFlux = np.size(self.dimDict['solarFlux']) 

275 

276 def _weighting(self, interpPoints, values): 

277 

278 results = np.zeros((interpPoints.size, np.size(values[0])), dtype=float) 

279 # Only interpolate point that lie in the model grid 

280 inRange = np.where((interpPoints['airmass'] <= np.max(self.dimDict['airmass'])) & 

281 (interpPoints['airmass'] >= np.min(self.dimDict['airmass'])) & 

282 (interpPoints['solarFlux'] >= np.min(self.dimDict['solarFlux'])) & 

283 (interpPoints['solarFlux'] <= np.max(self.dimDict['solarFlux']))) 

284 usePoints = interpPoints[inRange] 

285 amRightIndex, amLeftIndex, amRightW, amLeftW = self.indxAndWeights(usePoints['airmass'], 

286 self.dimDict['airmass']) 

287 

288 sfRightIndex, sfLeftIndex, sfRightW, sfLeftW = self.indxAndWeights(usePoints['solarFlux'], 

289 self.dimDict['solarFlux']) 

290 

291 for amIndex, amW in zip([amRightIndex, amLeftIndex], [amRightW, amLeftW]): 

292 for sfIndex, sfW in zip([sfRightIndex, sfLeftIndex], [sfRightW, sfLeftW]): 

293 results[inRange] += amW[:, np.newaxis]*sfW[:, np.newaxis] * \ 

294 values[amIndex*self.nSolarFlux+sfIndex] 

295 return results 

296 

297 

298class TwilightInterp(object): 

299 

300 def __init__(self, mags=False, darkSkyMags=None, fitResults=None): 

301 """ 

302 Read the Solar spectrum into a handy object and compute mags in different filters 

303 mags: If true, only return the LSST filter magnitudes, otherwise return the full spectrum 

304 

305 darkSkyMags = dict of the zenith dark sky values to be assumed. The twilight fits are 

306 done relative to the dark sky level. 

307 fitResults = dict of twilight parameters based on twilightFunc. Keys should be filter names. 

308 """ 

309 

310 if darkSkyMags is None: 

311 darkSkyMags = {'u': 22.8, 'g': 22.3, 'r': 21.2, 

312 'i': 20.3, 'z': 19.3, 'y': 18.0, 

313 'B': 22.35, 'G': 21.71, 'R': 21.3} 

314 

315 self.mags = mags 

316 

317 dataDir = getPackageDir('sims_skybrightness_data') 

318 

319 solarSaved = np.load(os.path.join(dataDir, 'solarSpec/solarSpec.npz')) 

320 self.solarSpec = Sed(wavelen=solarSaved['wave'], flambda=solarSaved['spec']) 

321 solarSaved.close() 

322 

323 canonFilters = {} 

324 fnames = ['blue_canon.csv', 'green_canon.csv', 'red_canon.csv'] 

325 

326 # Filter names, from bluest to reddest. 

327 self.filterNames = ['B', 'G', 'R'] 

328 

329 for fname, filterName in zip(fnames, self.filterNames): 

330 bpdata = np.genfromtxt(os.path.join(dataDir, 'Canon/', fname), delimiter=', ', 

331 dtype=list(zip(['wave', 'through'], [float]*2))) 

332 bpTemp = Bandpass() 

333 bpTemp.setBandpass(bpdata['wave'], bpdata['through']) 

334 canonFilters[filterName] = bpTemp 

335 

336 # Tack on the LSST filters 

337 throughPath = os.path.join(getPackageDir('throughputs'), 'baseline') 

338 lsstKeys = ['u', 'g', 'r', 'i', 'z', 'y'] 

339 for key in lsstKeys: 

340 bp = np.loadtxt(os.path.join(throughPath, 'total_'+key+'.dat'), 

341 dtype=list(zip(['wave', 'trans'], [float]*2))) 

342 tempB = Bandpass() 

343 tempB.setBandpass(bp['wave'], bp['trans']) 

344 canonFilters[key] = tempB 

345 self.filterNames.append(key) 

346 

347 # MAGIC NUMBERS from fitting the all-sky camera: 

348 # Code to generate values in sims_skybrightness/examples/fitTwiSlopesSimul.py 

349 # Which in turn uses twilight maps from sims_skybrightness/examples/buildTwilMaps.py 

350 # values are of the form: 

351 # 0: ratio of f^z_12 to f_dark^z 

352 # 1: slope of curve wrt sun alt 

353 # 2: airmass term (10^(arg[2]*(X-1))) 

354 # 3: azimuth term. 

355 # 4: zenith dark sky flux (erg/s/cm^2) 

356 

357 # For z and y, just assuming the shape parameter fits are similar to the other bands. 

358 # Looks like the diode is not sensitive enough to detect faint sky. 

359 # Using the Patat et al 2006 I-band values for z and modeified a little for y as a temp fix. 

360 if fitResults is None: 

361 self.fitResults = {'B': [7.56765633e+00, 2.29798055e+01, 2.86879956e-01, 

362 3.01162143e-01, 2.58462036e-04], 

363 'G': [2.38561156e+00, 2.29310648e+01, 2.97733083e-01, 

364 3.16403197e-01, 7.29660095e-04], 

365 'R': [1.75498017e+00, 2.22011802e+01, 2.98619033e-01, 

366 3.28880254e-01, 3.24411056e-04], 

367 'z': [2.29, 24.08, 0.3, 0.3, -666], 

368 'y': [2.0, 24.08, 0.3, 0.3, -666]} 

369 

370 # XXX-completely arbitrary fudge factor to make things brighter in the blue 

371 # Just copy the blue and say it's brighter. 

372 self.fitResults['u'] = [16., 2.29622121e+01, 2.85862729e-01, 

373 2.99902574e-01, 2.32325117e-04] 

374 else: 

375 self.fitResults = fitResults 

376 

377 # Take out any filters that don't have fit results 

378 self.filterNames = [key for key in self.filterNames if key in self.fitResults] 

379 

380 self.effWave = [] 

381 self.solarMag = [] 

382 for filterName in self.filterNames: 

383 self.effWave.append(canonFilters[filterName].calcEffWavelen()[0]) 

384 self.solarMag.append(self.solarSpec.calcMag(canonFilters[filterName])) 

385 

386 ord = np.argsort(self.effWave) 

387 self.filterNames = np.array(self.filterNames)[ord] 

388 self.effWave = np.array(self.effWave)[ord] 

389 self.solarMag = np.array(self.solarMag)[ord] 

390 

391 # update the fit results to be zeropointed properly 

392 for key in self.fitResults: 

393 f0 = 10.**(-0.4*(darkSkyMags[key]-np.log10(3631.))) 

394 self.fitResults[key][-1] = f0 

395 

396 self.solarWave = self.solarSpec.wavelen 

397 self.solarFlux = self.solarSpec.flambda 

398 # This one isn't as bad as the model grids, maybe we could get away with computing the magnitudes 

399 # in the __call__ each time. 

400 if mags: 

401 # Load up the LSST filters and convert the solarSpec.flabda and solarSpec.wavelen to fluxes 

402 throughPath = os.path.join(getPackageDir('throughputs'), 'baseline') 

403 self.lsstFilterNames = ['u', 'g', 'r', 'i', 'z', 'y'] 

404 self.lsstEquations = np.zeros((np.size(self.lsstFilterNames), 

405 np.size(self.fitResults['B'])), dtype=float) 

406 self.lsstEffWave = [] 

407 

408 fits = np.empty((np.size(self.effWave), np.size(self.fitResults['B'])), dtype=float) 

409 for i, fn in enumerate(self.filterNames): 

410 fits[i, :] = self.fitResults[fn] 

411 

412 for filtername in self.lsstFilterNames: 

413 bp = np.loadtxt(os.path.join(throughPath, 'total_'+filtername+'.dat'), 

414 dtype=list(zip(['wave', 'trans'], [float]*2))) 

415 tempB = Bandpass() 

416 tempB.setBandpass(bp['wave'], bp['trans']) 

417 self.lsstEffWave.append(tempB.calcEffWavelen()[0]) 

418 # Loop through the parameters and interpolate to new eff wavelengths 

419 for i in np.arange(self.lsstEquations[0, :].size): 

420 interp = InterpolatedUnivariateSpline(self.effWave, fits[:, i]) 

421 self.lsstEquations[:, i] = interp(self.lsstEffWave) 

422 # Set the dark sky flux 

423 for i, filterName in enumerate(self.lsstFilterNames): 

424 self.lsstEquations[i, -1] = 10.**(-0.4*(darkSkyMags[filterName]-np.log10(3631.))) 

425 

426 self.filterNameDict = {'u': 0, 'g': 1, 'r': 2, 'i': 3, 'z': 4, 'y': 5} 

427 

428 def printFitsUsed(self): 

429 """ 

430 Print out the fit parameters being used 

431 """ 

432 print('\\tablehead{\colhead{Filter} & \colhead{$r_{12/z}$} & \colhead{$a$ (1/radians)} & \colhead{$b$ (1/airmass)} & \colhead{$c$ (az term/airmass)} & \colhead{$f_z_dark$ (erg/s/cm$^2$)$\\times 10^8$} & \colhead{m$_z_dark$}}') 

433 for key in self.fitResults: 

434 numbers = '' 

435 for num in self.fitResults[key]: 

436 if num > .001: 

437 numbers += ' & %.2f' % num 

438 else: 

439 numbers += ' & %.2f' % (num*1e8) 

440 print(key, numbers, ' & ', '%.2f' % (-2.5*np.log10(self.fitResults[key][-1])+np.log10(3631.))) 

441 

442 def __call__(self, intepPoints, filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

443 if self.mags: 

444 return self.interpMag(intepPoints, filterNames=filterNames) 

445 else: 

446 return self.interpSpec(intepPoints) 

447 

448 def interpMag(self, interpPoints, maxAM=2.5, 

449 limits=[np.radians(-5.), np.radians(-20.)], 

450 filterNames=['u', 'g', 'r', 'i', 'z', 'y']): 

451 """ 

452 Originally fit the twilight with a cutoff of sun altitude of -11 degrees. I think it can be safely 

453 extrapolated farther, but be warned you may be entering a regime where it breaks down. 

454 """ 

455 npts = len(filterNames) 

456 result = np.zeros((np.size(interpPoints), npts), dtype=float) 

457 

458 good = np.where((interpPoints['sunAlt'] >= np.min(limits)) & 

459 (interpPoints['sunAlt'] <= np.max(limits)) & 

460 (interpPoints['airmass'] <= maxAM) & 

461 (interpPoints['airmass'] >= 1.))[0] 

462 

463 for i, filterName in enumerate(filterNames): 

464 result[good, i] = twilightFunc(interpPoints[good], 

465 *self.lsstEquations[self.filterNameDict[filterName], :].tolist()) 

466 

467 return {'spec': result, 'wave': self.lsstEffWave} 

468 

469 def interpSpec(self, interpPoints, maxAM=2.5, 

470 limits=[np.radians(-5.), np.radians(-20.)]): 

471 """ 

472 interpPoints should have airmass, azRelSun, and sunAlt. 

473 """ 

474 

475 npts = np.size(self.solarWave) 

476 result = np.zeros((np.size(interpPoints), npts), dtype=float) 

477 

478 good = np.where((interpPoints['sunAlt'] >= np.min(limits)) & 

479 (interpPoints['sunAlt'] <= np.max(limits)) & 

480 (interpPoints['airmass'] <= maxAM) & 

481 (interpPoints['airmass'] >= 1.))[0] 

482 

483 # Compute the expected flux in each of the filters that we have fits for 

484 fluxes = [] 

485 for filterName in self.filterNames: 

486 fluxes.append(twilightFunc(interpPoints[good], *self.fitResults[filterName])) 

487 fluxes = np.array(fluxes) 

488 

489 # ratio of model flux to raw solar flux: 

490 yvals = fluxes.T/(10.**(-0.4*(self.solarMag-np.log10(3631.)))) 

491 

492 # Find wavelengths bluer than cutoff 

493 blueRegion = np.where(self.solarWave < np.min(self.effWave)) 

494 

495 for i, yval in enumerate(yvals): 

496 interpF = interp1d(self.effWave, yval, bounds_error=False, fill_value=yval[-1]) 

497 ratio = interpF(self.solarWave) 

498 interpBlue = InterpolatedUnivariateSpline(self.effWave, yval, k=1) 

499 ratio[blueRegion] = interpBlue(self.solarWave[blueRegion]) 

500 result[good[i]] = self.solarFlux*ratio 

501 

502 return {'spec': result, 'wave': self.solarWave} 

503 

504 

505class MoonInterp(BaseSingleInterp): 

506 """ 

507 Read in the saved Lunar spectra and interpolate. 

508 """ 

509 

510 def __init__(self, compName='Moon', sortedOrder=['moonSunSep', 'moonAltitude', 'hpid'], mags=False): 

511 super(MoonInterp, self).__init__(compName=compName, sortedOrder=sortedOrder, mags=mags) 

512 # Magic number from when the templates were generated 

513 self.nside = 4 

514 

515 def _weighting(self, interpPoints, values): 

516 """ 

517 Weighting for the scattered moonlight. 

518 """ 

519 

520 result = np.zeros((interpPoints.size, np.size(values[0])), dtype=float) 

521 

522 # Check that moonAltitude is in range, otherwise return zero array 

523 if np.max(interpPoints['moonAltitude']) < np.min(self.dimDict['moonAltitude']): 

524 return result 

525 

526 # Find the neighboring healpixels 

527 hpids, hweights = get_neighbours(self.nside, np.pi/2.-interpPoints['alt'], 

528 interpPoints['azRelMoon']) 

529 

530 badhp = np.in1d(hpids.ravel(), self.dimDict['hpid'], invert=True).reshape(hpids.shape) 

531 hweights[badhp] = 0. 

532 

533 norm = np.sum(hweights, axis=0) 

534 good = np.where(norm != 0.)[0] 

535 hweights[:, good] = hweights[:, good]/norm[good] 

536 

537 # Find the neighboring moonAltitude points in the grid 

538 rightMAs, leftMAs, maRightW, maLeftW = self.indxAndWeights(interpPoints['moonAltitude'], 

539 self.dimDict['moonAltitude']) 

540 

541 # Find the neighboring moonSunSep points in the grid 

542 rightMss, leftMss, mssRightW, mssLeftW = self.indxAndWeights(interpPoints['moonSunSep'], 

543 self.dimDict['moonSunSep']) 

544 

545 nhpid = self.dimDict['hpid'].size 

546 nMA = self.dimDict['moonAltitude'].size 

547 # Convert the hpid to an index. 

548 tmp = intid2id(hpids.ravel(), self.dimDict['hpid'], 

549 np.arange(self.dimDict['hpid'].size)) 

550 hpindx = tmp.reshape(hpids.shape) 

551 # loop though the hweights and the moonAltitude weights 

552 

553 for hpid, hweight in zip(hpindx, hweights): 

554 for maid, maW in zip([rightMAs, leftMAs], [maRightW, maLeftW]): 

555 for mssid, mssW in zip([rightMss, leftMss], [mssRightW, mssLeftW]): 

556 weight = hweight*maW*mssW 

557 result += weight[:, np.newaxis]*values[mssid*nhpid*nMA+maid*nhpid+hpid] 

558 

559 return result 

560 

561 

562class ZodiacalInterp(BaseSingleInterp): 

563 """ 

564 Interpolate the zodiacal light based on the airmass and the healpix ID where 

565 the healpixels are in ecliptic coordinates, with the sun at ecliptic longitude zero 

566 """ 

567 

568 def __init__(self, compName='Zodiacal', sortedOrder=['airmass', 'hpid'], mags=False): 

569 super(ZodiacalInterp, self).__init__(compName=compName, sortedOrder=sortedOrder, mags=mags) 

570 self.nside = hp.npix2nside(np.size(np.where(self.spec['airmass'] == 

571 np.unique(self.spec['airmass'])[0])[0])) 

572 

573 def _weighting(self, interpPoints, values): 

574 """ 

575 interpPoints is a numpy array where interpolation is desired 

576 values are the model values. 

577 """ 

578 result = np.zeros((interpPoints.size, np.size(values[0])), dtype=float) 

579 

580 inRange = np.where((interpPoints['airmass'] <= np.max(self.dimDict['airmass'])) & 

581 (interpPoints['airmass'] >= np.min(self.dimDict['airmass']))) 

582 usePoints = interpPoints[inRange] 

583 # Find the neighboring healpixels 

584 hpids, hweights = get_neighbours(self.nside, np.pi/2.-usePoints['altEclip'], 

585 usePoints['azEclipRelSun']) 

586 

587 badhp = np.in1d(hpids.ravel(), self.dimDict['hpid'], invert=True).reshape(hpids.shape) 

588 hweights[badhp] = 0. 

589 

590 norm = np.sum(hweights, axis=0) 

591 good = np.where(norm != 0.)[0] 

592 hweights[:, good] = hweights[:, good]/norm[good] 

593 

594 amRightIndex, amLeftIndex, amRightW, amLeftW = self.indxAndWeights(usePoints['airmass'], 

595 self.dimDict['airmass']) 

596 

597 nhpid = self.dimDict['hpid'].size 

598 # loop though the hweights and the airmass weights 

599 for hpid, hweight in zip(hpids, hweights): 

600 for amIndex, amW in zip([amRightIndex, amLeftIndex], [amRightW, amLeftW]): 

601 weight = hweight*amW 

602 result[inRange] += weight[:, np.newaxis]*values[amIndex*nhpid+hpid] 

603 

604 return result