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 zip 

2import numbers 

3import copy 

4import numpy as np 

5import warnings 

6import healpy as hp 

7from matplotlib import colors 

8from matplotlib import ticker 

9import matplotlib.pyplot as plt 

10import matplotlib.cm as cm 

11from matplotlib.ticker import FuncFormatter 

12import matplotlib as mpl 

13from matplotlib.patches import Ellipse 

14from matplotlib.collections import PatchCollection 

15 

16from lsst.sims.maf.utils import optimalBins, percentileClipping 

17from .plotHandler import BasePlotter, applyZPNorm 

18 

19from lsst.sims.utils import _equatorialFromGalactic, _healbin 

20from .perceptual_rainbow import makePRCmap 

21perceptual_rainbow = makePRCmap() 

22import numpy.ma as ma 

23 

24__all__ = ['setColorLims', 'setColorMap', 'HealpixSkyMap', 'HealpixPowerSpectrum', 

25 'HealpixHistogram', 'OpsimHistogram', 'BaseHistogram', 

26 'BaseSkyMap', 'HealpixSDSSSkyMap', 'LambertSkyMap'] 

27 

28baseDefaultPlotDict = {'title': None, 'xlabel': None, 'label': None, 

29 'logScale': False, 'percentileClip': None, 'normVal': None, 'zp': None, 

30 'cbarFormat': None, 'cmap': perceptual_rainbow, 'cbar_edge': True, 'nTicks': 10, 

31 'colorMin': None, 'colorMax': None, 

32 'xMin': None, 'xMax': None, 'yMin': None, 'yMax': None, 

33 'labelsize': None, 'fontsize': None, 'figsize': None, 'subplot': 111, 

34 'maskBelow': None} 

35 

36 

37def setColorLims(metricValue, plotDict): 

38 """Set up color bar limits.""" 

39 # Use plot dict if these values are set. 

40 colorMin = plotDict['colorMin'] 

41 colorMax = plotDict['colorMax'] 

42 # If not, try to use percentile clipping. 

43 if (plotDict['percentileClip'] is not None) & (np.size(metricValue.compressed()) > 0): 

44 pcMin, pcMax = percentileClipping(metricValue.compressed(), percentile=plotDict['percentileClip']) 

45 if colorMin is None: 

46 colorMin = pcMin 

47 if colorMax is None: 

48 colorMax = pcMax 

49 # If not, just use the data limits. 

50 if colorMin is None: 

51 colorMin = metricValue.compressed().min() 

52 if colorMax is None: 

53 colorMax = metricValue.compressed().max() 

54 # But make sure there is some range on the colorbar 

55 if colorMin == colorMax: 

56 colorMin = colorMin - 0.5 

57 colorMax = colorMax + 0.5 

58 return np.sort([colorMin, colorMax]) 

59 

60 

61def setColorMap(plotDict): 

62 cmap = plotDict['cmap'] 

63 if cmap is None: 

64 cmap = 'perceptual_rainbow' 

65 if type(cmap) == str: 

66 cmap = getattr(cm, cmap) 

67 # Set background and masked pixel colors default healpy white and gray. 

68 cmap = copy.copy(cmap) 

69 cmap.set_over(cmap(1.0)) 

70 cmap.set_under('w') 

71 cmap.set_bad('gray') 

72 return cmap 

73 

74 

75class HealpixSkyMap(BasePlotter): 

76 """ 

77 Generate a sky map of healpix metric values using healpy's mollweide view. 

78 """ 

79 def __init__(self): 

80 super(HealpixSkyMap, self).__init__() 

81 # Set the plotType 

82 self.plotType = 'SkyMap' 

83 self.objectPlotter = False 

84 # Set up the default plotting parameters. 

85 self.defaultPlotDict = {} 

86 self.defaultPlotDict.update(baseDefaultPlotDict) 

87 self.defaultPlotDict.update({'rot': (0, 0, 0), 'flip': 'astro', 'coord': 'C', 

88 'nside': 8, 'reduceFunc': np.mean}) 

89 # Note: for alt/az sky maps using the healpix plotter, you can use 

90 # {'rot': (90, 90, 90), 'flip': 'geo'} 

91 self.healpy_visufunc = hp.mollview 

92 self.healpy_visufunc_params = {} 

93 self.ax = None 

94 self.im = None 

95 

96 def __call__(self, metricValueIn, slicer, userPlotDict, fignum=None): 

97 """ 

98 Parameters 

99 ---------- 

100 metricValue : numpy.ma.MaskedArray 

101 slicer : lsst.sims.maf.slicers.HealpixSlicer 

102 userPlotDict: dict 

103 Dictionary of plot parameters set by user (overrides default values). 

104 fignum : int 

105 Matplotlib figure number to use (default = None, starts new figure). 

106 

107 Returns 

108 ------- 

109 int 

110 Matplotlib figure number used to create the plot. 

111 """ 

112 # Override the default plotting parameters with user specified values. 

113 plotDict = {} 

114 plotDict.update(self.defaultPlotDict) 

115 plotDict.update(userPlotDict) 

116 

