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