LSSTApplications  16.0+12,16.0+13,16.0-1-g16dcc35+1,16.0-1-g44ebaf5+1,16.0-1-g4515a79+10,16.0-1-g928b041+1,16.0-1-g943889d+1,16.0-1-g9884240+1,16.0-1-g98efed3+9,16.0-1-gce273f5+1,16.0-19-g2da375352+1,16.0-2-g568a347+8,16.0-2-g839ba83+2,16.0-2-g8827b2c+1,16.0-2-g934b51e+2,16.0-2-g9d5294e+1,16.0-3-g06ec41f+2,16.0-3-g30bbe86+2,16.0-3-g3806c63+1,16.0-3-g5dc86b7+2,16.0-3-gc6a11d1+2,16.0-4-g4dc7b79+2,16.0-4-g50d071e+1,16.0-4-gd6aeece+1,16.0-4-gedb2a3a+2,16.0-5-g0fa37f0+1,16.0-5-g188214f+1,16.0-5-g81851deb+3,16.0-5-gb866ac2+2,16.0-6-g44ca919+1,16.0-6-gce2b4a7+2,16.0-6-gf89f2d9+2,16.0-7-g1411911+2,16.0-7-g2664ab2+1,16.0-7-ge62bbd2a+2,16.0-8-g2ce35ff+2,16.0-8-g7a94e69+2,w.2018.29
LSSTDataManagementBasePackage
utils.py
Go to the documentation of this file.
1 # See COPYRIGHT file at the top of the source tree.
2 """
3 Statistics of jointcal vs. single-frame procesing and diagnostic plots.
4 
5 Notes
6 -----
7 Some of the algorithms and data structures in this code are temporary
8 kludges and will no longer be necessary once the following are available:
9 
10 - a composite data structure that contains all ccds from a single visit
11 - an n-way matching system that preserves the separations between sources
12 """
13 import collections
14 import os
15 
16 import numpy as np
17 from astropy import units as u
18 
19 import lsst.log
20 import lsst.afw.table
21 import lsst.afw.image
22 from lsst.afw.image import abMagFromFlux
23 from lsst.afw.geom import arcseconds
24 
25 __all__ = ['JointcalStatistics']
26 
27 MatchDict = collections.namedtuple('MatchDict', ['relative', 'absolute'])
28 
29 
31  """
32  Compute statistics on jointcal-processed data, and optionally generate plots.
33 
34  Notes
35  -----
36  Instantiate JointcalStatistics and call compute_rms() to get the relevant
37  statistics for e.g. unittests, and call make_plots() to generate a suite of
38  diagnostic plots.
39  """
40 
41  def __init__(self, match_radius=0.1*arcseconds, flux_limit=100.0,
42  do_photometry=True, do_astrometry=True,
43  verbose=False):
44  """
45  Parameters
46  ----------
47  match_radius : lsst.afw.geom.Angle
48  match sources within this radius for RMS statistics
49  flux_limit : float
50  Signal/Noise (flux/fluxSigma) for sources to be included in the RMS cross-match.
51  100 is a balance between good centroids and enough sources.
52  do_photometry : bool, optional
53  Perform calculations/make plots for photometric metrics.
54  do_astrometry : bool, optional
55  Perform calculations/make plots for astrometric metrics.
56  verbose : bool, optional
57  Print extra things
58  """
59  self.match_radius = match_radius
60  self.flux_limit = flux_limit
61  self.do_photometry = do_photometry
62  self.do_astrometry = do_astrometry
63  self.verbose = verbose
64  self.log = lsst.log.Log.getLogger('JointcalStatistics')
65 
66  def compute_rms(self, data_refs, reference):
67  """
68  Match all data_refs to compute the RMS, for all detections above self.flux_limit.
69 
70  Parameters
71  ----------
72  data_refs : list of lsst.daf.persistence.butlerSubset.ButlerDataRef
73  A list of data refs to do the calculations between.
74  reference : lsst reference catalog
75  reference catalog to do absolute matching against.
76 
77  Returns
78  -------
79  namedtuple:
80  astropy.Quantity
81  Post-jointcal relative RMS of the matched sources.
82  astropy.Quantity
83  Post-jointcal absolute RMS of matched sources.
84  float
85  Post-jointcal photometric repeatability (PA1 from the SRD).
86  """
87 
88  # DECAM doesn't have "filter" in its registry, so we have to get the filter names directly.
89  self.filters = [ref.get('calexp_filter').getName() for ref in data_refs]
90  self.visits_per_dataRef = [ref.dataId['visit'] for ref in data_refs]
91 
92  def compute(catalogs, photoCalibs):
93  """Compute the relative and absolute matches in distance and flux."""
94  visit_catalogs = self._make_visit_catalogs(catalogs, self.visits_per_dataRef)
95  catalogs = [visit_catalogs[x] for x in self.visits_per_dataRef]
96  # use the first catalog as the relative reference catalog
97  # NOTE: The "first" catalog depends on the original ordering of the data_refs.
98  # NOTE: Thus, because I'm doing a many-1 match in _make_match_dict,
99  # the number of matches (and thus the details of the match statistics)
100  # will change if the data_refs are ordered differently.
101  # All the more reason to use a proper n-way matcher here. See: DM-8664
102  refcat = catalogs[0]
103  refcalib = photoCalibs[0] if photoCalibs != [] else None
104  dist_rel, flux_rel, ref_flux_rel, source_rel = self._make_match_dict(refcat,
105  catalogs[1:],
106  photoCalibs[1:],
107  refcalib=refcalib)
108  dist_abs, flux_abs, ref_flux_abs, source_abs = self._make_match_dict(reference,
109  catalogs,
110  photoCalibs)
111  dist = MatchDict(dist_rel, dist_abs)
112  flux = MatchDict(flux_rel, flux_abs)
113  ref_flux = MatchDict(ref_flux_rel, ref_flux_abs)
114  source = MatchDict(source_rel, source_abs)
115  return dist, flux, ref_flux, source
116 
117  old_cats = [ref.get('src') for ref in data_refs]
118  # NOTE: build photoCalibs from existing old Calib objects.
119  # TODO: we can make this a listcomp again once DM-10153 is finished.
120  old_calibs = []
121  if self.do_photometry:
122  for ref in data_refs:
123  calib = ref.get('calexp_calib')
124  fluxMag0 = calib.getFluxMag0()
125  old_calibs.append(lsst.afw.image.PhotoCalib(1.0/fluxMag0[0], fluxMag0[1]/fluxMag0[0]**2))
126 
127  self.old_dist, self.old_flux, self.old_ref_flux, self.old_source = compute(old_cats, old_calibs)
128 
129  # Update coordinates with the new wcs, and get the new photoCalibs.
130  new_cats = [ref.get('src') for ref in data_refs]
131  new_wcss = []
132  if self.do_astrometry:
133  new_wcss = [ref.get('jointcal_wcs') for ref in data_refs]
134  new_calibs = []
135  if self.do_photometry:
136  new_calibs = [ref.get('jointcal_photoCalib') for ref in data_refs]
137  if self.do_astrometry:
138  for wcs, cat in zip(new_wcss, new_cats):
139  # update in-place the object coordinates based on the new wcs
141 
142  self.new_dist, self.new_flux, self.new_ref_flux, self.new_source = compute(new_cats, new_calibs)
143 
144  if self.verbose:
145  print('old, new relative distance matches:',
146  len(self.old_dist.relative), len(self.new_dist.relative))
147  print('old, new absolute distance matches:',
148  len(self.old_dist.absolute), len(self.new_dist.absolute))
149  print('old, new relative flux matches:',
150  len(self.old_flux.relative), len(self.new_flux.relative))
151  print('old, new absolute flux matches:',
152  len(self.old_flux.absolute), len(self.new_flux.absolute))
153 
154  if self.do_photometry:
155  self._photometric_rms()
156  else:
157  self.new_PA1 = None
158 
159  def rms_total(data):
160  """Compute the total rms across all sources."""
161  total = sum(sum(dd**2) for dd in data.values())
162  n = sum(len(dd) for dd in data.values())
163  return np.sqrt(total/n)
164 
165  if self.do_astrometry:
166  self.old_dist_total = MatchDict(*(tuple(map(rms_total, self.old_dist))*u.radian).to(u.arcsecond))
167  self.new_dist_total = MatchDict(*(tuple(map(rms_total, self.new_dist))*u.radian).to(u.arcsecond))
168  else:
169  self.old_dist_total = MatchDict(None, None)
170  self.new_dist_total = MatchDict(None, None)
171 
172  Rms_result = collections.namedtuple("Rms_result", ["dist_relative", "dist_absolute", "pa1"])
173  return Rms_result(self.new_dist_total.relative, self.new_dist_total.absolute, self.new_PA1)
174 
175  def make_plots(self, data_refs, old_wcs_list,
176  name='', interactive=False, per_ccd_plot=False, outdir='.plots'):
177  """
178  Make plots of various quantites to help with debugging.
179  Requires that `compute_rms()` was run first.
180 
181  Parameters
182  ----------
183  data_refs : list of lsst.daf.persistence.butlerSubset.ButlerDataRef
184  A list of data refs to do the calculations between.
185  old_wcs_list : list of lsst.afw.image.wcs.Wcs
186  A list of the old (pre-jointcal) WCSs, one-to-one corresponding to data_refs.
187  name : str
188  Name to include in plot titles and save files.
189  interactive : bool
190  Turn on matplotlib interactive mode and drop into a debugger when
191  plotting is finished. Otherwise, use a non-interactive backend.
192  per_ccd_plot : bool
193  Plot the WCS per CCD (takes longer and generates many plots for a large camera)
194  outdir : str
195  directory to save plots to.
196  """
197  import matplotlib
198 
199  if not interactive:
200  # Use a non-interactive backend for faster plotting.
201  matplotlib.use('pdf')
202 
203  import matplotlib.pyplot as plt
204  import astropy.visualization
205  # make quantities behave nicely when plotted.
206  astropy.visualization.quantity_support()
207  if interactive:
208  plt.ion()
209 
210  self.log.info("N data_refs: %d", len(data_refs))
211 
212  if self.do_photometry:
213  plot_flux_distributions(plt, self.old_mag, self.new_mag,
215  self.faint, self.bright, self.old_PA1, self.new_PA1,
216  name=name, outdir=outdir)
217 
218  def rms_per_source(data):
219  """Each element of data must already be the "delta" of whatever measurement."""
220  return (np.sqrt([np.mean(dd**2) for dd in data.values()])*u.radian).to(u.arcsecond)
221 
222  if self.do_astrometry:
223  old_dist_rms = MatchDict(*(tuple(map(rms_per_source, self.old_dist))))
224  new_dist_rms = MatchDict(*(tuple(map(rms_per_source, self.new_dist))))
225 
226  self.log.info("relative RMS (old, new): {:.2e} {:.2e}".format(self.old_dist_total.relative,
227  self.new_dist_total.relative))
228  self.log.info("absolute RMS (old, new): {:.2e} {:.2e}".format(self.old_dist_total.absolute,
229  self.new_dist_total.absolute))
230  plot_rms_histogram(plt, old_dist_rms.relative, old_dist_rms.absolute,
231  new_dist_rms.relative, new_dist_rms.absolute,
232  self.old_dist_total.relative, self.old_dist_total.absolute,
233  self.new_dist_total.relative, self.new_dist_total.absolute,
234  name=name, outdir=outdir)
235 
236  plot_all_wcs_deltas(plt, data_refs, self.visits_per_dataRef, old_wcs_list,
237  per_ccd_plot=per_ccd_plot,
238  name=name, outdir=outdir)
239 
240  if interactive:
241  plt.show()
242  import pdb
243  pdb.set_trace()
244 
245  def _photometric_rms(self, sn_cut=300, magnitude_range=3):
246  """
247  Compute the photometric RMS and the photometric repeatablity values (PA1).
248 
249  Parameters
250  ----------
251  sn_cut : float
252  The minimum signal/noise for sources to be included in the PA1 calculation.
253  magnitude_range : float
254  The range of magnitudes above sn_cut to include in the PA1 calculation.
255  """
256  def rms(flux, ref_flux):
257  return np.sqrt([np.mean((ref_flux[dd] - flux[dd])**2) for dd in flux])
258 
259  self.old_rms = MatchDict(*map(rms, self.old_flux, self.old_ref_flux))
260  self.new_rms = MatchDict(*map(rms, self.new_flux, self.new_ref_flux))
261 
262  # we want to use the absolute fluxes for all of these calculations.
263  self.old_ref = np.fromiter(self.old_ref_flux.absolute.values(), dtype=float)
264  self.new_ref = np.fromiter(self.new_ref_flux.absolute.values(), dtype=float)
265  self.old_mag = np.fromiter((abMagFromFlux(r) for r in self.old_ref), dtype=float)
266  self.new_mag = np.fromiter((abMagFromFlux(r) for r in self.new_ref), dtype=float)
267 
268  def signal_to_noise(sources, flux_key='slot_PsfFlux_flux', sigma_key='slot_PsfFlux_fluxSigma'):
269  """Compute the mean signal/noise per source from a MatchDict of SourceRecords."""
270  result = np.empty(len(sources))
271  for i, src in enumerate(sources.values()):
272  result[i] = np.mean([x[flux_key]/x[sigma_key] for x in src])
273  return result
274 
275  old_sn = signal_to_noise(self.old_source.absolute)
276  # Find the faint/bright magnitude limits that are the "flat" part of the rms/magnitude relation.
277  self.faint = self.old_mag[old_sn > sn_cut].max()
278  self.bright = self.faint - magnitude_range
279  if self.verbose:
280  print("PA1 Magnitude range: {:.3f}, {:.3f}".format(self.bright, self.faint))
281  old_good = (self.old_mag < self.faint) & (self.old_mag > self.bright)
282  new_good = (self.new_mag < self.faint) & (self.new_mag > self.bright)
283  self.old_weighted_rms = self.old_rms.absolute/self.old_ref
284  self.new_weighted_rms = self.new_rms.absolute/self.new_ref
285  self.old_PA1 = np.median(self.old_weighted_rms[old_good])
286  self.new_PA1 = np.median(self.new_weighted_rms[new_good])
287 
288  def _make_match_dict(self, reference, visit_catalogs, photoCalibs, refcalib=None):
289  """
290  Return several dicts of sourceID:[values] over the catalogs, to be used in RMS calculations.
291 
292  Parameters
293  ----------
294  reference : lsst.afw.table.SourceCatalog
295  Catalog to do the matching against.
296  visit_catalogs : list of lsst.afw.table.SourceCatalog
297  Visit source catalogs (values() produced by _make_visit_catalogs)
298  to cross-match against reference.
299  photoCalibs : list of lsst.afw.image.PhotoCalib
300  Exposure PhotoCalibs, 1-1 coorespondent with visit_catalogs.
301  refcalib : lsst.afw.image.PhotoCalib or None
302  Pass a PhotoCalib here to use it to compute Janskys from the
303  reference catalog ADU slot_flux.
304 
305  Returns
306  -------
307  distances : dict
308  dict of sourceID: array(separation distances for that source)
309  fluxes : dict
310  dict of sourceID: array(fluxes (Jy) for that source)
311  ref_fluxes : dict
312  dict of sourceID: flux (Jy) of the reference object
313  sources : dict
314  dict of sourceID: list(each SourceRecord that was position-matched
315  to this sourceID)
316  """
317  # If we have no photoCalibs, make it the same length as the others for zipping.
318  if photoCalibs == []:
319  photoCalibs = [[]]*len(visit_catalogs)
320 
321  distances = collections.defaultdict(list)
322  fluxes = collections.defaultdict(list)
323  ref_fluxes = {}
324  sources = collections.defaultdict(list)
325  if 'slot_CalibFlux_flux' in reference.schema:
326  ref_flux_key = 'slot_CalibFlux'
327  else:
328  ref_flux_key = '{}_flux'
329 
330  def get_fluxes(photoCalib, match):
331  """Return (flux, ref_flux) or None if either is invalid."""
332  # NOTE: Protect against negative fluxes: ignore this match if we find one.
333  maggiesToJansky = 3631
334  flux = match[1]['slot_CalibFlux_flux']
335  if flux < 0:
336  return None
337  else:
338  flux = maggiesToJansky * photoCalib.instFluxToMaggies(match[1], "slot_CalibFlux").value
339 
340  # NOTE: Have to protect against negative reference fluxes too.
341  if 'slot' in ref_flux_key:
342  ref_flux = match[0][ref_flux_key+'_flux']
343  if ref_flux < 0:
344  return None
345  else:
346  ref_flux = maggiesToJansky * photoCalib.instFluxToMaggies(match[0], ref_flux_key).value
347  else:
348  # a.net fluxes are already in Janskys.
349  ref_flux = match[0][ref_flux_key.format(filt)]
350  if ref_flux < 0:
351  return None
352 
353  Flux = collections.namedtuple('Flux', ('flux', 'ref_flux'))
354  return Flux(flux, ref_flux)
355 
356  for cat, photoCalib, filt in zip(visit_catalogs, photoCalibs, self.filters):
357  good = (cat.get('base_PsfFlux_flux')/cat.get('base_PsfFlux_fluxSigma')) > self.flux_limit
358  # things the classifier called sources are not extended.
359  good &= (cat.get('base_ClassificationExtendedness_value') == 0)
360  matches = lsst.afw.table.matchRaDec(reference, cat[good], self.match_radius)
361  for m in matches:
362  if self.do_photometry:
363  flux = get_fluxes(photoCalib, m)
364  if flux is None:
365  continue
366  else:
367  fluxes[m[0].getId()].append(flux.flux)
368  # we can just use assignment here, since the value is always the same.
369  ref_fluxes[m[0].getId()] = flux.ref_flux
370 
371  if self.do_astrometry:
372  # Just use the computed separation distance directly.
373  distances[m[0].getId()].append(m[2])
374 
375  sources[m[0].getId()].append(m[1])
376  # Convert to numpy array for easier math
377  for source in distances:
378  distances[source] = np.array(distances[source])
379  for source in fluxes:
380  fluxes[source] = np.array(fluxes[source])
381 
382  return distances, fluxes, ref_fluxes, sources
383 
384  def _make_visit_catalogs(self, catalogs, visits):
385  """
386  Merge all catalogs from the each visit.
387  NOTE: creating this structure is somewhat slow, and will be unnecessary
388  once a full-visit composite dataset is available.
389 
390  Parameters
391  ----------
392  catalogs : list of lsst.afw.table.SourceCatalog
393  Catalogs to combine into per-visit catalogs.
394  visits : list of visit id (usually int)
395  list of visit identifiers, one-to-one correspondent with catalogs.
396 
397  Returns
398  -------
399  dict
400  dict of visit: catalog of all sources from all CCDs of that visit.
401  """
402  visit_dict = {v: lsst.afw.table.SourceCatalog(catalogs[0].schema) for v in visits}
403  for v, cat in zip(visits, catalogs):
404  visit_dict[v].extend(cat)
405  # We want catalog contiguity to do object selection later.
406  for v in visit_dict:
407  visit_dict[v] = visit_dict[v].copy(deep=True)
408 
409  return visit_dict
410 
411 
412 def plot_flux_distributions(plt, old_mag, new_mag, old_weighted_rms, new_weighted_rms,
413  faint, bright, old_PA1, new_PA1,
414  name='', outdir='.plots'):
415  """Plot various distributions of fluxes and magnitudes.
416 
417  Parameters
418  ----------
419  plt : matplotlib.pyplot instance
420  pyplot instance to plot with
421  old_mag : np.array
422  old magnitudes
423  new_mag : np.array
424  new magnitudes
425  old_weighted_rms : np.array
426  old rms weighted by the mean (rms(data)/mean(data))
427  new_weighted_rms : np.array
428  old rms weighted by the mean (rms(data)/mean(data))
429  faint : float
430  Faint end of range that PA1 was computed from.
431  bright : float
432  Bright end of range that PA1 was computed from.
433  old_PA1 : float
434  Old value of PA1, to plot as horizontal line.
435  new_PA1 : float
436  New value of PA1, to plot as horizontal line.
437  name : str
438  Name to include in plot titles and save files.
439  outdir : str, optional
440  Directory to write the saved plots to.
441  """
442 
443  import seaborn
444  seaborn.set_style('whitegrid')
445  import scipy.stats
446 
447  old_color = 'blue'
448  new_color = 'red'
449  plt.figure()
450  plt.plot(old_mag, old_weighted_rms, '.', color=old_color, label='old')
451  plt.plot(new_mag, new_weighted_rms, '.', color=new_color, label='new')
452  plt.axvline(faint, ls=':', color=old_color)
453  plt.axvline(bright, ls=':', color=old_color)
454  plt.axhline(old_PA1, ls='--', color=old_color)
455  plt.axhline(new_PA1, ls='--', color=new_color)
456  plt.legend(loc='upper left')
457  plt.title('Where is the systematic flux rms limit?')
458  plt.xlabel('magnitude')
459  plt.ylabel('rms/mean per source')
460  filename = os.path.join(outdir, '{}-photometry-PA1.pdf')
461  plt.savefig(filename.format(name))
462 
463  plt.figure()
464  seaborn.distplot(old_weighted_rms, fit=scipy.stats.lognorm, kde=False, label="old", color=old_color)
465  seaborn.distplot(new_weighted_rms, fit=scipy.stats.lognorm, kde=False, label="new", color=new_color)
466  plt.title('Source RMS pre/post-jointcal')
467  plt.xlabel('rms(flux)/mean(flux)')
468  plt.ylabel('number')
469  plt.legend(loc='upper right')
470  filename = os.path.join(outdir, '{}-photometry-rms.pdf')
471  plt.savefig(filename.format(name))
472 
473 
474 def plot_all_wcs_deltas(plt, data_refs, visits, old_wcs_list, per_ccd_plot=False,
475  name='', outdir='.plots'):
476  """
477  Various plots of the difference between old and new Wcs.
478 
479  Parameters
480  ----------
481  plt : matplotlib.pyplot instance
482  pyplot instance to plot with.
483  data_refs : list of lsst.daf.persistence.butlerSubset.ButlerDataRef
484  A list of data refs to plot.
485  visits : list of visit id (usually int)
486  list of visit identifiers, one-to-one correspondent with catalogs.
487  old_wcs_list : list of lsst.afw.image.wcs.Wcs
488  A list of the old (pre-jointcal) WCSs, one-to-one corresponding to data_refs.
489  per_ccd_plot : bool, optional
490  Make per-ccd plots of the "wcs different" (warning: slow!)
491  name : str
492  Name to include in plot titles and save files.
493  outdir : str, optional
494  Directory to write the saved plots to.
495  """
496 
497  plot_wcs_magnitude(plt, data_refs, visits, old_wcs_list, name, outdir=outdir)
498  plot_all_wcs_quivers(plt, data_refs, visits, old_wcs_list, name, outdir=outdir)
499 
500  if per_ccd_plot:
501  for i, ref in enumerate(data_refs):
502  md = ref.get('calexp_md')
503  dims = lsst.afw.image.bboxFromMetadata(md).getDimensions()
504  plot_wcs(plt, old_wcs_list[i], ref.get('jointcal_wcs'),
505  dims.getX(), dims.getY(),
506  center=(md.getScalar('CRVAL1'), md.getScalar('CRVAL2')), name='dataRef %d'%i,
507  outdir=outdir)
508 
509 
510 def make_xy_wcs_grid(x_dim, y_dim, wcs1, wcs2, num=50):
511  """Return num x/y grid coordinates for wcs1 and wcs2."""
512  x = np.linspace(0, x_dim, num)
513  y = np.linspace(0, y_dim, num)
514  x1, y1 = wcs_convert(x, y, wcs1)
515  x2, y2 = wcs_convert(x, y, wcs2)
516  return x1, y1, x2, y2
517 
518 
519 def wcs_convert(xv, yv, wcs):
520  """Convert two arrays of x/y points into an on-sky grid."""
521  xout = np.zeros((xv.shape[0], yv.shape[0]))
522  yout = np.zeros((xv.shape[0], yv.shape[0]))
523  for i, x in enumerate(xv):
524  for j, y in enumerate(yv):
525  sky = wcs.pixelToSky(x, y)
526  xout[i, j] = sky.getRa()
527  yout[i, j] = sky.getDec()
528  return xout, yout
529 
530 
531 def plot_all_wcs_quivers(plt, data_refs, visits, old_wcs_list, name, outdir='.plots'):
532  """
533  Make quiver plots of the WCS deltas for each CCD in each visit.
534 
535  Parameters
536  ----------
537  plt : matplotlib.pyplot instance
538  pyplot instance to plot with.
539  data_refs : list of lsst.daf.persistence.butlerSubset.ButlerDataRef
540  A list of data refs to plot.
541  visits : list of visit id (usually int)
542  list of visit identifiers, one-to-one correspondent with catalogs.
543  old_wcs_list : list of lsst.afw.image.wcs.Wcs
544  A list of the old (pre-jointcal) WCSs, one-to-one corresponding to data_refs.
545  name : str
546  Name to include in plot titles and save files.
547  outdir : str, optional
548  Directory to write the saved plots to.
549  """
550 
551  for visit in visits:
552  fig = plt.figure()
553  # fig.set_tight_layout(True)
554  ax = fig.add_subplot(111)
555  for old_wcs, ref in zip(old_wcs_list, data_refs):
556  if ref.dataId['visit'] != visit:
557  continue
558  md = ref.get('calexp_md')
559  dims = lsst.afw.image.bboxFromMetadata(md).getDimensions()
560  Q = plot_wcs_quivers(ax, old_wcs, ref.get('jointcal_wcs'),
561  dims.getX(), dims.getY())
562  # TODO: add CCD bounding boxes to plot once DM-5503 is finished.
563  # TODO: add a circle for the full focal plane.
564  length = (0.1*u.arcsecond).to(u.radian).value
565  ax.quiverkey(Q, 0.9, 0.95, length, '0.1 arcsec', coordinates='figure', labelpos='W')
566  plt.xlabel('RA')
567  plt.ylabel('Dec')
568  plt.title('visit: {}'.format(visit))
569  filename = os.path.join(outdir, '{}-{}-quivers.pdf')
570  plt.savefig(filename.format(name, visit))
571 
572 
573 def plot_wcs_quivers(ax, wcs1, wcs2, x_dim, y_dim):
574  """
575  Plot the delta between wcs1 and wcs2 as vector arrows.
576 
577  Parameters
578  ----------
579  ax : matplotlib.axis
580  Matplotlib axis instance to plot to.
581  wcs1 : lsst.afw.image.wcs.Wcs
582  First WCS to compare.
583  wcs2 : lsst.afw.image.wcs.Wcs
584  Second WCS to compare.
585  x_dim : int
586  Size of array in X-coordinate to make the grid over.
587  y_dim : int
588  Size of array in Y-coordinate to make the grid over.
589  """
590 
591  x1, y1, x2, y2 = make_xy_wcs_grid(x_dim, y_dim, wcs1, wcs2)
592  uu = x2 - x1
593  vv = y2 - y1
594  return ax.quiver(x1, y1, uu, vv, units='x', pivot='tail', scale=1e-3, width=1e-5)
595 
596 
597 def plot_wcs_magnitude(plt, data_refs, visits, old_wcs_list, name, outdir='.plots'):
598  """Plot the magnitude of the WCS change between old and new visits as a heat map.
599 
600  Parameters
601  ----------
602  plt : matplotlib.pyplot instance
603  pyplot instance to plot with.
604  data_refs : list of lsst.daf.persistence.butlerSubset.ButlerDataRef
605  A list of data refs to plot.
606  visits : list of visit id (usually int)
607  list of visit identifiers, one-to-one correspondent with catalogs.
608  old_wcs_list : list of lsst.afw.image.wcs.Wcs
609  A list of the old (pre-jointcal) WCSs, one-to-one corresponding to data_refs.
610  name : str
611  Name to include in plot titles and save files.
612  outdir : str, optional
613  Directory to write the saved plots to.
614  """
615  for visit in visits:
616  fig = plt.figure()
617  fig.set_tight_layout(True)
618  ax = fig.add_subplot(111)
619  # Start min/max at the "opposite" ends so they always get the first valid value.
620  xmin = np.inf
621  ymin = np.inf
622  xmax = -np.inf
623  ymax = -np.inf
624  for old_wcs, ref in zip(old_wcs_list, data_refs):
625  if ref.dataId['visit'] != visit:
626  continue
627  md = ref.get('calexp_md')
628  dims = lsst.afw.image.bboxFromMetadata(md).getDimensions()
629  x1, y1, x2, y2 = make_xy_wcs_grid(dims.getX(), dims.getY(),
630  old_wcs, ref.get('jointcal_wcs'))
631  uu = x2 - x1
632  vv = y2 - y1
633  extent = (x1[0, 0], x1[-1, -1], y1[0, 0], y1[-1, -1])
634  xmin = min(x1.min(), xmin)
635  ymin = min(y1.min(), ymin)
636  xmax = max(x1.max(), xmax)
637  ymax = max(y1.max(), ymax)
638  magnitude = (np.linalg.norm((uu, vv), axis=0)*u.radian).to(u.arcsecond).value
639  img = ax.imshow(magnitude, vmin=0, vmax=0.3,
640  aspect='auto', extent=extent, cmap=plt.get_cmap('magma'))
641  # TODO: add CCD bounding boxes to the plot once DM-5503 is finished.
642  # TODO: add a circle for the full focal plane.
643 
644  # We're reusing only one of the returned images here for colorbar scaling,
645  # but it doesn't matter because we set vmin/vmax so they are all scaled the same.
646  cbar = plt.colorbar(img)
647  cbar.ax.set_ylabel('distortion (arcseconds)')
648  plt.xlim(xmin, xmax)
649  plt.ylim(ymin, ymax)
650  plt.xlabel('RA')
651  plt.ylabel('Dec')
652  plt.title('visit: {}'.format(visit))
653  filename = os.path.join(outdir, '{}-{}-heatmap.pdf')
654  plt.savefig(filename.format(name, visit))
655 
656 
657 def plot_wcs(plt, wcs1, wcs2, x_dim, y_dim, center=(0, 0), name="", outdir='.plots'):
658  """Plot the "distortion map": wcs1-wcs2 delta of points in the CCD grid.
659 
660  Parameters
661  ----------
662  plt : matplotlib.pyplot instance
663  pyplot instance to plot with.
664  wcs1 : lsst.afw.image.wcs.Wcs
665  First WCS to compare.
666  wcs2 : lsst.afw.image.wcs.Wcs
667  Second WCS to compare.
668  x_dim : int
669  Size of array in X-coordinate to make the grid over.
670  y_dim : int
671  Size of array in Y-coordinate to make the grid over.
672  center : tuple, optional
673  Center of the data, in on-chip coordinates.
674  name : str
675  Name to include in plot titles and save files.
676  outdir : str, optional
677  Directory to write the saved plots to.
678  """
679 
680  plt.figure()
681 
682  x1, y1, x2, y2 = make_xy_wcs_grid(x_dim, y_dim, wcs1, wcs2, num=50)
683  plt.plot((x1 - x2) + center[0], (y1 - y2) + center[1], '-')
684  plt.xlabel('delta RA (arcsec)')
685  plt.ylabel('delta Dec (arcsec)')
686  plt.title(name)
687  filename = os.path.join(outdir, '{}-wcs.pdf')
688  plt.savefig(filename.format(name))
689 
690 
691 def plot_rms_histogram(plt, old_rms_relative, old_rms_absolute,
692  new_rms_relative, new_rms_absolute,
693  old_rel_total, old_abs_total, new_rel_total, new_abs_total,
694  name="", outdir='.plots'):
695  """Plot histograms of the source separations and their RMS values.
696 
697  Parameters
698  ----------
699  plt : matplotlib.pyplot instance
700  pyplot instance to plot with.
701  old_rms_relative : np.array
702  old relative rms/star
703  old_rms_absolute : np.array
704  old absolute rms/star
705  new_rms_relative : np.array
706  new relative rms/star
707  new_rms_absolute : np.array
708  new absolute rms/star
709  old_rel_total : float
710  old relative rms over all stars
711  old_abs_total : float
712  old absolute rms over all stars
713  new_rel_total : float
714  new relative rms over all stars
715  new_abs_total : float
716  new absolute rms over all stars
717  name : str
718  Name to include in plot titles and save files.
719  outdir : str, optional
720  Directory to write the saved plots to.
721  """
722  plt.figure()
723 
724  color_rel = 'black'
725  ls_old = 'dotted'
726  color_abs = 'green'
727  ls_new = 'dashed'
728  plotOptions = {'lw': 2, 'range': (0, 0.1)*u.arcsecond, 'normed': True,
729  'bins': 30, 'histtype': 'step'}
730 
731  plt.title('relative vs. absolute: %d vs. %d'%(len(old_rms_relative), len(old_rms_absolute)))
732 
733  plt.hist(old_rms_absolute, color=color_abs, ls=ls_old, label='old abs', **plotOptions)
734  plt.hist(new_rms_absolute, color=color_abs, ls=ls_new, label='new abs', **plotOptions)
735 
736  plt.hist(old_rms_relative, color=color_rel, ls=ls_old, label='old rel', **plotOptions)
737  plt.hist(new_rms_relative, color=color_rel, ls=ls_new, label='new rel', **plotOptions)
738 
739  plt.axvline(x=old_abs_total.value, linewidth=1.5, color=color_abs, ls=ls_old)
740  plt.axvline(x=new_abs_total.value, linewidth=1.5, color=color_abs, ls=ls_new)
741  plt.axvline(x=old_rel_total.value, linewidth=1.5, color=color_rel, ls=ls_old)
742  plt.axvline(x=new_rel_total.value, linewidth=1.5, color=color_rel, ls=ls_new)
743 
744  plt.xlim(plotOptions['range'])
745  plt.xlabel('arcseconds')
746  plt.legend(loc='best')
747  filename = os.path.join(outdir, '{}-histogram.pdf')
748  plt.savefig(filename.format(name))
def _photometric_rms(self, sn_cut=300, magnitude_range=3)
Definition: utils.py:245
def make_plots(self, data_refs, old_wcs_list, name='', interactive=False, per_ccd_plot=False, outdir='.plots')
Definition: utils.py:176
def __init__(self, match_radius=0.1 *arcseconds, flux_limit=100.0, do_photometry=True, do_astrometry=True, verbose=False)
Definition: utils.py:43
The photometric calibration of an exposure.
Definition: PhotoCalib.h:81
def wcs_convert(xv, yv, wcs)
Definition: utils.py:519
int min
Definition: Coord.cc:82
def plot_rms_histogram(plt, old_rms_relative, old_rms_absolute, new_rms_relative, new_rms_absolute, old_rel_total, old_abs_total, new_rel_total, new_abs_total, name="", outdir='.plots')
Definition: utils.py:694
def _make_match_dict(self, reference, visit_catalogs, photoCalibs, refcalib=None)
Definition: utils.py:288
std::shared_ptr< FrameSet > append(FrameSet const &first, FrameSet const &second)
Construct a FrameSet that performs two transformations in series.
Definition: functional.cc:33
STL class.
double abMagFromFlux(double flux)
Compute AB magnitude from flux in Janskys.
Definition: Calib.h:58
template SourceMatchVector matchRaDec(SourceCatalog const &, lsst::geom::Angle, MatchControl const &)
void updateSourceCoords(geom::SkyWcs const &wcs, SourceCollection &sourceList)
Update sky coordinates in a collection of source objects.
Definition: wcsUtils.cc:96
Definition: Log.h:716
def plot_all_wcs_deltas(plt, data_refs, visits, old_wcs_list, per_ccd_plot=False, name='', outdir='.plots')
Definition: utils.py:475
def plot_wcs(plt, wcs1, wcs2, x_dim, y_dim, center=(0, 0), name="", outdir='.plots')
Definition: utils.py:657
def plot_wcs_quivers(ax, wcs1, wcs2, x_dim, y_dim)
Definition: utils.py:573
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:134
int max
Definition: BoundedField.cc:98
def make_xy_wcs_grid(x_dim, y_dim, wcs1, wcs2, num=50)
Definition: utils.py:510
def plot_flux_distributions(plt, old_mag, new_mag, old_weighted_rms, new_weighted_rms, faint, bright, old_PA1, new_PA1, name='', outdir='.plots')
Definition: utils.py:414
def compute_rms(self, data_refs, reference)
Definition: utils.py:66
static Log getLogger(Log const &logger)
Definition: Log.h:772
def plot_wcs_magnitude(plt, data_refs, visits, old_wcs_list, name, outdir='.plots')
Definition: utils.py:597
def plot_all_wcs_quivers(plt, data_refs, visits, old_wcs_list, name, outdir='.plots')
Definition: utils.py:531
lsst::geom::Box2I bboxFromMetadata(daf::base::PropertySet &metadata)
Determine the image bounding box from its metadata (FITS header)
Definition: Image.cc:703
def _make_visit_catalogs(self, catalogs, visits)
Definition: utils.py:384