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