117 # Check if we have a valid HEALpix slicer 

118 if 'Heal' in slicer.slicerName: 

119 # Update the metric data with zeropoint or normalization. 

120 metricValue = applyZPNorm(metricValueIn, plotDict) 

121 else: 

122 # Bin the values up on a healpix grid. 

123 metricValue = _healbin(slicer.slicePoints['ra'], slicer.slicePoints['dec'], 

124 metricValueIn.filled(slicer.badval), nside=plotDict['nside'], 

125 reduceFunc=plotDict['reduceFunc'], fillVal=slicer.badval) 

126 mask = np.zeros(metricValue.size) 

127 mask[np.where(metricValue == slicer.badval)] = 1 

128 metricValue = ma.array(metricValue, mask=mask) 

129 metricValue = applyZPNorm(metricValue, plotDict) 

130 

131 if plotDict['maskBelow'] is not None: 

132 toMask = np.where(metricValue <= plotDict['maskBelow'])[0] 

133 metricValue.mask[toMask] = True 

134 badval = hp.UNSEEN 

135 else: 

136 badval = slicer.badval 

137 

138 # Generate a Mollweide full-sky plot. 

139 fig = plt.figure(fignum, figsize=plotDict['figsize']) 

140 # Set up color bar limits. 

141 clims = setColorLims(metricValue, plotDict) 

142 cmap = setColorMap(plotDict) 

143 # Set log scale? 

144 norm = None 

145 if plotDict['logScale']: 

146 norm = 'log' 

147 # Avoid trying to log scale when zero is in the range. 

148 if (norm == 'log') & ((clims[0] <= 0 <= clims[1]) or (clims[0] >= 0 >= clims[1])): 

149 # Try something simple 

150 above = metricValue[np.where(metricValue > 0)] 

151 if len(above) > 0: 

152 clims[0] = above.max() 

153 # If still bad, give up and turn off norm 

154 if ((clims[0] <= 0 <= clims[1]) or (clims[0] >= 0 >= clims[1])): 

155 norm = None 

156 warnings.warn("Using norm was set to log, but color limits pass through 0. " 

157 "Adjusting so plotting doesn't fail") 

158 if plotDict['coord'] == 'C': 

159 notext = True 

160 else: 

161 notext = False 

162 

163 visufunc_params = {'title': plotDict['title'], 

164 'cbar': False, 

165 'min': clims[0], 

166 'max': clims[1], 

167 'rot': plotDict['rot'], 

168 'flip': plotDict['flip'], 

169 'coord': plotDict['coord'], 

170 'cmap': cmap, 

171 'norm': norm, 

172 'sub': plotDict['subplot'], 

173 'fig':fig.number, 

174 'notext': notext} 

175 visufunc_params.update(self.healpy_visufunc_params) 

176 self.healpy_visufunc(metricValue.filled(badval), **visufunc_params) 

177 

178 # Add a graticule (grid) over the globe. 

179 hp.graticule(dpar=30, dmer=30, verbose=False) 

180 # Add colorbar (not using healpy default colorbar because we want more tickmarks). 

181 self.ax = plt.gca() 

182 im = self.ax.get_images()[0] 

183 # Add label. 

184 if plotDict['label'] is not None: 

185 plt.figtext(0.8, 0.8, '%s' % (plotDict['label'])) 

186 # Make a color bar. Supress silly colorbar warnings. 

187 with warnings.catch_warnings(): 

188 warnings.simplefilter("ignore") 

189 cb = plt.colorbar(im, shrink=0.75, aspect=25, pad=0.1, orientation='horizontal', 

190 format=plotDict['cbarFormat'], extendrect=True) 

191 cb.set_label(plotDict['xlabel'], fontsize=plotDict['fontsize']) 

192 if plotDict['labelsize'] is not None: 

193 cb.ax.tick_params(labelsize=plotDict['labelsize']) 

194 if norm == 'log': 

195 tick_locator = ticker.LogLocator(numticks=plotDict['nTicks']) 

196 cb.locator = tick_locator 

197 cb.update_ticks() 

198 if (plotDict['nTicks'] is not None) & (norm != 'log'): 

199 tick_locator = ticker.MaxNLocator(nbins=plotDict['nTicks']) 

200 cb.locator = tick_locator 

201 cb.update_ticks() 

202 # If outputing to PDF, this fixes the colorbar white stripes 

203 if plotDict['cbar_edge']: 

204 cb.solids.set_edgecolor("face") 

205 return fig.number 

206 

207 

208class HealpixPowerSpectrum(BasePlotter): 

209 def __init__(self): 

210 self.plotType = 'PowerSpectrum' 

211 self.objectPlotter = False 

212 self.defaultPlotDict = {} 

213 self.defaultPlotDict.update(baseDefaultPlotDict) 

214 self.defaultPlotDict.update({'maxl': None, 'removeDipole': True, 'linestyle': '-'}) 

215 

216 def __call__(self, metricValue, slicer, userPlotDict, fignum=None): 

217 """ 

218 Generate and plot the power spectrum of metricValue (calculated on a healpix grid). 

219 """ 

