lsst.pipe.tasks  13.0-66-gfbf2f2ce+5
matchBackgrounds.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 # Copyright 2008, 2009, 2010 LSST Corporation.
4 #
5 # This product includes software developed by the
6 # LSST Project (http://www.lsst.org/).
7 #
8 # This program is free software: you can redistribute it and/or modify
9 # it under the terms of the GNU General Public License as published by
10 # the Free Software Foundation, either version 3 of the License, or
11 # (at your option) any later version.
12 #
13 # This program is distributed in the hope that it will be useful,
14 # but WITHOUT ANY WARRANTY; without even the implied warranty of
15 # MERCHANTABILIY or FITNESS FOR A PARTICULAR PURPOSE. See the
16 # GNU General Public License for more details.
17 #
18 # You should have received a copy of the LSST License Statement and
19 # the GNU General Public License along with this program. If not,
20 # see <http://www.lsstcorp.org/LegalNotices/>.
21 
22 from __future__ import absolute_import, division, print_function
23 from builtins import zip
24 from builtins import range
25 from builtins import object
26 import numpy
27 import lsst.afw.image as afwImage
28 import lsst.afw.math as afwMath
29 import lsst.afw.geom as afwGeom
30 import lsst.pex.config as pexConfig
31 import lsst.pipe.base as pipeBase
32 import lsstDebug
33 
34 
35 class MatchBackgroundsConfig(pexConfig.Config):
36 
37  usePolynomial = pexConfig.Field(
38  dtype=bool,
39  doc="Fit background difference with Chebychev polynomial interpolation "
40  "(using afw.math.Approximate)? If False, fit with spline interpolation using afw.math.Background",
41  default=False
42  )
43  order = pexConfig.Field(
44  dtype=int,
45  doc="Order of Chebyshev polynomial background model. Ignored if usePolynomial False",
46  default=8
47  )
48  badMaskPlanes = pexConfig.ListField(
49  doc="Names of mask planes to ignore while estimating the background",
50  dtype=str, default=["NO_DATA", "DETECTED", "DETECTED_NEGATIVE", "SAT", "BAD", "INTRP", "CR"],
51  itemCheck=lambda x: x in afwImage.Mask().getMaskPlaneDict(),
52  )
53  gridStatistic = pexConfig.ChoiceField(
54  dtype=str,
55  doc="Type of statistic to estimate pixel value for the grid points",
56  default="MEAN",
57  allowed={
58  "MEAN": "mean",
59  "MEDIAN": "median",
60  "MEANCLIP": "clipped mean"
61  }
62  )
63  undersampleStyle = pexConfig.ChoiceField(
64  doc="Behaviour if there are too few points in grid for requested interpolation style. "
65  "Note: INCREASE_NXNYSAMPLE only allowed for usePolynomial=True.",
66  dtype=str,
67  default="REDUCE_INTERP_ORDER",
68  allowed={
69  "THROW_EXCEPTION": "throw an exception if there are too few points",
70  "REDUCE_INTERP_ORDER": "use an interpolation style with a lower order.",
71  "INCREASE_NXNYSAMPLE": "Increase the number of samples used to make the interpolation grid.",
72  }
73  )
74  binSize = pexConfig.Field(
75  doc="Bin size for gridding the difference image and fitting a spatial model",
76  dtype=int,
77  default=256
78  )
79  interpStyle = pexConfig.ChoiceField(
80  dtype=str,
81  doc="Algorithm to interpolate the background values; ignored if usePolynomial is True"
82  "Maps to an enum; see afw.math.Background",
83  default="AKIMA_SPLINE",
84  allowed={
85  "CONSTANT": "Use a single constant value",
86  "LINEAR": "Use linear interpolation",
87  "NATURAL_SPLINE": "cubic spline with zero second derivative at endpoints",
88  "AKIMA_SPLINE": "higher-level nonlinear spline that is more robust to outliers",
89  "NONE": "No background estimation is to be attempted",
90  }
91  )
92  numSigmaClip = pexConfig.Field(
93  dtype=int,
94  doc="Sigma for outlier rejection; ignored if gridStatistic != 'MEANCLIP'.",
95  default=3
96  )
97  numIter = pexConfig.Field(
98  dtype=int,
99  doc="Number of iterations of outlier rejection; ignored if gridStatistic != 'MEANCLIP'.",
100  default=2
101  )
102  bestRefWeightCoverage = pexConfig.RangeField(
103  dtype=float,
104  doc="Weight given to coverage (number of pixels that overlap with patch), "
105  "when calculating best reference exposure. Higher weight prefers exposures with high coverage."
106  "Ignored when reference visit is supplied",
107  default=0.4,
108  min=0., max=1.
109  )
110  bestRefWeightVariance = pexConfig.RangeField(
111  dtype=float,
112  doc="Weight given to image variance when calculating best reference exposure. "
113  "Higher weight prefers exposures with low image variance. Ignored when reference visit is supplied",
114  default=0.4,
115  min=0., max=1.
116  )
117  bestRefWeightLevel = pexConfig.RangeField(
118  dtype=float,
119  doc="Weight given to mean background level when calculating best reference exposure. "
120  "Higher weight prefers exposures with low mean background level. "
121  "Ignored when reference visit is supplied.",
122  default=0.2,
123  min=0., max=1.
124  )
125  approxWeighting = pexConfig.Field(
126  dtype=bool,
127  doc=("Use inverse-variance weighting when approximating background offset model? " +
128  "This will fail when the background offset is constant " +
129  "(this is usually only the case in testing with artificial images)." +
130  "(usePolynomial=True)"),
131  default=True,
132  )
133  gridStdevEpsilon = pexConfig.RangeField(
134  dtype=float,
135  doc="Tolerance on almost zero standard deviation in a background-offset grid bin. "
136  "If all bins have a standard deviation below this value, the background offset model "
137  "is approximated without inverse-variance weighting. (usePolynomial=True)",
138  default=1e-8,
139  min=0.
140  )
141 
142 
143 class MatchBackgroundsTask(pipeBase.Task):
144  ConfigClass = MatchBackgroundsConfig
145  _DefaultName = "matchBackgrounds"
146 
147  def __init__(self, *args, **kwargs):
148  pipeBase.Task.__init__(self, *args, **kwargs)
149 
150  self.sctrl = afwMath.StatisticsControl()
151  self.sctrl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.badMaskPlanes))
152  self.sctrl.setNanSafe(True)
153 
154  @pipeBase.timeMethod
155  def run(self, expRefList, expDatasetType, imageScalerList=None, refExpDataRef=None, refImageScaler=None):
156  """Match the backgrounds of a list of coadd temp exposures to a reference coadd temp exposure.
157 
158  Choose a refExpDataRef automatically if none supplied.
159 
160  @param[in] expRefList: list of data references to science exposures to be background-matched;
161  all exposures must exist.
162  @param[in] expDatasetType: dataset type of exposures, e.g. 'goodSeeingCoadd_tempExp'
163  @param[in] imageScalerList: list of image scalers (coaddUtils.ImageScaler);
164  if None then the images are not scaled
165  @param[in] refExpDataRef: data reference for the reference exposure.
166  If None, then this task selects the best exposures from expRefList.
167  if not None then must be one of the exposures in expRefList.
168  @param[in] refImageScaler: image scaler for reference image;
169  ignored if refExpDataRef is None, else scaling is not performed if None
170 
171  @return: a pipBase.Struct containing these fields:
172  - backgroundInfoList: a list of pipeBase.Struct, one per exposure in expRefList,
173  each of which contains these fields:
174  - isReference: this is the reference exposure (only one returned Struct will
175  contain True for this value, unless the ref exposure is listed multiple times)
176  - backgroundModel: differential background model (afw.Math.Background or afw.Math.Approximate).
177  Add this to the science exposure to match the reference exposure.
178  - fitRMS: rms of the fit. This is the sqrt(mean(residuals**2)).
179  - matchedMSE: the MSE of the reference and matched images: mean((refImage - matchedSciImage)**2);
180  should be comparable to difference image's mean variance.
181  - diffImVar: the mean variance of the difference image.
182  All fields except isReference will be None if isReference True or the fit failed.
183 
184  @warning: all exposures must exist on disk
185  """
186 
187  numExp = len(expRefList)
188  if numExp < 1:
189  raise pipeBase.TaskError("No exposures to match")
190 
191  if expDatasetType is None:
192  raise pipeBase.TaskError("Must specify expDatasetType")
193 
194  if imageScalerList is None:
195  self.log.info("imageScalerList is None; no scaling will be performed")
196  imageScalerList = [None] * numExp
197 
198  if len(expRefList) != len(imageScalerList):
199  raise RuntimeError("len(expRefList) = %s != %s = len(imageScalerList)" %
200  (len(expRefList), len(imageScalerList)))
201 
202  refInd = None
203  if refExpDataRef is None:
204  # select the best reference exposure from expRefList
205  refInd = self.selectRefExposure(
206  expRefList=expRefList,
207  imageScalerList=imageScalerList,
208  expDatasetType=expDatasetType,
209  )
210  refExpDataRef = expRefList[refInd]
211  refImageScaler = imageScalerList[refInd]
212 
213  # refIndSet is the index of all exposures in expDataList that match the reference.
214  # It is used to avoid background-matching an exposure to itself. It is a list
215  # because it is possible (though unlikely) that expDataList will contain duplicates.
216  expKeyList = refExpDataRef.butlerSubset.butler.getKeys(expDatasetType)
217  refMatcher = DataRefMatcher(refExpDataRef.butlerSubset.butler, expDatasetType)
218  refIndSet = set(refMatcher.matchList(ref0=refExpDataRef, refList=expRefList))
219 
220  if refInd is not None and refInd not in refIndSet:
221  raise RuntimeError("Internal error: selected reference %s not found in expRefList")
222 
223  refExposure = refExpDataRef.get(expDatasetType, immediate=True)
224  if refImageScaler is not None:
225  refMI = refExposure.getMaskedImage()
226  refImageScaler.scaleMaskedImage(refMI)
227 
228  debugIdKeyList = tuple(set(expKeyList) - set(['tract', 'patch']))
229 
230  self.log.info("Matching %d Exposures" % (numExp))
231 
232  backgroundInfoList = []
233  for ind, (toMatchRef, imageScaler) in enumerate(zip(expRefList, imageScalerList)):
234  if ind in refIndSet:
235  backgroundInfoStruct = pipeBase.Struct(
236  isReference=True,
237  backgroundModel=None,
238  fitRMS=0.0,
239  matchedMSE=None,
240  diffImVar=None,
241  )
242  else:
243  self.log.info("Matching background of %s to %s" % (toMatchRef.dataId, refExpDataRef.dataId))
244  try:
245  toMatchExposure = toMatchRef.get(expDatasetType, immediate=True)
246  if imageScaler is not None:
247  toMatchMI = toMatchExposure.getMaskedImage()
248  imageScaler.scaleMaskedImage(toMatchMI)
249  # store a string specifying the visit to label debug plot
250  self.debugDataIdString = ''.join([str(toMatchRef.dataId[vk]) for vk in debugIdKeyList])
251  backgroundInfoStruct = self.matchBackgrounds(
252  refExposure=refExposure,
253  sciExposure=toMatchExposure,
254  )
255  backgroundInfoStruct.isReference = False
256  except Exception as e:
257  self.log.warn("Failed to fit background %s: %s" % (toMatchRef.dataId, e))
258  backgroundInfoStruct = pipeBase.Struct(
259  isReference=False,
260  backgroundModel=None,
261  fitRMS=None,
262  matchedMSE=None,
263  diffImVar=None,
264  )
265 
266  backgroundInfoList.append(backgroundInfoStruct)
267 
268  return pipeBase.Struct(
269  backgroundInfoList=backgroundInfoList)
270 
271  @pipeBase.timeMethod
272  def selectRefExposure(self, expRefList, imageScalerList, expDatasetType):
273  """Find best exposure to use as the reference exposure
274 
275  Calculate an appropriate reference exposure by minimizing a cost function that penalizes
276  high variance, high background level, and low coverage. Use the following config parameters:
277  - bestRefWeightCoverage
278  - bestRefWeightVariance
279  - bestRefWeightLevel
280 
281  @param[in] expRefList: list of data references to exposures.
282  Retrieves dataset type specified by expDatasetType.
283  If an exposure is not found, it is skipped with a warning.
284  @param[in] imageScalerList: list of image scalers (coaddUtils.ImageScaler);
285  must be the same length as expRefList
286  @param[in] expDatasetType: dataset type of exposure: e.g. 'goodSeeingCoadd_tempExp'
287 
288  @return: index of best exposure
289 
290  @raise pipeBase.TaskError if none of the exposures in expRefList are found.
291  """
292  self.log.info("Calculating best reference visit")
293  varList = []
294  meanBkgdLevelList = []
295  coverageList = []
296 
297  if len(expRefList) != len(imageScalerList):
298  raise RuntimeError("len(expRefList) = %s != %s = len(imageScalerList)" %
299  (len(expRefList), len(imageScalerList)))
300 
301  for expRef, imageScaler in zip(expRefList, imageScalerList):
302  exposure = expRef.get(expDatasetType, immediate=True)
303  maskedImage = exposure.getMaskedImage()
304  if imageScaler is not None:
305  try:
306  imageScaler.scaleMaskedImage(maskedImage)
307  except:
308  # need to put a place holder in Arr
309  varList.append(numpy.nan)
310  meanBkgdLevelList.append(numpy.nan)
311  coverageList.append(numpy.nan)
312  continue
313  statObjIm = afwMath.makeStatistics(maskedImage.getImage(), maskedImage.getMask(),
314  afwMath.MEAN | afwMath.NPOINT | afwMath.VARIANCE, self.sctrl)
315  meanVar, meanVarErr = statObjIm.getResult(afwMath.VARIANCE)
316  meanBkgdLevel, meanBkgdLevelErr = statObjIm.getResult(afwMath.MEAN)
317  npoints, npointsErr = statObjIm.getResult(afwMath.NPOINT)
318  varList.append(meanVar)
319  meanBkgdLevelList.append(meanBkgdLevel)
320  coverageList.append(npoints)
321  if not coverageList:
322  raise pipeBase.TaskError(
323  "None of the candidate %s exist; cannot select best reference exposure" % (expDatasetType,))
324 
325  # Normalize metrics to range from 0 to 1
326  varArr = numpy.array(varList)/numpy.nanmax(varList)
327  meanBkgdLevelArr = numpy.array(meanBkgdLevelList)/numpy.nanmax(meanBkgdLevelList)
328  coverageArr = numpy.nanmin(coverageList)/numpy.array(coverageList)
329 
330  costFunctionArr = self.config.bestRefWeightVariance * varArr
331  costFunctionArr += self.config.bestRefWeightLevel * meanBkgdLevelArr
332  costFunctionArr += self.config.bestRefWeightCoverage * coverageArr
333  return numpy.nanargmin(costFunctionArr)
334 
335  @pipeBase.timeMethod
336  def matchBackgrounds(self, refExposure, sciExposure):
337  """
338  Match science exposure's background level to that of reference exposure.
339 
340  Process creates a difference image of the reference exposure minus the science exposure, and then
341  generates an afw.math.Background object. It assumes (but does not require/check) that the mask plane
342  already has detections set. If detections have not been set/masked, sources will bias the
343  background estimation.
344  The 'background' of the difference image is smoothed by spline interpolation (by the Background class)
345  or by polynomial interpolation by the Approximate class. This model of difference image is added to the
346  science exposure in memory.
347  Fit diagnostics are also calculated and returned.
348 
349  @param[in] refExposure: reference exposure
350  @param[in,out] sciExposure: science exposure; modified by changing the background level
351  to match that of the reference exposure
352  @returns a pipBase.Struct with fields:
353  - backgroundModel: an afw.math.Approximate or an afw.math.Background.
354  - fitRMS: rms of the fit. This is the sqrt(mean(residuals**2)).
355  - matchedMSE: the MSE of the reference and matched images: mean((refImage - matchedSciImage)**2);
356  should be comparable to difference image's mean variance.
357  - diffImVar: the mean variance of the difference image.
358  """
359 
360  if lsstDebug.Info(__name__).savefits:
361  refExposure.writeFits(lsstDebug.Info(__name__).figpath + 'refExposure.fits')
362  sciExposure.writeFits(lsstDebug.Info(__name__).figpath + 'sciExposure.fits')
363 
364  # Check Configs for polynomials:
365  if self.config.usePolynomial:
366  x, y = sciExposure.getDimensions()
367  shortSideLength = min(x, y)
368  if shortSideLength < self.config.binSize:
369  raise ValueError("%d = config.binSize > shorter dimension = %d" % (self.config.binSize,
370  shortSideLength))
371  npoints = shortSideLength // self.config.binSize
372  if shortSideLength % self.config.binSize != 0:
373  npoints += 1
374 
375  if self.config.order > npoints - 1:
376  raise ValueError("%d = config.order > npoints - 1 = %d" % (self.config.order, npoints - 1))
377 
378  # Check that exposures are same shape
379  if (sciExposure.getDimensions() != refExposure.getDimensions()):
380  wSci, hSci = sciExposure.getDimensions()
381  wRef, hRef = refExposure.getDimensions()
382  raise RuntimeError(
383  "Exposures are different dimensions. sci:(%i, %i) vs. ref:(%i, %i)" %
384  (wSci, hSci, wRef, hRef))
385 
386  statsFlag = getattr(afwMath, self.config.gridStatistic)
387  self.sctrl.setNumSigmaClip(self.config.numSigmaClip)
388  self.sctrl.setNumIter(self.config.numIter)
389 
390  im = refExposure.getMaskedImage()
391  diffMI = im.Factory(im, True)
392  diffMI -= sciExposure.getMaskedImage()
393 
394  width = diffMI.getWidth()
395  height = diffMI.getHeight()
396  nx = width // self.config.binSize
397  if width % self.config.binSize != 0:
398  nx += 1
399  ny = height // self.config.binSize
400  if height % self.config.binSize != 0:
401  ny += 1
402 
403  bctrl = afwMath.BackgroundControl(nx, ny, self.sctrl, statsFlag)
404  bctrl.setUndersampleStyle(self.config.undersampleStyle)
405  bctrl.setInterpStyle(self.config.interpStyle)
406 
407  bkgd = afwMath.makeBackground(diffMI, bctrl)
408 
409  # Some config and input checks if config.usePolynomial:
410  # 1) Check that order/bin size make sense:
411  # 2) Change binsize or order if underconstrained.
412  if self.config.usePolynomial:
413  order = self.config.order
414  bgX, bgY, bgZ, bgdZ = self._gridImage(diffMI, self.config.binSize, statsFlag)
415  minNumberGridPoints = min(len(set(bgX)), len(set(bgY)))
416  if len(bgZ) == 0:
417  raise ValueError("No overlap with reference. Nothing to match")
418  elif minNumberGridPoints <= self.config.order:
419  # must either lower order or raise number of bins or throw exception
420  if self.config.undersampleStyle == "THROW_EXCEPTION":
421  raise ValueError("Image does not cover enough of ref image for order and binsize")
422  elif self.config.undersampleStyle == "REDUCE_INTERP_ORDER":
423  self.log.warn("Reducing order to %d"%(minNumberGridPoints - 1))
424  order = minNumberGridPoints - 1
425  elif self.config.undersampleStyle == "INCREASE_NXNYSAMPLE":
426  newBinSize = (minNumberGridPoints*self.config.binSize) // (self.config.order + 1)
427  bctrl.setNxSample(newBinSize)
428  bctrl.setNySample(newBinSize)
429  bkgd = afwMath.makeBackground(diffMI, bctrl) # do over
430  self.log.warn("Decreasing binsize to %d"%(newBinSize))
431 
432  # If there is no variance in any image pixels, do not weight bins by inverse variance
433  isUniformImageDiff = not numpy.any(bgdZ > self.config.gridStdevEpsilon)
434  weightByInverseVariance = False if isUniformImageDiff else self.config.approxWeighting
435 
436  # Add offset to sciExposure
437  try:
438  if self.config.usePolynomial:
439  actrl = afwMath.ApproximateControl(afwMath.ApproximateControl.CHEBYSHEV,
440  order, order, weightByInverseVariance)
441  undersampleStyle = getattr(afwMath, self.config.undersampleStyle)
442  approx = bkgd.getApproximate(actrl, undersampleStyle)
443  bkgdImage = approx.getImage()
444  else:
445  bkgdImage = bkgd.getImageF()
446  except Exception as e:
447  raise RuntimeError("Background/Approximation failed to interp image %s: %s" % (
448  self.debugDataIdString, e))
449 
450  sciMI = sciExposure.getMaskedImage()
451  sciMI += bkgdImage
452  del sciMI
453 
454  # Need RMS from fit: 2895 will replace this:
455  rms = 0.0
456  X, Y, Z, dZ = self._gridImage(diffMI, self.config.binSize, statsFlag)
457  x0, y0 = diffMI.getXY0()
458  modelValueArr = numpy.empty(len(Z))
459  for i in range(len(X)):
460  modelValueArr[i] = bkgdImage.get(int(X[i]-x0), int(Y[i]-y0))
461  resids = Z - modelValueArr
462  rms = numpy.sqrt(numpy.mean(resids[~numpy.isnan(resids)]**2))
463 
464  if lsstDebug.Info(__name__).savefits:
465  sciExposure.writeFits(lsstDebug.Info(__name__).figpath + 'sciMatchedExposure.fits')
466 
467  if lsstDebug.Info(__name__).savefig:
468  bbox = afwGeom.Box2D(refExposure.getMaskedImage().getBBox())
469  try:
470  self._debugPlot(X, Y, Z, dZ, bkgdImage, bbox, modelValueArr, resids)
471  except Exception as e:
472  self.log.warn('Debug plot not generated: %s'%(e))
473 
474  meanVar = afwMath.makeStatistics(diffMI.getVariance(), diffMI.getMask(),
475  afwMath.MEANCLIP, self.sctrl).getValue()
476 
477  diffIm = diffMI.getImage()
478  diffIm -= bkgdImage # diffMI should now have a mean ~ 0
479  del diffIm
480  mse = afwMath.makeStatistics(diffMI, afwMath.MEANSQUARE, self.sctrl).getValue()
481 
482  outBkgd = approx if self.config.usePolynomial else bkgd
483 
484  return pipeBase.Struct(
485  backgroundModel=outBkgd,
486  fitRMS=rms,
487  matchedMSE=mse,
488  diffImVar=meanVar)
489 
490  def _debugPlot(self, X, Y, Z, dZ, modelImage, bbox, model, resids):
491  """Generate a plot showing the background fit and residuals.
492 
493  It is called when lsstDebug.Info(__name__).savefig = True
494  Saves the fig to lsstDebug.Info(__name__).figpath
495  Displays on screen if lsstDebug.Info(__name__).display = True
496 
497  @param X: array of x positions
498  @param Y: array of y positions
499  @param Z: array of the grid values that were interpolated
500  @param dZ: array of the error on the grid values
501  @param modelImage: image ofthe model of the fit
502  @param model: array of len(Z) containing the grid values predicted by the model
503  @param resids: Z - model
504  """
505  import matplotlib.pyplot as plt
506  import matplotlib.colors
507  from mpl_toolkits.axes_grid1 import ImageGrid
508  zeroIm = afwImage.MaskedImageF(afwGeom.Box2I(bbox))
509  zeroIm += modelImage
510  x0, y0 = zeroIm.getXY0()
511  dx, dy = zeroIm.getDimensions()
512  if len(Z) == 0:
513  self.log.warn("No grid. Skipping plot generation.")
514  else:
515  max, min = numpy.max(Z), numpy.min(Z)
516  norm = matplotlib.colors.normalize(vmax=max, vmin=min)
517  maxdiff = numpy.max(numpy.abs(resids))
518  diffnorm = matplotlib.colors.normalize(vmax=maxdiff, vmin=-maxdiff)
519  rms = numpy.sqrt(numpy.mean(resids**2))
520  fig = plt.figure(1, (8, 6))
521  meanDz = numpy.mean(dZ)
522  grid = ImageGrid(fig, 111, nrows_ncols=(1, 2), axes_pad=0.1,
523  share_all=True, label_mode="L", cbar_mode="each",
524  cbar_size="7%", cbar_pad="2%", cbar_location="top")
525  im = grid[0].imshow(zeroIm.getImage().getArray(),
526  extent=(x0, x0+dx, y0+dy, y0), norm=norm,
527  cmap='Spectral')
528  im = grid[0].scatter(X, Y, c=Z, s=15.*meanDz/dZ, edgecolor='none', norm=norm,
529  marker='o', cmap='Spectral')
530  im2 = grid[1].scatter(X, Y, c=resids, edgecolor='none', norm=diffnorm,
531  marker='s', cmap='seismic')
532  grid.cbar_axes[0].colorbar(im)
533  grid.cbar_axes[1].colorbar(im2)
534  grid[0].axis([x0, x0+dx, y0+dy, y0])
535  grid[1].axis([x0, x0+dx, y0+dy, y0])
536  grid[0].set_xlabel("model and grid")
537  grid[1].set_xlabel("residuals. rms = %0.3f"%(rms))
538  if lsstDebug.Info(__name__).savefig:
539  fig.savefig(lsstDebug.Info(__name__).figpath + self.debugDataIdString + '.png')
540  if lsstDebug.Info(__name__).display:
541  plt.show()
542  plt.clf()
543 
544  def _gridImage(self, maskedImage, binsize, statsFlag):
545  """Private method to grid an image for debugging"""
546  width, height = maskedImage.getDimensions()
547  x0, y0 = maskedImage.getXY0()
548  xedges = numpy.arange(0, width, binsize)
549  yedges = numpy.arange(0, height, binsize)
550  xedges = numpy.hstack((xedges, width)) # add final edge
551  yedges = numpy.hstack((yedges, height)) # add final edge
552 
553  # Use lists/append to protect against the case where
554  # a bin has no valid pixels and should not be included in the fit
555  bgX = []
556  bgY = []
557  bgZ = []
558  bgdZ = []
559 
560  for ymin, ymax in zip(yedges[0:-1], yedges[1:]):
561  for xmin, xmax in zip(xedges[0:-1], xedges[1:]):
562  subBBox = afwGeom.Box2I(afwGeom.PointI(int(x0 + xmin), int(y0 + ymin)),
563  afwGeom.PointI(int(x0 + xmax-1), int(y0 + ymax-1)))
564  subIm = afwImage.MaskedImageF(maskedImage, subBBox, afwImage.PARENT, False)
565  stats = afwMath.makeStatistics(subIm,
566  afwMath.MEAN | afwMath.MEANCLIP | afwMath.MEDIAN |
567  afwMath.NPOINT | afwMath.STDEV,
568  self.sctrl)
569  npoints, _ = stats.getResult(afwMath.NPOINT)
570  if npoints >= 2:
571  stdev, _ = stats.getResult(afwMath.STDEV)
572  if stdev < self.config.gridStdevEpsilon:
573  stdev = self.config.gridStdevEpsilon
574  bgX.append(0.5 * (x0 + xmin + x0 + xmax))
575  bgY.append(0.5 * (y0 + ymin + y0 + ymax))
576  bgdZ.append(stdev/numpy.sqrt(npoints))
577  est, _ = stats.getResult(statsFlag)
578  bgZ.append(est)
579 
580  return numpy.array(bgX), numpy.array(bgY), numpy.array(bgZ), numpy.array(bgdZ)
581 
582 
583 class DataRefMatcher(object):
584  """Match data references for a specified dataset type
585 
586  Note that this is not exact, but should suffice for this task
587  until there is better support for this kind of thing in the butler.
588  """
589 
590  def __init__(self, butler, datasetType):
591  """Construct a DataRefMatcher
592 
593  @param[in] butler
594  @param[in] datasetType: dataset type to match
595  """
596  self._datasetType = datasetType # for diagnostics
597  self._keyNames = butler.getKeys(datasetType)
598 
599  def _makeKey(self, ref):
600  """Return a tuple of values for the specified keyNames
601 
602  @param[in] ref: data reference
603 
604  @raise KeyError if ref.dataId is missing a key in keyNames
605  """
606  return tuple(ref.dataId[key] for key in self._keyNames)
607 
608  def isMatch(self, ref0, ref1):
609  """Return True if ref0 == ref1
610 
611  @param[in] ref0: data ref 0
612  @param[in] ref1: data ref 1
613 
614  @raise KeyError if either ID is missing a key in keyNames
615  """
616  return self._makeKey(ref0) == self._makeKey(ref1)
617 
618  def matchList(self, ref0, refList):
619  """Return a list of indices of matches
620 
621  @return tuple of indices of matches
622 
623  @raise KeyError if any ID is missing a key in keyNames
624  """
625  key0 = self._makeKey(ref0)
626  return tuple(ind for ind, ref in enumerate(refList) if self._makeKey(ref) == key0)
def selectRefExposure(self, expRefList, imageScalerList, expDatasetType)
def matchBackgrounds(self, refExposure, sciExposure)
def _debugPlot(self, X, Y, Z, dZ, modelImage, bbox, model, resids)
def run(self, expRefList, expDatasetType, imageScalerList=None, refExpDataRef=None, refImageScaler=None)
def _gridImage(self, maskedImage, binsize, statsFlag)