lsst.ip.diffim  17.0.1-11-ge9de802+8
utils.py
Go to the documentation of this file.
1 # This file is part of ip_diffim.
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 """Support utilities for Measuring sources"""
23 
24 # Export DipoleTestImage to expose fake image generating funcs
25 __all__ = ["DipoleTestImage"]
26 
27 import numpy as np
28 
29 import lsst.afw.detection as afwDet
30 import lsst.afw.display as afwDisplay
31 import lsst.afw.geom as afwGeom
32 import lsst.afw.image as afwImage
33 import lsst.afw.math as afwMath
34 import lsst.afw.table as afwTable
35 from lsst.log import Log
36 import lsst.meas.algorithms as measAlg
37 import lsst.meas.base as measBase
38 from .dipoleFitTask import DipoleFitAlgorithm
39 from . import diffimLib
40 from . import diffimTools
41 
42 afwDisplay.setDefaultMaskTransparency(75)
43 keptPlots = False # Have we arranged to keep spatial plots open?
44 
45 
46 def showSourceSet(sSet, xy0=(0, 0), frame=0, ctype=afwDisplay.GREEN, symb="+", size=2):
47  """Draw the (XAstrom, YAstrom) positions of a set of Sources.
48 
49  Image has the given XY0.
50  """
51  disp = afwDisplay.afwDisplay(frame=frame)
52  with disp.Buffering():
53  for s in sSet:
54  xc, yc = s.getXAstrom() - xy0[0], s.getYAstrom() - xy0[1]
55 
56  if symb == "id":
57  disp.dot(str(s.getId()), xc, yc, ctype=ctype, size=size)
58  else:
59  disp.dot(symb, xc, yc, ctype=ctype, size=size)
60 
61 
62 # Kernel display utilities
63 #
64 
65 
66 def showKernelSpatialCells(maskedIm, kernelCellSet, showChi2=False, symb="o",
67  ctype=None, ctypeUnused=None, ctypeBad=None, size=3,
68  frame=None, title="Spatial Cells"):
69  """Show the SpatialCells.
70 
71  If symb is something that display.dot understands (e.g. "o"), the top
72  nMaxPerCell candidates will be indicated with that symbol, using ctype
73  and size.
74  """
75  disp = afwDisplay.Display(frame=frame)
76  disp.mtv(maskedIm, title=title)
77  with disp.Buffering():
78  origin = [-maskedIm.getX0(), -maskedIm.getY0()]
79  for cell in kernelCellSet.getCellList():
80  afwDisplay.utils.drawBBox(cell.getBBox(), origin=origin, display=disp)
81 
82  goodies = ctypeBad is None
83  for cand in cell.begin(goodies):
84  xc, yc = cand.getXCenter() + origin[0], cand.getYCenter() + origin[1]
85  if cand.getStatus() == afwMath.SpatialCellCandidate.BAD:
86  color = ctypeBad
87  elif cand.getStatus() == afwMath.SpatialCellCandidate.GOOD:
88  color = ctype
89  elif cand.getStatus() == afwMath.SpatialCellCandidate.UNKNOWN:
90  color = ctypeUnused
91  else:
92  continue
93 
94  if color:
95  disp.dot(symb, xc, yc, ctype=color, size=size)
96 
97  if showChi2:
98  rchi2 = cand.getChi2()
99  if rchi2 > 1e100:
100  rchi2 = np.nan
101  disp.dot("%d %.1f" % (cand.getId(), rchi2),
102  xc - size, yc - size - 4, ctype=color, size=size)
103 
104 
105 def showDiaSources(sources, exposure, isFlagged, isDipole, frame=None):
106  """Display Dia Sources.
107  """
108  #
109  # Show us the ccandidates
110  #
111  # Too many mask planes in diffims
112  disp = afwDisplay.Display(frame=frame)
113  for plane in ("BAD", "CR", "EDGE", "INTERPOlATED", "INTRP", "SAT", "SATURATED"):
114  disp.setMaskPlaneColor(plane, color="ignore")
115 
116  mos = afwDisplay.utils.Mosaic()
117  for i in range(len(sources)):
118  source = sources[i]
119  badFlag = isFlagged[i]
120  dipoleFlag = isDipole[i]
121  bbox = source.getFootprint().getBBox()
122  stamp = exposure.Factory(exposure, bbox, True)
123  im = afwDisplay.utils.Mosaic(gutter=1, background=0, mode="x")
124  im.append(stamp.getMaskedImage())
125  lab = "%.1f,%.1f:" % (source.getX(), source.getY())
126  if badFlag:
127  ctype = afwDisplay.RED
128  lab += "BAD"
129  if dipoleFlag:
130  ctype = afwDisplay.YELLOW
131  lab += "DIPOLE"
132  if not badFlag and not dipoleFlag:
133  ctype = afwDisplay.GREEN
134  lab += "OK"
135  mos.append(im.makeMosaic(), lab, ctype)
136  title = "Dia Sources"
137  mosaicImage = mos.makeMosaic(frame=frame, title=title)
138  return mosaicImage
139 
140 
141 def showKernelCandidates(kernelCellSet, kernel, background, frame=None, showBadCandidates=True,
142  resids=False, kernels=False):
143  """Display the Kernel candidates.
144 
145  If kernel is provided include spatial model and residuals;
146  If chi is True, generate a plot of residuals/sqrt(variance), i.e. chi.
147  """
148  #
149  # Show us the ccandidates
150  #
151  if kernels:
152  mos = afwDisplay.utils.Mosaic(gutter=5, background=0)
153  else:
154  mos = afwDisplay.utils.Mosaic(gutter=5, background=-1)
155  #
156  candidateCenters = []
157  candidateCentersBad = []
158  candidateIndex = 0
159  for cell in kernelCellSet.getCellList():
160  for cand in cell.begin(False): # include bad candidates
161  # Original difference image; if does not exist, skip candidate
162  try:
163  resid = cand.getDifferenceImage(diffimLib.KernelCandidateF.ORIG)
164  except Exception:
165  continue
166 
167  rchi2 = cand.getChi2()
168  if rchi2 > 1e100:
169  rchi2 = np.nan
170 
171  if not showBadCandidates and cand.isBad():
172  continue
173 
174  im_resid = afwDisplay.utils.Mosaic(gutter=1, background=-0.5, mode="x")
175 
176  try:
177  im = cand.getScienceMaskedImage()
178  im = im.Factory(im, True)
179  im.setXY0(cand.getScienceMaskedImage().getXY0())
180  except Exception:
181  continue
182  if (not resids and not kernels):
183  im_resid.append(im.Factory(im, True))
184  try:
185  im = cand.getTemplateMaskedImage()
186  im = im.Factory(im, True)
187  im.setXY0(cand.getTemplateMaskedImage().getXY0())
188  except Exception:
189  continue
190  if (not resids and not kernels):
191  im_resid.append(im.Factory(im, True))
192 
193  # Difference image with original basis
194  if resids:
195  var = resid.getVariance()
196  var = var.Factory(var, True)
197  np.sqrt(var.getArray(), var.getArray()) # inplace sqrt
198  resid = resid.getImage()
199  resid /= var
200  bbox = kernel.shrinkBBox(resid.getBBox())
201  resid = resid.Factory(resid, bbox, deep=True)
202  elif kernels:
203  kim = cand.getKernelImage(diffimLib.KernelCandidateF.ORIG).convertF()
204  resid = kim.Factory(kim, True)
205  im_resid.append(resid)
206 
207  # residuals using spatial model
208  ski = afwImage.ImageD(kernel.getDimensions())
209  kernel.computeImage(ski, False, int(cand.getXCenter()), int(cand.getYCenter()))
210  sk = afwMath.FixedKernel(ski)
211  sbg = 0.0
212  if background:
213  sbg = background(int(cand.getXCenter()), int(cand.getYCenter()))
214  sresid = cand.getDifferenceImage(sk, sbg)
215  resid = sresid
216  if resids:
217  resid = sresid.getImage()
218  resid /= var
219  bbox = kernel.shrinkBBox(resid.getBBox())
220  resid = resid.Factory(resid, bbox, deep=True)
221  elif kernels:
222  kim = ski.convertF()
223  resid = kim.Factory(kim, True)
224  im_resid.append(resid)
225 
226  im = im_resid.makeMosaic()
227 
228  lab = "%d chi^2 %.1f" % (cand.getId(), rchi2)
229  ctype = afwDisplay.RED if cand.isBad() else afwDisplay.GREEN
230 
231  mos.append(im, lab, ctype)
232 
233  if False and np.isnan(rchi2):
234  disp = afwDisplay.Display(frame=1)
235  disp.mtv(cand.getScienceMaskedImage.getImage(), title="candidate")
236  print("rating", cand.getCandidateRating())
237 
238  im = cand.getScienceMaskedImage()
239  center = (candidateIndex, cand.getXCenter() - im.getX0(), cand.getYCenter() - im.getY0())
240  candidateIndex += 1
241  if cand.isBad():
242  candidateCentersBad.append(center)
243  else:
244  candidateCenters.append(center)
245 
246  if resids:
247  title = "chi Diffim"
248  elif kernels:
249  title = "Kernels"
250  else:
251  title = "Candidates & residuals"
252  mosaicImage = mos.makeMosaic(frame=frame, title=title)
253 
254  return mosaicImage
255 
256 
257 def showKernelBasis(kernel, frame=None):
258  """Display a Kernel's basis images.
259  """
260  mos = afwDisplay.utils.Mosaic()
261 
262  for k in kernel.getKernelList():
263  im = afwImage.ImageD(k.getDimensions())
264  k.computeImage(im, False)
265  mos.append(im)
266  mos.makeMosaic(frame=frame, title="Kernel Basis Images")
267 
268  return mos
269 
270 
271 
272 
273 def plotKernelSpatialModel(kernel, kernelCellSet, showBadCandidates=True,
274  numSample=128, keepPlots=True, maxCoeff=10):
275  """Plot the Kernel spatial model.
276  """
277  try:
278  import matplotlib.pyplot as plt
279  import matplotlib.colors
280  except ImportError as e:
281  print("Unable to import numpy and matplotlib: %s" % e)
282  return
283 
284  x0 = kernelCellSet.getBBox().getBeginX()
285  y0 = kernelCellSet.getBBox().getBeginY()
286 
287  candPos = list()
288  candFits = list()
289  badPos = list()
290  badFits = list()
291  candAmps = list()
292  badAmps = list()
293  for cell in kernelCellSet.getCellList():
294  for cand in cell.begin(False):
295  if not showBadCandidates and cand.isBad():
296  continue
297  candCenter = afwGeom.PointD(cand.getXCenter(), cand.getYCenter())
298  try:
299  im = cand.getTemplateMaskedImage()
300  except Exception:
301  continue
302 
303  targetFits = badFits if cand.isBad() else candFits
304  targetPos = badPos if cand.isBad() else candPos
305  targetAmps = badAmps if cand.isBad() else candAmps
306 
307  # compare original and spatial kernel coefficients
308  kp0 = np.array(cand.getKernel(diffimLib.KernelCandidateF.ORIG).getKernelParameters())
309  amp = cand.getCandidateRating()
310 
311  targetFits = badFits if cand.isBad() else candFits
312  targetPos = badPos if cand.isBad() else candPos
313  targetAmps = badAmps if cand.isBad() else candAmps
314 
315  targetFits.append(kp0)
316  targetPos.append(candCenter)
317  targetAmps.append(amp)
318 
319  xGood = np.array([pos.getX() for pos in candPos]) - x0
320  yGood = np.array([pos.getY() for pos in candPos]) - y0
321  zGood = np.array(candFits)
322 
323  xBad = np.array([pos.getX() for pos in badPos]) - x0
324  yBad = np.array([pos.getY() for pos in badPos]) - y0
325  zBad = np.array(badFits)
326  numBad = len(badPos)
327 
328  xRange = np.linspace(0, kernelCellSet.getBBox().getWidth(), num=numSample)
329  yRange = np.linspace(0, kernelCellSet.getBBox().getHeight(), num=numSample)
330 
331  if maxCoeff:
332  maxCoeff = min(maxCoeff, kernel.getNKernelParameters())
333  else:
334  maxCoeff = kernel.getNKernelParameters()
335 
336  for k in range(maxCoeff):
337  func = kernel.getSpatialFunction(k)
338  dfGood = zGood[:, k] - np.array([func(pos.getX(), pos.getY()) for pos in candPos])
339  yMin = dfGood.min()
340  yMax = dfGood.max()
341  if numBad > 0:
342  dfBad = zBad[:, k] - np.array([func(pos.getX(), pos.getY()) for pos in badPos])
343  # Can really screw up the range...
344  yMin = min([yMin, dfBad.min()])
345  yMax = max([yMax, dfBad.max()])
346  yMin -= 0.05*(yMax - yMin)
347  yMax += 0.05*(yMax - yMin)
348 
349  fRange = np.ndarray((len(xRange), len(yRange)))
350  for j, yVal in enumerate(yRange):
351  for i, xVal in enumerate(xRange):
352  fRange[j][i] = func(xVal, yVal)
353 
354  fig = plt.figure(k)
355 
356  fig.clf()
357  try:
358  fig.canvas._tkcanvas._root().lift() # == Tk's raise, but raise is a python reserved word
359  except Exception: # protect against API changes
360  pass
361 
362  fig.suptitle('Kernel component %d' % k)
363 
364  # LL
365  ax = fig.add_axes((0.1, 0.05, 0.35, 0.35))
366  vmin = fRange.min() # - 0.05*np.fabs(fRange.min())
367  vmax = fRange.max() # + 0.05*np.fabs(fRange.max())
368  norm = matplotlib.colors.Normalize(vmin=vmin, vmax=vmax)
369  im = ax.imshow(fRange, aspect='auto', norm=norm,
370  extent=[0, kernelCellSet.getBBox().getWidth() - 1,
371  0, kernelCellSet.getBBox().getHeight() - 1])
372  ax.set_title('Spatial polynomial')
373  plt.colorbar(im, orientation='horizontal', ticks=[vmin, vmax])
374 
375  # UL
376  ax = fig.add_axes((0.1, 0.55, 0.35, 0.35))
377  ax.plot(-2.5*np.log10(candAmps), zGood[:, k], 'b+')
378  if numBad > 0:
379  ax.plot(-2.5*np.log10(badAmps), zBad[:, k], 'r+')
380  ax.set_title("Basis Coefficients")
381  ax.set_xlabel("Instr mag")
382  ax.set_ylabel("Coeff")
383 
384  # LR
385  ax = fig.add_axes((0.55, 0.05, 0.35, 0.35))
386  ax.set_autoscale_on(False)
387  ax.set_xbound(lower=0, upper=kernelCellSet.getBBox().getHeight())
388  ax.set_ybound(lower=yMin, upper=yMax)
389  ax.plot(yGood, dfGood, 'b+')
390  if numBad > 0:
391  ax.plot(yBad, dfBad, 'r+')
392  ax.axhline(0.0)
393  ax.set_title('dCoeff (indiv-spatial) vs. y')
394 
395  # UR
396  ax = fig.add_axes((0.55, 0.55, 0.35, 0.35))
397  ax.set_autoscale_on(False)
398  ax.set_xbound(lower=0, upper=kernelCellSet.getBBox().getWidth())
399  ax.set_ybound(lower=yMin, upper=yMax)
400  ax.plot(xGood, dfGood, 'b+')
401  if numBad > 0:
402  ax.plot(xBad, dfBad, 'r+')
403  ax.axhline(0.0)
404  ax.set_title('dCoeff (indiv-spatial) vs. x')
405 
406  fig.show()
407 
408  global keptPlots
409  if keepPlots and not keptPlots:
410  # Keep plots open when done
411  def show():
412  print("%s: Please close plots when done." % __name__)
413  try:
414  plt.show()
415  except Exception:
416  pass
417  print("Plots closed, exiting...")
418  import atexit
419  atexit.register(show)
420  keptPlots = True
421 
422 
423 def showKernelMosaic(bbox, kernel, nx=7, ny=None, frame=None, title=None,
424  showCenter=True, showEllipticity=True):
425  """Show a mosaic of Kernel images.
426  """
427  mos = afwDisplay.utils.Mosaic()
428 
429  x0 = bbox.getBeginX()
430  y0 = bbox.getBeginY()
431  width = bbox.getWidth()
432  height = bbox.getHeight()
433 
434  if not ny:
435  ny = int(nx*float(height)/width + 0.5)
436  if not ny:
437  ny = 1
438 
439  schema = afwTable.SourceTable.makeMinimalSchema()
440  centroidName = "base_SdssCentroid"
441  shapeName = "base_SdssShape"
442  control = measBase.SdssCentroidControl()
443  schema.getAliasMap().set("slot_Centroid", centroidName)
444  schema.getAliasMap().set("slot_Centroid_flag", centroidName + "_flag")
445  centroider = measBase.SdssCentroidAlgorithm(control, centroidName, schema)
446  sdssShape = measBase.SdssShapeControl()
447  shaper = measBase.SdssShapeAlgorithm(sdssShape, shapeName, schema)
448  table = afwTable.SourceTable.make(schema)
449  table.defineCentroid(centroidName)
450  table.defineShape(shapeName)
451 
452  centers = []
453  shapes = []
454  for iy in range(ny):
455  for ix in range(nx):
456  x = int(ix*(width - 1)/(nx - 1)) + x0
457  y = int(iy*(height - 1)/(ny - 1)) + y0
458 
459  im = afwImage.ImageD(kernel.getDimensions())
460  ksum = kernel.computeImage(im, False, x, y)
461  lab = "Kernel(%d,%d)=%.2f" % (x, y, ksum) if False else ""
462  mos.append(im, lab)
463 
464  # SdssCentroidAlgorithm.measure requires an exposure of floats
465  exp = afwImage.makeExposure(afwImage.makeMaskedImage(im.convertF()))
466 
467  w, h = im.getWidth(), im.getHeight()
468  centerX = im.getX0() + w//2
469  centerY = im.getY0() + h//2
470  src = table.makeRecord()
471  spans = afwGeom.SpanSet(exp.getBBox())
472  foot = afwDet.Footprint(spans)
473  foot.addPeak(centerX, centerY, 1)
474  src.setFootprint(foot)
475 
476  try: # The centroider requires a psf, so this will fail if none is attached to exp
477  centroider.measure(src, exp)
478  centers.append((src.getX(), src.getY()))
479 
480  shaper.measure(src, exp)
481  shapes.append((src.getIxx(), src.getIxy(), src.getIyy()))
482  except Exception:
483  pass
484 
485  mos.makeMosaic(frame=frame, title=title if title else "Model Kernel", mode=nx)
486 
487  if centers and frame is not None:
488  disp = afwDisplay.Display(frame=frame)
489  i = 0
490  with disp.Buffering():
491  for cen, shape in zip(centers, shapes):
492  bbox = mos.getBBox(i)
493  i += 1
494  xc, yc = cen[0] + bbox.getMinX(), cen[1] + bbox.getMinY()
495  if showCenter:
496  disp.dot("+", xc, yc, ctype=afwDisplay.BLUE)
497 
498  if showEllipticity:
499  ixx, ixy, iyy = shape
500  disp.dot("@:%g,%g,%g" % (ixx, ixy, iyy), xc, yc, ctype=afwDisplay.RED)
501 
502  return mos
503 
504 
505 def plotPixelResiduals(exposure, warpedTemplateExposure, diffExposure, kernelCellSet,
506  kernel, background, testSources, config,
507  origVariance=False, nptsFull=1e6, keepPlots=True, titleFs=14):
508  """Plot diffim residuals for LOCAL and SPATIAL models.
509  """
510  candidateResids = []
511  spatialResids = []
512  nonfitResids = []
513 
514  for cell in kernelCellSet.getCellList():
515  for cand in cell.begin(True): # only look at good ones
516  # Be sure
517  if not (cand.getStatus() == afwMath.SpatialCellCandidate.GOOD):
518  continue
519 
520  diffim = cand.getDifferenceImage(diffimLib.KernelCandidateF.ORIG)
521  orig = cand.getScienceMaskedImage()
522 
523  ski = afwImage.ImageD(kernel.getDimensions())
524  kernel.computeImage(ski, False, int(cand.getXCenter()), int(cand.getYCenter()))
525  sk = afwMath.FixedKernel(ski)
526  sbg = background(int(cand.getXCenter()), int(cand.getYCenter()))
527  sdiffim = cand.getDifferenceImage(sk, sbg)
528 
529  # trim edgs due to convolution
530  bbox = kernel.shrinkBBox(diffim.getBBox())
531  tdiffim = diffim.Factory(diffim, bbox)
532  torig = orig.Factory(orig, bbox)
533  tsdiffim = sdiffim.Factory(sdiffim, bbox)
534 
535  if origVariance:
536  candidateResids.append(np.ravel(tdiffim.getImage().getArray() /
537  np.sqrt(torig.getVariance().getArray())))
538  spatialResids.append(np.ravel(tsdiffim.getImage().getArray() /
539  np.sqrt(torig.getVariance().getArray())))
540  else:
541  candidateResids.append(np.ravel(tdiffim.getImage().getArray() /
542  np.sqrt(tdiffim.getVariance().getArray())))
543  spatialResids.append(np.ravel(tsdiffim.getImage().getArray() /
544  np.sqrt(tsdiffim.getVariance().getArray())))
545 
546  fullIm = diffExposure.getMaskedImage().getImage().getArray()
547  fullMask = diffExposure.getMaskedImage().getMask().getArray()
548  if origVariance:
549  fullVar = exposure.getMaskedImage().getVariance().getArray()
550  else:
551  fullVar = diffExposure.getMaskedImage().getVariance().getArray()
552 
553  bitmaskBad = 0
554  bitmaskBad |= afwImage.Mask.getPlaneBitMask('NO_DATA')
555  bitmaskBad |= afwImage.Mask.getPlaneBitMask('SAT')
556  idx = np.where((fullMask & bitmaskBad) == 0)
557  stride = int(len(idx[0])//nptsFull)
558  sidx = idx[0][::stride], idx[1][::stride]
559  allResids = fullIm[sidx]/np.sqrt(fullVar[sidx])
560 
561  testFootprints = diffimTools.sourceToFootprintList(testSources, warpedTemplateExposure,
562  exposure, config, Log.getDefaultLogger())
563  for fp in testFootprints:
564  subexp = diffExposure.Factory(diffExposure, fp["footprint"].getBBox())
565  subim = subexp.getMaskedImage().getImage()
566  if origVariance:
567  subvar = afwImage.ExposureF(exposure, fp["footprint"].getBBox()).getMaskedImage().getVariance()
568  else:
569  subvar = subexp.getMaskedImage().getVariance()
570  nonfitResids.append(np.ravel(subim.getArray()/np.sqrt(subvar.getArray())))
571 
572  candidateResids = np.ravel(np.array(candidateResids))
573  spatialResids = np.ravel(np.array(spatialResids))
574  nonfitResids = np.ravel(np.array(nonfitResids))
575 
576  try:
577  import pylab
578  from matplotlib.font_manager import FontProperties
579  except ImportError as e:
580  print("Unable to import pylab: %s" % e)
581  return
582 
583  fig = pylab.figure()
584  fig.clf()
585  try:
586  fig.canvas._tkcanvas._root().lift() # == Tk's raise, but raise is a python reserved word
587  except Exception: # protect against API changes
588  pass
589  if origVariance:
590  fig.suptitle("Diffim residuals: Normalized by sqrt(input variance)", fontsize=titleFs)
591  else:
592  fig.suptitle("Diffim residuals: Normalized by sqrt(diffim variance)", fontsize=titleFs)
593 
594  sp1 = pylab.subplot(221)
595  sp2 = pylab.subplot(222, sharex=sp1, sharey=sp1)
596  sp3 = pylab.subplot(223, sharex=sp1, sharey=sp1)
597  sp4 = pylab.subplot(224, sharex=sp1, sharey=sp1)
598  xs = np.arange(-5, 5.05, 0.1)
599  ys = 1./np.sqrt(2*np.pi)*np.exp(-0.5*xs**2)
600 
601  sp1.hist(candidateResids, bins=xs, normed=True, alpha=0.5, label="N(%.2f, %.2f)"
602  % (np.mean(candidateResids), np.var(candidateResids)))
603  sp1.plot(xs, ys, "r-", lw=2, label="N(0,1)")
604  sp1.set_title("Candidates: basis fit", fontsize=titleFs - 2)
605  sp1.legend(loc=1, fancybox=True, shadow=True, prop=FontProperties(size=titleFs - 6))
606 
607  sp2.hist(spatialResids, bins=xs, normed=True, alpha=0.5, label="N(%.2f, %.2f)"
608  % (np.mean(spatialResids), np.var(spatialResids)))
609  sp2.plot(xs, ys, "r-", lw=2, label="N(0,1)")
610  sp2.set_title("Candidates: spatial fit", fontsize=titleFs - 2)
611  sp2.legend(loc=1, fancybox=True, shadow=True, prop=FontProperties(size=titleFs - 6))
612 
613  sp3.hist(nonfitResids, bins=xs, normed=True, alpha=0.5, label="N(%.2f, %.2f)"
614  % (np.mean(nonfitResids), np.var(nonfitResids)))
615  sp3.plot(xs, ys, "r-", lw=2, label="N(0,1)")
616  sp3.set_title("Control sample: spatial fit", fontsize=titleFs - 2)
617  sp3.legend(loc=1, fancybox=True, shadow=True, prop=FontProperties(size=titleFs - 6))
618 
619  sp4.hist(allResids, bins=xs, normed=True, alpha=0.5, label="N(%.2f, %.2f)"
620  % (np.mean(allResids), np.var(allResids)))
621  sp4.plot(xs, ys, "r-", lw=2, label="N(0,1)")
622  sp4.set_title("Full image (subsampled)", fontsize=titleFs - 2)
623  sp4.legend(loc=1, fancybox=True, shadow=True, prop=FontProperties(size=titleFs - 6))
624 
625  pylab.setp(sp1.get_xticklabels() + sp1.get_yticklabels(), fontsize=titleFs - 4)
626  pylab.setp(sp2.get_xticklabels() + sp2.get_yticklabels(), fontsize=titleFs - 4)
627  pylab.setp(sp3.get_xticklabels() + sp3.get_yticklabels(), fontsize=titleFs - 4)
628  pylab.setp(sp4.get_xticklabels() + sp4.get_yticklabels(), fontsize=titleFs - 4)
629 
630  sp1.set_xlim(-5, 5)
631  sp1.set_ylim(0, 0.5)
632  fig.show()
633 
634  global keptPlots
635  if keepPlots and not keptPlots:
636  # Keep plots open when done
637  def show():
638  print("%s: Please close plots when done." % __name__)
639  try:
640  pylab.show()
641  except Exception:
642  pass
643  print("Plots closed, exiting...")
644  import atexit
645  atexit.register(show)
646  keptPlots = True
647 
648 
649 def calcCentroid(arr):
650  """Calculate first moment of a (kernel) image.
651  """
652  y, x = arr.shape
653  sarr = arr*arr
654  xarr = np.asarray([[el for el in range(x)] for el2 in range(y)])
655  yarr = np.asarray([[el2 for el in range(x)] for el2 in range(y)])
656  narr = xarr*sarr
657  sarrSum = sarr.sum()
658  centx = narr.sum()/sarrSum
659  narr = yarr*sarr
660  centy = narr.sum()/sarrSum
661  return centx, centy
662 
663 
664 def calcWidth(arr, centx, centy):
665  """Calculate second moment of a (kernel) image.
666  """
667  y, x = arr.shape
668  # Square the flux so we don't have to deal with negatives
669  sarr = arr*arr
670  xarr = np.asarray([[el for el in range(x)] for el2 in range(y)])
671  yarr = np.asarray([[el2 for el in range(x)] for el2 in range(y)])
672  narr = sarr*np.power((xarr - centx), 2.)
673  sarrSum = sarr.sum()
674  xstd = np.sqrt(narr.sum()/sarrSum)
675  narr = sarr*np.power((yarr - centy), 2.)
676  ystd = np.sqrt(narr.sum()/sarrSum)
677  return xstd, ystd
678 
679 
680 def printSkyDiffs(sources, wcs):
681  """Print differences in sky coordinates.
682 
683  The difference is that between the source Position and its Centroid mapped
684  through Wcs.
685  """
686  for s in sources:
687  sCentroid = s.getCentroid()
688  sPosition = s.getCoord().getPosition(afwGeom.degrees)
689  dra = 3600*(sPosition.getX() - wcs.pixelToSky(sCentroid).getPosition(afwGeom.degrees).getX())/0.2
690  ddec = 3600*(sPosition.getY() - wcs.pixelToSky(sCentroid).getPosition(afwGeom.degrees).getY())/0.2
691  if np.isfinite(dra) and np.isfinite(ddec):
692  print(dra, ddec)
693 
694 
695 def makeRegions(sources, outfilename, wcs=None):
696  """Create regions file for display from input source list.
697  """
698  fh = open(outfilename, "w")
699  fh.write("global color=red font=\"helvetica 10 normal\" "
700  "select=1 highlite=1 edit=1 move=1 delete=1 include=1 fixed=0 source\nfk5\n")
701  for s in sources:
702  if wcs:
703  (ra, dec) = wcs.pixelToSky(s.getCentroid()).getPosition(afwGeom.degrees)
704  else:
705  (ra, dec) = s.getCoord().getPosition(afwGeom.degrees)
706  if np.isfinite(ra) and np.isfinite(dec):
707  fh.write("circle(%f,%f,2\")\n"%(ra, dec))
708  fh.flush()
709  fh.close()
710 
711 
712 def showSourceSetSky(sSet, wcs, xy0, frame=0, ctype=afwDisplay.GREEN, symb="+", size=2):
713  """Draw the (RA, Dec) positions of a set of Sources. Image has the XY0.
714  """
715  disp = afwDisplay.Display(frame=frame)
716  with disp.Buffering():
717  for s in sSet:
718  (xc, yc) = wcs.skyToPixel(s.getCoord().getRa(), s.getCoord().getDec())
719  xc -= xy0[0]
720  yc -= xy0[1]
721  disp.dot(symb, xc, yc, ctype=ctype, size=size)
722 
723 
724 def plotWhisker(results, newWcs):
725  """Plot whisker diagram of astromeric offsets between results.matches.
726  """
727  refCoordKey = results.matches[0].first.getTable().getCoordKey()
728  inCentroidKey = results.matches[0].second.getTable().getCentroidKey()
729  positions = [m.first.get(refCoordKey) for m in results.matches]
730  residuals = [m.first.get(refCoordKey).getOffsetFrom(
731  newWcs.pixelToSky(m.second.get(inCentroidKey))) for
732  m in results.matches]
733  import matplotlib.pyplot as plt
734  fig = plt.figure()
735  sp = fig.add_subplot(1, 1, 0)
736  xpos = [x[0].asDegrees() for x in positions]
737  ypos = [x[1].asDegrees() for x in positions]
738  xpos.append(0.02*(max(xpos) - min(xpos)) + min(xpos))
739  ypos.append(0.98*(max(ypos) - min(ypos)) + min(ypos))
740  xidxs = np.isfinite(xpos)
741  yidxs = np.isfinite(ypos)
742  X = np.asarray(xpos)[xidxs]
743  Y = np.asarray(ypos)[yidxs]
744  distance = [x[1].asArcseconds() for x in residuals]
745  distance.append(0.2)
746  distance = np.asarray(distance)[xidxs]
747  # NOTE: This assumes that the bearing is measured positive from +RA through North.
748  # From the documentation this is not clear.
749  bearing = [x[0].asRadians() for x in residuals]
750  bearing.append(0)
751  bearing = np.asarray(bearing)[xidxs]
752  U = (distance*np.cos(bearing))
753  V = (distance*np.sin(bearing))
754  sp.quiver(X, Y, U, V)
755  sp.set_title("WCS Residual")
756  plt.show()
757 
758 
759 class DipoleTestImage(object):
760  """Utility class for dipole measurement testing.
761 
762  Generate an image with simulated dipoles and noise; store the original
763  "pre-subtraction" images and catalogs as well.
764  Used to generate test data for DMTN-007 (http://dmtn-007.lsst.io).
765  """
766 
767  def __init__(self, w=101, h=101, xcenPos=[27.], ycenPos=[25.], xcenNeg=[23.], ycenNeg=[25.],
768  psfSigma=2., flux=[30000.], fluxNeg=None, noise=10., gradientParams=None):
769  self.w = w
770  self.h = h
771  self.xcenPos = xcenPos
772  self.ycenPos = ycenPos
773  self.xcenNeg = xcenNeg
774  self.ycenNeg = ycenNeg
775  self.psfSigma = psfSigma
776  self.flux = flux
777  self.fluxNeg = fluxNeg
778  if fluxNeg is None:
779  self.fluxNeg = self.flux
780  self.noise = noise
781  self.gradientParams = gradientParams
782  self._makeDipoleImage()
783 
784  def _makeDipoleImage(self):
785  """Generate an exposure and catalog with the given dipole source(s).
786  """
787  # Must seed the pos/neg images with different values to ensure they get different noise realizations
788  posImage, posCatalog = self._makeStarImage(
789  xc=self.xcenPos, yc=self.ycenPos, flux=self.flux, randomSeed=111)
790 
791  negImage, negCatalog = self._makeStarImage(
792  xc=self.xcenNeg, yc=self.ycenNeg, flux=self.fluxNeg, randomSeed=222)
793 
794  dipole = posImage.clone()
795  di = dipole.getMaskedImage()
796  di -= negImage.getMaskedImage()
797 
798  # Carry through pos/neg detection masks to new planes in diffim
799  dm = di.getMask()
800  posDetectedBits = posImage.getMaskedImage().getMask().getArray() == dm.getPlaneBitMask("DETECTED")
801  negDetectedBits = negImage.getMaskedImage().getMask().getArray() == dm.getPlaneBitMask("DETECTED")
802  pos_det = dm.addMaskPlane("DETECTED_POS") # new mask plane -- different from "DETECTED"
803  neg_det = dm.addMaskPlane("DETECTED_NEG") # new mask plane -- different from "DETECTED_NEGATIVE"
804  dma = dm.getArray()
805  # set the two custom mask planes to these new masks
806  dma[:, :] = posDetectedBits*pos_det + negDetectedBits*neg_det
807  self.diffim, self.posImage, self.posCatalog, self.negImage, self.negCatalog \
808  = dipole, posImage, posCatalog, negImage, negCatalog
809 
810  def _makeStarImage(self, xc=[15.3], yc=[18.6], flux=[2500], schema=None, randomSeed=None):
811  """Generate an exposure and catalog with the given stellar source(s).
812  """
813  from lsst.meas.base.tests import TestDataset
814  bbox = afwGeom.Box2I(afwGeom.Point2I(0, 0), afwGeom.Point2I(self.w - 1, self.h - 1))
815  dataset = TestDataset(bbox, psfSigma=self.psfSigma, threshold=1.)
816 
817  for i in range(len(xc)):
818  dataset.addSource(instFlux=flux[i], centroid=afwGeom.Point2D(xc[i], yc[i]))
819 
820  if schema is None:
821  schema = TestDataset.makeMinimalSchema()
822  exposure, catalog = dataset.realize(noise=self.noise, schema=schema, randomSeed=randomSeed)
823 
824  if self.gradientParams is not None:
825  y, x = np.mgrid[:self.w, :self.h]
826  gp = self.gradientParams
827  gradient = gp[0] + gp[1]*x + gp[2]*y
828  if len(self.gradientParams) > 3: # it includes a set of 2nd-order polynomial params
829  gradient += gp[3]*x*y + gp[4]*x*x + gp[5]*y*y
830  imgArr = exposure.getMaskedImage().getArrays()[0]
831  imgArr += gradient
832 
833  return exposure, catalog
834 
835  def fitDipoleSource(self, source, **kwds):
836  alg = DipoleFitAlgorithm(self.diffim, self.posImage, self.negImage)
837  fitResult = alg.fitDipole(source, **kwds)
838  return fitResult
839 
840  def detectDipoleSources(self, doMerge=True, diffim=None, detectSigma=5.5, grow=3, minBinSize=32):
841  """Utility function for detecting dipoles.
842 
843  Detect pos/neg sources in the diffim, then merge them. A
844  bigger "grow" parameter leads to a larger footprint which
845  helps with dipole measurement for faint dipoles.
846 
847  Parameters
848  ----------
849  doMerge : `bool`
850  Whether to merge the positive and negagive detections into a single
851  source table.
852  diffim : `lsst.afw.image.exposure.exposure.ExposureF`
853  Difference image on which to perform detection.
854  detectSigma : `float`
855  Threshold for object detection.
856  grow : `int`
857  Number of pixels to grow the footprints before merging.
858  minBinSize : `int`
859  Minimum bin size for the background (re)estimation (only applies if
860  the default leads to min(nBinX, nBinY) < fit order so the default
861  config parameter needs to be decreased, but not to a value smaller
862  than ``minBinSize``, in which case the fitting algorithm will take
863  over and decrease the fit order appropriately.)
864 
865  Returns
866  -------
867  sources : `lsst.afw.table.SourceCatalog`
868  If doMerge=True, the merged source catalog is returned OR
869  detectTask : `lsst.meas.algorithms.SourceDetectionTask`
870  schema : `lsst.afw.table.Schema`
871  If doMerge=False, the source detection task and its schema are
872  returned.
873  """
874  if diffim is None:
875  diffim = self.diffim
876 
877  # Start with a minimal schema - only the fields all SourceCatalogs need
878  schema = afwTable.SourceTable.makeMinimalSchema()
879 
880  # Customize the detection task a bit (optional)
881  detectConfig = measAlg.SourceDetectionConfig()
882  detectConfig.returnOriginalFootprints = False # should be the default
883 
884  psfSigma = diffim.getPsf().computeShape().getDeterminantRadius()
885 
886  # code from imageDifference.py:
887  detectConfig.thresholdPolarity = "both"
888  detectConfig.thresholdValue = detectSigma
889  # detectConfig.nSigmaToGrow = psfSigma
890  detectConfig.reEstimateBackground = True # if False, will fail often for faint sources on gradients?
891  detectConfig.thresholdType = "pixel_stdev"
892  # Test images are often quite small, so may need to adjust background binSize
893  while ((min(diffim.getWidth(), diffim.getHeight()))//detectConfig.background.binSize <
894  detectConfig.background.approxOrderX and detectConfig.background.binSize > minBinSize):
895  detectConfig.background.binSize = max(minBinSize, detectConfig.background.binSize//2)
896 
897  # Create the detection task. We pass the schema so the task can declare a few flag fields
898  detectTask = measAlg.SourceDetectionTask(schema, config=detectConfig)
899 
900  table = afwTable.SourceTable.make(schema)
901  catalog = detectTask.makeSourceCatalog(table, diffim, sigma=psfSigma)
902 
903  # Now do the merge.
904  if doMerge:
905  fpSet = catalog.fpSets.positive
906  fpSet.merge(catalog.fpSets.negative, grow, grow, False)
907  sources = afwTable.SourceCatalog(table)
908  fpSet.makeSources(sources)
909 
910  return sources
911 
912  else:
913  return detectTask, schema
def showDiaSources(sources, exposure, isFlagged, isDipole, frame=None)
Definition: utils.py:105
def makeRegions(sources, outfilename, wcs=None)
Definition: utils.py:695
def __init__(self, w=101, h=101, xcenPos=[27.], ycenPos=[25.], xcenNeg=[23.], ycenNeg=[25.], psfSigma=2., flux=[30000.], fluxNeg=None, noise=10., gradientParams=None)
Definition: utils.py:768
MaskedImage< ImagePixelT, MaskPixelT, VariancePixelT > * makeMaskedImage(typename std::shared_ptr< Image< ImagePixelT >> image, typename std::shared_ptr< Mask< MaskPixelT >> mask=Mask< MaskPixelT >(), typename std::shared_ptr< Image< VariancePixelT >> variance=Image< VariancePixelT >())
def showKernelSpatialCells(maskedIm, kernelCellSet, showChi2=False, symb="o", ctype=None, ctypeUnused=None, ctypeBad=None, size=3, frame=None, title="Spatial Cells")
Definition: utils.py:68
std::shared_ptr< Exposure< ImagePixelT, MaskPixelT, VariancePixelT > > makeExposure(MaskedImage< ImagePixelT, MaskPixelT, VariancePixelT > &mimage, std::shared_ptr< geom::SkyWcs const > wcs=std::shared_ptr< geom::SkyWcs const >())
def _makeStarImage(self, xc=[15.3], yc=[18.6], flux=[2500], schema=None, randomSeed=None)
Definition: utils.py:810
def plotKernelSpatialModel(kernel, kernelCellSet, showBadCandidates=True, numSample=128, keepPlots=True, maxCoeff=10)
Definition: utils.py:274
def plotPixelResiduals(exposure, warpedTemplateExposure, diffExposure, kernelCellSet, kernel, background, testSources, config, origVariance=False, nptsFull=1e6, keepPlots=True, titleFs=14)
Definition: utils.py:507
def printSkyDiffs(sources, wcs)
Definition: utils.py:680
def showKernelCandidates(kernelCellSet, kernel, background, frame=None, showBadCandidates=True, resids=False, kernels=False)
Definition: utils.py:142
def detectDipoleSources(self, doMerge=True, diffim=None, detectSigma=5.5, grow=3, minBinSize=32)
Definition: utils.py:840
def calcCentroid(arr)
Definition: utils.py:649
def calcWidth(arr, centx, centy)
Definition: utils.py:664
def showSourceSet(sSet, xy0=(0, 0), frame=0, ctype=afwDisplay.GREEN, symb="+", size=2)
Definition: utils.py:46
def plotWhisker(results, newWcs)
Definition: utils.py:724
def showKernelMosaic(bbox, kernel, nx=7, ny=None, frame=None, title=None, showCenter=True, showEllipticity=True)
Definition: utils.py:424
def showKernelBasis(kernel, frame=None)
Definition: utils.py:257
def fitDipoleSource(self, source, kwds)
Definition: utils.py:835
def showSourceSetSky(sSet, wcs, xy0, frame=0, ctype=afwDisplay.GREEN, symb="+", size=2)
Definition: utils.py:712