220 if 'Healpix' not in slicer.slicerName: 

221 raise ValueError('HealpixPowerSpectrum for use with healpix metricBundles.') 

222 plotDict = {} 

223 plotDict.update(self.defaultPlotDict) 

224 plotDict.update(userPlotDict) 

225 

226 fig = plt.figure(fignum, figsize=plotDict['figsize']) 

227 ax = fig.add_subplot(plotDict['subplot']) 

228 # If the mask is True everywhere (no data), just plot zeros 

229 if False not in metricValue.mask: 

230 return None 

231 if plotDict['removeDipole']: 

232 cl = hp.anafast(hp.remove_dipole(metricValue.filled(slicer.badval)), lmax=plotDict['maxl']) 

233 else: 

234 cl = hp.anafast(metricValue.filled(slicer.badval), lmax=plotDict['maxl']) 

235 ell = np.arange(np.size(cl)) 

236 if plotDict['removeDipole']: 

237 condition = (ell > 1) 

238 else: 

239 condition = (ell > 0) 

240 ell = ell[condition] 

241 cl = cl[condition] 

242 # Plot the results. 

243 plt.plot(ell, (cl * ell * (ell + 1)) / 2.0 / np.pi, 

244 color=plotDict['color'], linestyle=plotDict['linestyle'], label=plotDict['label']) 

245 if cl.max() > 0 and plotDict['logScale']: 

246 plt.yscale('log') 

247 plt.xlabel(r'$l$', fontsize=plotDict['fontsize']) 

248 plt.ylabel(r'$l(l+1)C_l/(2\pi)$', fontsize=plotDict['fontsize']) 

249 if plotDict['labelsize'] is not None: 

250 plt.tick_params(axis='x', labelsize=plotDict['labelsize']) 

251 plt.tick_params(axis='y', labelsize=plotDict['labelsize']) 

252 if plotDict['title'] is not None: 

253 plt.title(plotDict['title']) 

254 # Return figure number (so we can reuse/add onto/save this figure if desired). 

255 return fig.number 

256 

257 

258class HealpixHistogram(BasePlotter): 

259 def __init__(self): 

260 self.plotType = 'Histogram' 

261 self.objectPlotter = False 

262 self.defaultPlotDict = {} 

263 self.defaultPlotDict.update(baseDefaultPlotDict) 

264 self.defaultPlotDict.update({'ylabel': 'Area (1000s of square degrees)', 

265 'bins': None, 'binsize': None, 'cumulative': False, 

266 'scale': None, 'linestyle': '-'}) 

267 self.baseHist = BaseHistogram() 

268 

269 def __call__(self, metricValue, slicer, userPlotDict, fignum=None): 

270 """ 

271 Histogram metricValue for all healpix points. 

272 """ 

273 if 'Healpix' not in slicer.slicerName: 

274 raise ValueError('HealpixHistogram is for use with healpix slicer.') 

275 plotDict = {} 

276 plotDict.update(self.defaultPlotDict) 

277 plotDict.update(userPlotDict) 

278 if plotDict['scale'] is None: 

279 plotDict['scale'] = (hp.nside2pixarea(slicer.nside, degrees=True) / 1000.0) 

280 fignum = self.baseHist(metricValue, slicer, plotDict, fignum=fignum) 

281 return fignum 

282 

283 

284class OpsimHistogram(BasePlotter): 

285 def __init__(self): 

286 self.plotType = 'Histogram' 

287 self.objectPlotter = False 

288 self.defaultPlotDict = {} 

289 self.defaultPlotDict.update(baseDefaultPlotDict) 

290 self.defaultPlotDict.update({'ylabel': 'Number of Fields', 'yaxisformat': '%d', 

291 'bins': None, 'binsize': None, 'cumulative': False, 

292 'scale': 1.0, 'linestyle': '-'}) 

293 self.baseHist = BaseHistogram() 

294 

295 def __call__(self, metricValue, slicer, userPlotDict, fignum=None): 

296 """ 

297 Histogram metricValue for all healpix points. 

298 """ 

299 if slicer.slicerName != 'OpsimFieldSlicer': 

300 raise ValueError('OpsimHistogram is for use with OpsimFieldSlicer.') 

301 plotDict = {} 

302 plotDict.update(self.defaultPlotDict) 

303 plotDict.update(userPlotDict) 

304 fignum = self.baseHist(metricValue, slicer, plotDict, fignum=fignum) 

305 return fignum 

306 

307 

308class BaseHistogram(BasePlotter): 

309 def __init__(self): 

310 self.plotType = 'Histogram' 

311 self.objectPlotter = False 

312 self.defaultPlotDict = {} 

313 self.defaultPlotDict.update(baseDefaultPlotDict) 

314 self.defaultPlotDict.update({'ylabel': 'Count', 'bins': None, 'binsize': None, 'cumulative': False, 

315 'scale': 1.0, 'yaxisformat': '%.3f', 'linestyle': '-'}) 

316 

317 def __call__(self, metricValueIn, slicer, userPlotDict, fignum=None): 

318 """ 

319 Plot a histogram of metricValues (such as would come from a spatial slicer). 

320 """ 

