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