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