321 # Adjust metric values by zeropoint or normVal, and use 'compressed' version of masked array. 

322 plotDict = {} 

323 plotDict.update(self.defaultPlotDict) 

324 plotDict.update(userPlotDict) 

325 metricValue = applyZPNorm(metricValueIn, plotDict) 

326 metricValue = metricValue.compressed() 

327 # Toss any NaNs or infs 

328 metricValue = metricValue[np.isfinite(metricValue)] 

329 # Determine percentile clipped X range, if set. (and xmin/max not set). 

330 if plotDict['xMin'] is None and plotDict['xMax'] is None: 

331 if plotDict['percentileClip']: 

332 plotDict['xMin'], plotDict['xMax'] = percentileClipping(metricValue, 

333 percentile=plotDict['percentileClip']) 

334 # Set the histogram range values, to avoid cases of trying to histogram single-valued data. 

335 # First we try to use the range specified by a user, if there is one. Then use the data if not. 

336 # all of this only works if plotDict is not cumulative. 

337 histRange = [plotDict['xMin'], plotDict['xMax']] 

338 if histRange[0] is None: 

339 histRange[0] = metricValue.min() 

340 if histRange[1] is None: 

341 histRange[1] = metricValue.max() 

342 # Need to have some range of values on the histogram, or it will fail. 

343 if histRange[0] == histRange[1]: 

344 warnings.warn('Histogram range was single-valued; expanding default range.') 

345 histRange[1] = histRange[0] + 1.0 

346 # Set up the bins for the histogram. User specified 'bins' overrides 'binsize'. 

347 # Note that 'bins' could be a single number or an array, simply passed to plt.histogram. 

348 if plotDict['bins'] is not None: 

349 bins = plotDict['bins'] 

350 elif plotDict['binsize'] is not None: 

351 # If generating a cumulative histogram, want to use full range of data (but with given binsize). 

352 # .. but if user set histRange to be wider than full range of data, then 

353 # extend bins to cover this range, so we can make prettier plots. 

354 if plotDict['cumulative']: 

355 if plotDict['xMin'] is not None: 

356 # Potentially, expand the range for the cumulative histogram. 

357 bmin = np.min([metricValue.min(), plotDict['xMin']]) 

358 else: 

359 bmin = metricValue.min() 

360 if plotDict['xMax'] is not None: 

361 bmax = np.max([metricValue.max(), plotDict['xMax']]) 

362 else: 

363 bmax = metricValue.max() 

364 bins = np.arange(bmin, bmax + plotDict['binsize'] / 2.0, plotDict['binsize']) 

365 # Otherwise, not cumulative so just use metric values, without potential expansion. 

366 else: 

367 bins = np.arange(histRange[0], histRange[1] + plotDict['binsize'] / 2.0, plotDict['binsize']) 

368 # Catch edge-case where there is only 1 bin value 

369 if bins.size < 2: 

370 bins = np.arange(bins.min() - plotDict['binsize'] * 2.0, 

371 bins.max() + plotDict['binsize'] * 2.0, plotDict['binsize']) 

372 else: 

373 # If user did not specify bins or binsize, then we try to figure out a good number of bins. 

374 bins = optimalBins(metricValue) 

375 # Generate plots. 

376 fig = plt.figure(fignum, figsize=plotDict['figsize']) 

377 ax = fig.add_subplot(plotDict['subplot']) 

378 # Check if any data falls within histRange, because otherwise histogram generation will fail. 

379 if isinstance(bins, np.ndarray): 

380 condition = ((metricValue >= bins.min()) & (metricValue <= bins.max())) 

381 else: 

382 condition = ((metricValue >= histRange[0]) & (metricValue <= histRange[1])) 

383 plotValue = metricValue[condition] 

384 if len(plotValue) == 0: 

385 # No data is within histRange/bins. So let's just make a simple histogram anyway. 

386 n, b, p = plt.hist(metricValue, bins=50, histtype='step', cumulative=plotDict['cumulative'], 

387 log=plotDict['logScale'], label=plotDict['label'], 

388 color=plotDict['color']) 

389 else: 

390 # There is data to plot, and we've already ensured histRange/bins are more than single value. 

391 n, b, p = plt.hist(metricValue, bins=bins, range=histRange, 

392 histtype='step', log=plotDict['logScale'], 

393 cumulative=plotDict['cumulative'], 

394 label=plotDict['label'], color=plotDict['color']) 

395 hist_ylims = plt.ylim() 

396 if n.max() > hist_ylims[1]: 

397 plt.ylim(top = n.max()) 

398 if n.min() < hist_ylims[0] and not plotDict['logScale']: 

399 plt.ylim(bottom = n.min()) 

400 # Fill in axes labels and limits. 

401 # Option to use 'scale' to turn y axis into area or other value. 

402 

403 def mjrFormatter(y, pos): 

404 if not isinstance(plotDict['scale'], numbers.Number): 

405 raise ValueError('plotDict["scale"] must be a number to scale the y axis.') 

406 return plotDict['yaxisformat'] % (y * plotDict['scale']) 

407 

408 ax.yaxis.set_major_formatter(FuncFormatter(mjrFormatter)) 

409 # Set optional x, y limits. 

410 if 'xMin' in plotDict: 

411 plt.xlim(left=plotDict['xMin']) 

412 if 'xMax' in plotDict: 

413 plt.xlim(right=plotDict['xMax']) 

414 if 'yMin' in plotDict: 

415 plt.ylim(bottom=plotDict['yMin']) 

416 if 'yMax' in plotDict: 

417 plt.ylim(top=plotDict['yMax']) 

418 # Set/Add various labels. 

419 plt.xlabel(plotDict['xlabel'], fontsize=plotDict['fontsize']) 

420 plt.ylabel(plotDict['ylabel'], fontsize=plotDict['fontsize']) 

421 plt.title(plotDict['title']) 

422 if plotDict['labelsize'] is not None: 

423 plt.tick_params(axis='x', labelsize=plotDict['labelsize']) 

424 plt.tick_params(axis='y', labelsize=plotDict['labelsize']) 

425 # Return figure number 

426 return fig.number 

427 

428 

429class BaseSkyMap(BasePlotter): 

430 def __init__(self): 

431 self.plotType = 'SkyMap' 

432 self.objectPlotter = False # unless 'metricIsColor' is true.. 

433 self.defaultPlotDict = {} 

434 self.defaultPlotDict.update(baseDefaultPlotDict) 

435 self.defaultPlotDict.update({'projection': 'aitoff', 'radius': np.radians(1.75), 'alpha': 1.0, 

436 'plotMask': False, 'metricIsColor': False, 'cbar': True, 

437 'raCen': 0.0, 'mwZone': True, 'bgcolor': 'gray'}) 

438 

439 def _plot_tissot_ellipse(self, lon, lat, radius, ax=None, **kwargs): 

440 """Plot Tissot Ellipse/Tissot Indicatrix 

441 

442 Parameters 

443 ---------- 

444 lon : float or array_like 

445 longitude-like of ellipse centers (radians) 

446 lat : float or array_like 

447 latitude-like of ellipse centers (radians) 

448 radius : float or array_like 

449 radius of ellipses (radians) 

450 ax : Axes object (optional) 

451 matplotlib axes instance on which to draw ellipses. 

452 

453 Other Parameters 

454 ---------------- 

455 other keyword arguments will be passed to matplotlib.patches.Ellipse. 

456 

457 # The code in this method adapted from astroML, which is BSD-licensed. 

458 # See http: //github.com/astroML/astroML for details. 

459 """ 

460 # Code adapted from astroML, which is BSD-licensed. 

461 # See http: //github.com/astroML/astroML for details. 

462 ellipses = [] 

463 if ax is None: 

464 ax = plt.gca() 

465 for l, b, diam in np.broadcast(lon, lat, radius * 2.0): 

466 el = Ellipse((l, b), diam / np.cos(b), diam, **kwargs) 

467 ellipses.append(el) 

468 return ellipses 

469 

470 def _plot_ecliptic(self, raCen=0, ax=None): 

471 """ 

472 Plot a red line at location of ecliptic. 

473 """ 

474 if ax is None: 

475 ax = plt.gca() 

476 ecinc = 23.439291 * (np.pi / 180.0) 

477 ra_ec = np.arange(0, np.pi * 2., (np.pi * 2. / 360.)) 

478 dec_ec = np.sin(ra_ec) * ecinc 

479 lon = -(ra_ec - raCen - np.pi) % (np.pi * 2) - np.pi 

480 ax.plot(lon, dec_ec, 'r.', markersize=1.8, alpha=0.4) 

481 

482 def _plot_mwZone(self, raCen=0, peakWidth=np.radians(10.), taperLength=np.radians(80.), ax=None): 

483 """ 

484 Plot blue lines to mark the milky way galactic exclusion zone. 

485 """ 

486 if ax is None: 

487 ax = plt.gca() 

488 # Calculate galactic coordinates for mw location. 

489 step = 0.02 

490 galL = np.arange(-np.pi, np.pi + step / 2., step) 

491 val = peakWidth * np.cos(galL / taperLength * np.pi / 2.) 

492 galB1 = np.where(np.abs(galL) <= taperLength, val, 0) 

493 galB2 = np.where(np.abs(galL) <= taperLength, -val, 0) 

494 # Convert to ra/dec. 

495 # Convert to lon/lat and plot. 

496 ra, dec = _equatorialFromGalactic(galL, galB1) 

497 lon = -(ra - raCen - np.pi) % (np.pi * 2) - np.pi 

498 ax.plot(lon, dec, 'b.', markersize=1.8, alpha=0.4) 

499 ra, dec = _equatorialFromGalactic(galL, galB2) 

500 lon = -(ra - raCen - np.pi) % (np.pi * 2) - np.pi 

501 ax.plot(lon, dec, 'b.', markersize=1.8, alpha=0.4) 

502 

503 def __call__(self, metricValueIn, slicer, userPlotDict, fignum=None): 

504 """ 

505 Plot the sky map of metricValue for a generic spatial slicer. 

506 """ 

507 if 'ra' not in slicer.slicePoints or 'dec' not in slicer.slicePoints: 

508 errMessage = 'SpatialSlicer must contain "ra" and "dec" in slicePoints metadata.' 

509 errMessage += ' SlicePoints only contains keys %s.' % (slicer.slicePoints.keys()) 

510 raise ValueError(errMessage) 

511 plotDict = {} 

512 plotDict.update(self.defaultPlotDict) 

513 plotDict.update(userPlotDict) 

514 metricValue = applyZPNorm(metricValueIn, plotDict) 

515 

516 fig = plt.figure(fignum, figsize=plotDict['figsize']) 

517 # other projections available include 

518 # ['aitoff', 'hammer', 'lambert', 'mollweide', 'polar', 'rectilinear'] 

519 ax = fig.add_subplot(plotDict['subplot'], projection=plotDict['projection']) 

520 # Set up valid datapoints and colormin/max values. 

521 if plotDict['plotMask']: 

522 # Plot all data points. 

523 mask = np.ones(len(metricValue), dtype='bool') 

524 else: 

525 # Only plot points which are not masked. Flip numpy ma mask where 'False' == 'good'. 

526 good = ~metricValue.mask 

527 

528 # Add ellipses at RA/Dec locations - but don't add colors yet. 

529 lon = -(slicer.slicePoints['ra'][good] - plotDict['raCen'] - np.pi) % (np.pi * 2) - np.pi 

530 ellipses = self._plot_tissot_ellipse(lon, slicer.slicePoints['dec'][good], 

531 plotDict['radius'], rasterized=True, ax=ax) 

532 if plotDict['metricIsColor']: 

533 current = None 

534 for ellipse, mVal in zip(ellipses, metricValue.data[good]): 

535 if mVal[3] > 1: 

536 ellipse.set_alpha(1.0) 

537 ellipse.set_facecolor((mVal[0], mVal[1], mVal[2])) 

538 ellipse.set_edgecolor('k') 

539 current = ellipse 

540 else: 

541 ellipse.set_alpha(mVal[3]) 

542 ellipse.set_color((mVal[0], mVal[1], mVal[2])) 

543 ax.add_patch(ellipse) 

544 if current: 

545 ax.add_patch(current) 

546 else: 

547 # Determine color min/max values. metricValue.compressed = non-masked points. 

548 clims = setColorLims(metricValue, plotDict) 

549 # Determine whether or not to use auto-log scale. 

550 if plotDict['logScale'] == 'auto': 

551 if clims[0] > 0: 

552 if np.log10(clims[1]) - np.log10(clims[0]) > 3: 

553 plotDict['logScale'] = True 

554 else: 

555 plotDict['logScale'] = False 

556 else: 

557 plotDict['logScale'] = False 

558 if plotDict['logScale']: 

559 # Move min/max values to things that can be marked on the colorbar. 

560 #clims[0] = 10 ** (int(np.log10(clims[0]))) 

561 #clims[1] = 10 ** (int(np.log10(clims[1]))) 

562 norml = colors.LogNorm() 

563 p = PatchCollection(ellipses, cmap=plotDict['cmap'], alpha=plotDict['alpha'], 

564 linewidth=0, edgecolor=None, norm=norml, rasterized=True) 

565 else: 

566 p = PatchCollection(ellipses, cmap=plotDict['cmap'], alpha=plotDict['alpha'], 

567 linewidth=0, edgecolor=None, rasterized=True) 

568 p.set_array(metricValue.data[good]) 

569 p.set_clim(clims) 

570 ax.add_collection(p) 

571 # Add color bar (with optional setting of limits) 

572 if plotDict['cbar']: 

573 cb = plt.colorbar(p, aspect=25, extendrect=True, orientation='horizontal', 

574 format=plotDict['cbarFormat']) 

575 # If outputing to PDF, this fixes the colorbar white stripes 

576 if plotDict['cbar_edge']: 

577 cb.solids.set_edgecolor("face") 

578 cb.set_label(plotDict['xlabel'], fontsize=plotDict['fontsize']) 

579 cb.ax.tick_params(labelsize=plotDict['labelsize']) 

580 tick_locator = ticker.MaxNLocator(nbins=plotDict['nTicks']) 

581 cb.locator = tick_locator 

582 cb.update_ticks() 

583 # Add ecliptic 

584 self._plot_ecliptic(plotDict['raCen'], ax=ax) 

585 if plotDict['mwZone']: 

586 self._plot_mwZone(plotDict['raCen'], ax=ax) 

587 ax.grid(True, zorder=1) 

588 ax.xaxis.set_ticklabels([]) 

589 if plotDict['bgcolor'] is not None: 

590 ax.set_facecolor(plotDict['bgcolor']) 

591 # Add label. 

592 if plotDict['label'] is not None: 

593 plt.figtext(0.75, 0.9, '%s' % plotDict['label'], fontsize=plotDict['fontsize']) 

594 if plotDict['title'] is not None: 

595 plt.text(0.5, 1.09, plotDict['title'], horizontalalignment='center', 

596 transform=ax.transAxes, fontsize=plotDict['fontsize']) 

597 return fig.number 

598 

599 

600class HealpixSDSSSkyMap(BasePlotter): 

601 def __init__(self): 

602 self.plotType = 'SkyMap' 

603 self.objectPlotter = False 

604 self.defaultPlotDict = {} 

605 self.defaultPlotDict.update(baseDefaultPlotDict) 

606 self.defaultPlotDict.update({'cbarFormat': '%.2f', 

607 'raMin': -90, 'raMax': 90, 'raLen': 45, 

608 'decMin': -2., 'decMax': 2.}) 

609 

610 def __call__(self, metricValueIn, slicer, userPlotDict, fignum=None): 

611 

612 """ 

613 Plot the sky map of metricValue using healpy cartview plots in thin strips. 

614 raMin: Minimum RA to plot (deg) 

615 raMax: Max RA to plot (deg). Note raMin/raMax define the centers that will be plotted. 

616 raLen: Length of the plotted strips in degrees 

617 decMin: minimum dec value to plot 

618 decMax: max dec value to plot 

619 metricValueIn: metric values 

620 """ 

621 plotDict = {} 

622 plotDict.update(self.defaultPlotDict) 

623 plotDict.update(userPlotDict) 

624 metricValue = applyZPNorm(metricValueIn, plotDict) 

625 norm = None 

626 if plotDict['logScale']: 

627 norm = 'log' 

628 clims = setColorLims(metricValue, plotDict) 

629 cmap = setColorMap(plotDict) 

630 racenters = np.arange(plotDict['raMin'], plotDict['raMax'], plotDict['raLen']) 

631 nframes = racenters.size 

632 fig = plt.figure(fignum) 

633 # Do not specify or use plotDict['subplot'] because this is done in each call to hp.cartview. 

634 for i, racenter in enumerate(racenters): 

635 if i == 0: 

636 useTitle = plotDict['title'] + ' /n' + '%i < RA < %i' % (racenter - plotDict['raLen'], 

637 racenter + plotDict['raLen']) 

638 else: 

639 useTitle = '%i < RA < %i' % (racenter - plotDict['raLen'], racenter + plotDict['raLen']) 

640 hp.cartview(metricValue.filled(slicer.badval), title=useTitle, cbar=False, 

641 min=clims[0], max=clims[1], flip='astro', rot=(racenter, 0, 0), 

642 cmap=cmap, norm=norm, lonra=[-plotDict['raLen'], plotDict['raLen']], 

643 latra=[plotDict['decMin'], plotDict['decMax']], sub=(nframes + 1, 1, i + 1), fig=fig) 

644 hp.graticule(dpar=20, dmer=20, verbose=False) 

645 # Add colorbar (not using healpy default colorbar because want more tickmarks). 

646 ax = fig.add_axes([0.1, .15, .8, .075]) # left, bottom, width, height 

647 # Add label. 

648 if plotDict['label'] is not None: 

649 plt.figtext(0.8, 0.9, '%s' % plotDict['label']) 

650 # Make the colorbar as a seperate figure, 

651 # from http: //matplotlib.org/examples/api/colorbar_only.html 

652 cnorm = colors.Normalize(vmin=clims[0], vmax=clims[1]) 

653 # supress silly colorbar warnings 

654 with warnings.catch_warnings(): 

655 warnings.simplefilter("ignore") 

656 cb = mpl.colorbar.ColorbarBase(ax, cmap=cmap, norm=cnorm, 

657 orientation='horizontal', format=plotDict['cbarFormat']) 

658 cb.set_label(plotDict['xlabel']) 

659 cb.ax.tick_params(labelsize=plotDict['labelsize']) 

660 if norm == 'log': 

661 tick_locator = ticker.LogLocator(numticks=plotDict['nTicks']) 

662 cb.locator = tick_locator 

663 cb.update_ticks() 

664 if (plotDict['nTicks'] is not None) & (norm != 'log'): 

665 tick_locator = ticker.MaxNLocator(nbins=plotDict['nTicks']) 

666 cb.locator = tick_locator 

667 cb.update_ticks() 

668 # If outputing to PDF, this fixes the colorbar white stripes 

669 if plotDict['cbar_edge']: 

670 cb.solids.set_edgecolor("face") 

671 fig = plt.gcf() 

672 return fig.number 

673 

674 

675def project_lambert(longitude, latitude): 

676 """Project from RA,dec to plane 

677 https://en.wikipedia.org/wiki/Lambert_azimuthal_equal-area_projection 

678 """ 

679 

680 # flipping the sign on latitude goes north pole or south pole centered 

681 r_polar = 2*np.cos((np.pi/2+latitude)/2.) 

682 # Add pi/2 so north is up 

683 theta_polar = longitude + np.pi/2 

684 

685 x = r_polar * np.cos(theta_polar) 

686 y = r_polar * np.sin(theta_polar) 

687 return x, y 

688 

689 

690def draw_grat(ax): 

691 """Draw some graticule lines on an axis 

692 """ 

693 decs = np.radians(90.-np.array([20, 40, 60, 80])) 

694 ra = np.radians(np.arange(0, 361, 1)) 

695 for dec in decs: 

696 temp_dec = ra*0+dec 

697 x, y = project_lambert(ra, temp_dec) 

698 ax.plot(x, y, 'k--', alpha=0.5) 

699 

700 ras = np.radians(np.arange(0, 360+45, 45)) 

701 dec = np.radians(90.-np.arange(0, 81, 1)) 

702 for ra in ras: 

703 temp_ra = dec*0 + ra 

704 x, y = project_lambert(temp_ra, dec) 

705 ax.plot(x, y, 'k--', alpha=0.5) 

706 

707 for dec in decs: 

708 x, y = project_lambert(np.radians(45.), dec) 

709 ax.text(x, y, '%i' % np.round(np.degrees(dec))) 

710 

711 return ax 

712 

713class LambertSkyMap(BasePlotter): 

714 """ 

715 Use basemap and contour to make a Lambertian projection. 

716 Note that the plotDict can include a 'basemap' key with a dictionary of 

717 arbitrary kwargs to use with the call to Basemap. 

718 """ 

719 

720 def __init__(self): 

721 self.plotType = 'SkyMap' 

722 self.objectPlotter = False 

723 self.defaultPlotDict = {} 

724 self.defaultPlotDict.update(baseDefaultPlotDict) 

725 self.defaultPlotDict.update({'basemap': {'projection': 'nplaea', 'boundinglat': 1, 'lon_0': 180, 

726 'resolution': None, 'celestial': False, 'round': False}, 

727 'levels': 200, 'cbarFormat': '%i', 'norm': None}) 

728 

729 def __call__(self, metricValueIn, slicer, userPlotDict, fignum=None): 

730 

731 if 'ra' not in slicer.slicePoints or 'dec' not in slicer.slicePoints: 

732 errMessage = 'SpatialSlicer must contain "ra" and "dec" in slicePoints metadata.' 

733 errMessage += ' SlicePoints only contains keys %s.' % (slicer.slicePoints.keys()) 

734 raise ValueError(errMessage) 

735 

736 plotDict = {} 

737 plotDict.update(self.defaultPlotDict) 

738 plotDict.update(userPlotDict) 

739 

740 metricValue = applyZPNorm(metricValueIn, plotDict) 

741 clims = setColorLims(metricValue, plotDict) 

742 # Calculate the levels to use for the contour 

743 if np.size(plotDict['levels']) > 1: 

744 levels = plotDict['levels'] 

745 else: 

746 step = (clims[1] - clims[0]) / plotDict['levels'] 

747 levels = np.arange(clims[0], clims[1] + step, step) 

748 

749 fig = plt.figure(fignum, figsize=plotDict['figsize']) 

750 ax = fig.add_subplot(plotDict['subplot']) 

751 

752 x, y = project_lambert(slicer.slicePoints['ra'], slicer.slicePoints['dec']) 

753 # Contour the plot first to remove any anti-aliasing artifacts. Doesn't seem to work though. See: 

754 # http: //stackoverflow.com/questions/15822159/aliasing-when-saving-matplotlib\ 

755 # -filled-contour-plot-to-pdf-or-eps 

756 # tmpContour = m.contour(np.degrees(slicer.slicePoints['ra']), 

757 # np.degrees(slicer.slicePoints['dec']), 

758 # metricValue.filled(np.min(clims)-1), levels, tri=True, 

759 # cmap=plotDict['cmap'], ax=ax, latlon=True, 

760 # lw=1) 

761 

762 # Set masked values to be below the lowest contour level. 

763 if plotDict['norm'] == 'log': 

764 z_val = metricValue.filled(np.min(clims)-0.9) 

765 norm = colors.LogNorm(vmin=z_val.min(), vmax=z_val.max()) 

766 else: 

767 norm = plotDict['norm'] 

768 tcf = ax.tricontourf(x, y, metricValue.filled(np.min(clims)-0.9), levels, 

769 cmap=plotDict['cmap'], norm=norm) 

770 

771 ax = draw_grat(ax) 

772 

773 ax.set_xticks([]) 

774 ax.set_yticks([]) 

775 alt_limit = 10. 

776 x, y = project_lambert(0, np.radians(alt_limit)) 

777 max_val = np.max(np.abs([x, y])) 

778 ax.set_xlim([-max_val, max_val]) 

779 ax.set_ylim([-max_val, max_val]) 

780 

781 # Try to fix the ugly pdf contour problem 

782 for c in tcf.collections: 

783 c.set_edgecolor("face") 

784 

785 cb = plt.colorbar(tcf, format=plotDict['cbarFormat']) 

786 cb.set_label(plotDict['xlabel']) 

787 if plotDict['labelsize'] is not None: 

788 cb.ax.tick_params(labelsize=plotDict['labelsize']) 

789 # Pop in an extra line to raise the title a bit 

790 ax.set_title(plotDict['title']+'\n ') 

791 # If outputing to PDF, this fixes the colorbar white stripes 

792 if plotDict['cbar_edge']: 

793 cb.solids.set_edgecolor("face") 

794 return fig.number