lsst.ip.diffim  16.0-22-g847a80f+2
imageDecorrelation.py
Go to the documentation of this file.
1 #
2 # LSST Data Management System
3 # Copyright 2016 AURA/LSST.
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 # MERCHANTABILITY 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 <https://www.lsstcorp.org/LegalNotices/>.
21 #
22 
23 import numpy as np
24 
25 import lsst.afw.image as afwImage
26 import lsst.afw.geom as afwGeom
27 import lsst.meas.algorithms as measAlg
28 import lsst.afw.math as afwMath
29 import lsst.pex.config as pexConfig
30 import lsst.pipe.base as pipeBase
31 import lsst.log
32 
33 from .imageMapReduce import (ImageMapReduceConfig, ImageMapReduceTask,
34  ImageMapper)
35 
36 __all__ = ("DecorrelateALKernelTask", "DecorrelateALKernelConfig",
37  "DecorrelateALKernelMapper", "DecorrelateALKernelMapReduceConfig",
38  "DecorrelateALKernelSpatialConfig", "DecorrelateALKernelSpatialTask")
39 
40 
41 class DecorrelateALKernelConfig(pexConfig.Config):
42  """Configuration parameters for the DecorrelateALKernelTask
43  """
44 
45  ignoreMaskPlanes = pexConfig.ListField(
46  dtype=str,
47  doc="""Mask planes to ignore for sigma-clipped statistics""",
48  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
49  )
50 
51 
52 class DecorrelateALKernelTask(pipeBase.Task):
53  """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference
54 
55  Notes
56  -----
57 
58  Pipe-task that removes the neighboring-pixel covariance in an
59  image difference that are added when the template image is
60  convolved with the Alard-Lupton PSF matching kernel.
61 
62  The image differencing pipeline task @link
63  ip.diffim.psfMatch.PsfMatchTask PSFMatchTask@endlink and @link
64  ip.diffim.psfMatch.PsfMatchConfigAL PSFMatchConfigAL@endlink uses
65  the Alard and Lupton (1998) method for matching the PSFs of the
66  template and science exposures prior to subtraction. The
67  Alard-Lupton method identifies a matching kernel, which is then
68  (typically) convolved with the template image to perform PSF
69  matching. This convolution has the effect of adding covariance
70  between neighboring pixels in the template image, which is then
71  added to the image difference by subtraction.
72 
73  The pixel covariance may be corrected by whitening the noise of
74  the image difference. This task performs such a decorrelation by
75  computing a decorrelation kernel (based upon the A&L matching
76  kernel and variances in the template and science images) and
77  convolving the image difference with it. This process is described
78  in detail in [DMTN-021](http://dmtn-021.lsst.io).
79 
80  This task has no standalone example, however it is applied as a
81  subtask of pipe.tasks.imageDifference.ImageDifferenceTask.
82  """
83  ConfigClass = DecorrelateALKernelConfig
84  _DefaultName = "ip_diffim_decorrelateALKernel"
85 
86  def __init__(self, *args, **kwargs):
87  """Create the image decorrelation Task
88 
89  Parameters
90  ----------
91  args :
92  arguments to be passed to ``lsst.pipe.base.task.Task.__init__``
93  kwargs :
94  keyword arguments to be passed to ``lsst.pipe.base.task.Task.__init__``
95  """
96  pipeBase.Task.__init__(self, *args, **kwargs)
97 
99  self.statsControl.setNumSigmaClip(3.)
100  self.statsControl.setNumIter(3)
101  self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes))
102 
103  def computeVarianceMean(self, exposure):
104  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
105  exposure.getMaskedImage().getMask(),
106  afwMath.MEANCLIP, self.statsControl)
107  var = statObj.getValue(afwMath.MEANCLIP)
108  return var
109 
110  @pipeBase.timeMethod
111  def run(self, exposure, templateExposure, subtractedExposure, psfMatchingKernel,
112  preConvKernel=None, xcen=None, ycen=None, svar=None, tvar=None):
113  """Perform decorrelation of an image difference exposure.
114 
115  Decorrelates the diffim due to the convolution of the templateExposure with the
116  A&L PSF matching kernel. Currently can accept a spatially varying matching kernel but in
117  this case it simply uses a static kernel from the center of the exposure. The decorrelation
118  is described in [DMTN-021, Equation 1](http://dmtn-021.lsst.io/#equation-1), where
119  `exposure` is I_1; templateExposure is I_2; `subtractedExposure` is D(k);
120  `psfMatchingKernel` is kappa; and svar and tvar are their respective
121  variances (see below).
122 
123  Parameters
124  ----------
125  exposure : `lsst.afw.image.Exposure`
126  The science afwImage.Exposure used for PSF matching
127  templateExposure : `lsst.afw.image.Exposure`
128  The template exposure used for PSF matching
129  subtractedExposure :
130  the subtracted exposure produced by
131  `ip_diffim.ImagePsfMatchTask.subtractExposures()`
132  psfMatchingKernel :
133  An (optionally spatially-varying) PSF matching kernel produced
134  by `ip_diffim.ImagePsfMatchTask.subtractExposures()`
135  preConvKernel :
136  if not None, then the `exposure` was pre-convolved with this kernel
137  xcen : `float`, optional
138  X-pixel coordinate to use for computing constant matching kernel to use
139  If `None` (default), then use the center of the image.
140  ycen : `float`, optional
141  Y-pixel coordinate to use for computing constant matching kernel to use
142  If `None` (default), then use the center of the image.
143  svar : `float`, optional
144  image variance for science image
145  If `None` (default) then compute the variance over the entire input science image.
146  tvar : `float`, optional
147  Image variance for template image
148  If `None` (default) then compute the variance over the entire input template image.
149 
150  Returns
151  -------
152  result : `Struct`
153  a `lsst.pipe.base.Struct` containing:
154 
155  - ``correctedExposure`` : the decorrelated diffim
156  - ``correctionKernel`` : the decorrelation correction kernel (which may be ignored)
157 
158  Notes
159  -----
160  The `subtractedExposure` is NOT updated
161 
162  The returned `correctedExposure` has an updated PSF as well.
163 
164  Here we currently convert a spatially-varying matching kernel into a constant kernel,
165  just by computing it at the center of the image (tickets DM-6243, DM-6244).
166 
167  We are also using a constant accross-the-image measure of sigma (sqrt(variance)) to compute
168  the decorrelation kernel.
169 
170  Still TBD (ticket DM-6580): understand whether the convolution is correctly modifying
171  the variance plane of the new subtractedExposure.
172  """
173  spatialKernel = psfMatchingKernel
174  kimg = afwImage.ImageD(spatialKernel.getDimensions())
175  bbox = subtractedExposure.getBBox()
176  if xcen is None:
177  xcen = (bbox.getBeginX() + bbox.getEndX()) / 2.
178  if ycen is None:
179  ycen = (bbox.getBeginY() + bbox.getEndY()) / 2.
180  self.log.info("Using matching kernel computed at (%d, %d)", xcen, ycen)
181  spatialKernel.computeImage(kimg, True, xcen, ycen)
182 
183  if svar is None:
184  svar = self.computeVarianceMean(exposure)
185  if tvar is None:
186  tvar = self.computeVarianceMean(templateExposure)
187  self.log.info("Variance (science, template): (%f, %f)", svar, tvar)
188 
189  # Should not happen unless entire image has been masked, which could happen
190  # if this is a small subimage of the main exposure. In this case, just return a full NaN
191  # exposure
192  if np.isnan(svar) or np.isnan(tvar):
193  # Double check that one of the exposures is all NaNs
194  if (np.all(np.isnan(exposure.getMaskedImage().getImage().getArray())) or
195  np.all(np.isnan(templateExposure.getMaskedImage().getImage().getArray()))):
196  self.log.warn('Template or science image is entirely NaNs: skipping decorrelation.')
197  outExposure = subtractedExposure.clone()
198  return pipeBase.Struct(correctedExposure=outExposure, correctionKernel=None)
199 
200  var = self.computeVarianceMean(subtractedExposure)
201  self.log.info("Variance (uncorrected diffim): %f", var)
202 
203  pck = None
204  if preConvKernel is not None:
205  self.log.info('Using a pre-convolution kernel as part of decorrelation.')
206  kimg2 = afwImage.ImageD(preConvKernel.getDimensions())
207  preConvKernel.computeImage(kimg2, False)
208  pck = kimg2.getArray()
209  corrKernel = DecorrelateALKernelTask._computeDecorrelationKernel(kimg.getArray(), svar, tvar,
210  pck)
211  correctedExposure, corrKern = DecorrelateALKernelTask._doConvolve(subtractedExposure, corrKernel)
212 
213  # Compute the subtracted exposure's updated psf
214  psf = subtractedExposure.getPsf().computeKernelImage(afwGeom.Point2D(xcen, ycen)).getArray()
215  psfc = DecorrelateALKernelTask.computeCorrectedDiffimPsf(corrKernel, psf, svar=svar, tvar=tvar)
216  psfcI = afwImage.ImageD(psfc.shape[0], psfc.shape[1])
217  psfcI.getArray()[:, :] = psfc
218  psfcK = afwMath.FixedKernel(psfcI)
219  psfNew = measAlg.KernelPsf(psfcK)
220  correctedExposure.setPsf(psfNew)
221 
222  var = self.computeVarianceMean(correctedExposure)
223  self.log.info("Variance (corrected diffim): %f", var)
224 
225  return pipeBase.Struct(correctedExposure=correctedExposure, correctionKernel=corrKern)
226 
227  @staticmethod
228  def _computeDecorrelationKernel(kappa, svar=0.04, tvar=0.04, preConvKernel=None):
229  """Compute the Lupton decorrelation post-conv. kernel for decorrelating an
230  image difference, based on the PSF-matching kernel.
231 
232  Parameters
233  ----------
234  kappa : `numpy.ndarray`
235  A matching kernel 2-d numpy.array derived from Alard & Lupton PSF matching
236  svar : `float`, optional
237  Average variance of science image used for PSF matching
238  tvar : `float`, optional
239  Average variance of template image used for PSF matching
240  preConvKernel If not None, then pre-filtering was applied
241  to science exposure, and this is the pre-convolution kernel.
242 
243  Returns
244  -------
245  fkernel : `numpy.ndarray`
246  a 2-d numpy.array containing the correction kernel
247 
248  Notes
249  -----
250  As currently implemented, kappa is a static (single, non-spatially-varying) kernel.
251  """
252  # Psf should not be <= 0, and messes up denominator; set the minimum value to MIN_KERNEL
253  MIN_KERNEL = 1.0e-4
254 
255  kappa = DecorrelateALKernelTask._fixOddKernel(kappa)
256  if preConvKernel is not None:
257  mk = DecorrelateALKernelTask._fixOddKernel(preConvKernel)
258  # Need to make them the same size
259  if kappa.shape[0] < mk.shape[0]:
260  diff = (mk.shape[0] - kappa.shape[0]) // 2
261  kappa = np.pad(kappa, (diff, diff), mode='constant')
262  elif kappa.shape[0] > mk.shape[0]:
263  diff = (kappa.shape[0] - mk.shape[0]) // 2
264  mk = np.pad(mk, (diff, diff), mode='constant')
265 
266  kft = np.fft.fft2(kappa)
267  kft2 = np.conj(kft) * kft
268  kft2[np.abs(kft2) < MIN_KERNEL] = MIN_KERNEL
269  denom = svar + tvar * kft2
270  if preConvKernel is not None:
271  mk = np.fft.fft2(mk)
272  mk2 = np.conj(mk) * mk
273  mk2[np.abs(mk2) < MIN_KERNEL] = MIN_KERNEL
274  denom = svar * mk2 + tvar * kft2
275  denom[np.abs(denom) < MIN_KERNEL] = MIN_KERNEL
276  kft = np.sqrt((svar + tvar) / denom)
277  pck = np.fft.ifft2(kft)
278  pck = np.fft.ifftshift(pck.real)
279  fkernel = DecorrelateALKernelTask._fixEvenKernel(pck)
280  if preConvKernel is not None:
281  # This is not pretty but seems to be necessary as the preConvKernel term seems to lead
282  # to a kernel that amplifies the noise way too much.
283  fkernel[fkernel > -np.min(fkernel)] = -np.min(fkernel)
284 
285  # I think we may need to "reverse" the PSF, as in the ZOGY (and Kaiser) papers...
286  # This is the same as taking the complex conjugate in Fourier space before FFT-ing back to real space.
287  if False: # TBD: figure this out. For now, we are turning it off.
288  fkernel = fkernel[::-1, :]
289 
290  return fkernel
291 
292  @staticmethod
293  def computeCorrectedDiffimPsf(kappa, psf, svar=0.04, tvar=0.04):
294  """Compute the (decorrelated) difference image's new PSF.
295  new_psf = psf(k) * sqrt((svar + tvar) / (svar + tvar * kappa_ft(k)**2))
296 
297  Parameters
298  ----------
299  kappa : `numpy.ndarray`
300  A matching kernel array derived from Alard & Lupton PSF matching
301  psf : `numpy.ndarray`
302  The uncorrected psf array of the science image (and also of the diffim)
303  svar : `float`, optional
304  Average variance of science image used for PSF matching
305  tvar : `float`, optional
306  Average variance of template image used for PSF matching
307 
308  Returns
309  -------
310  pcf : `numpy.ndarray`
311  a 2-d numpy.array containing the new PSF
312  """
313  def post_conv_psf_ft2(psf, kernel, svar, tvar):
314  # Pad psf or kernel symmetrically to make them the same size!
315  # Note this assumes they are both square (width == height)
316  if psf.shape[0] < kernel.shape[0]:
317  diff = (kernel.shape[0] - psf.shape[0]) // 2
318  psf = np.pad(psf, (diff, diff), mode='constant')
319  elif psf.shape[0] > kernel.shape[0]:
320  diff = (psf.shape[0] - kernel.shape[0]) // 2
321  kernel = np.pad(kernel, (diff, diff), mode='constant')
322  psf_ft = np.fft.fft2(psf)
323  kft = np.fft.fft2(kernel)
324  out = psf_ft * np.sqrt((svar + tvar) / (svar + tvar * kft**2))
325  return out
326 
327  def post_conv_psf(psf, kernel, svar, tvar):
328  kft = post_conv_psf_ft2(psf, kernel, svar, tvar)
329  out = np.fft.ifft2(kft)
330  return out
331 
332  pcf = post_conv_psf(psf=psf, kernel=kappa, svar=svar, tvar=tvar)
333  pcf = pcf.real / pcf.real.sum()
334  return pcf
335 
336  @staticmethod
337  def _fixOddKernel(kernel):
338  """Take a kernel with odd dimensions and make them even for FFT
339 
340  Parameters
341  ----------
342  kernel : `numpy.array`
343  a numpy.array
344 
345  Returns
346  -------
347  out : `numpy.array`
348  a fixed kernel numpy.array. Returns a copy if the dimensions needed to change;
349  otherwise just return the input kernel.
350  """
351  # Note this works best for the FFT if we left-pad
352  out = kernel
353  changed = False
354  if (out.shape[0] % 2) == 1:
355  out = np.pad(out, ((1, 0), (0, 0)), mode='constant')
356  changed = True
357  if (out.shape[1] % 2) == 1:
358  out = np.pad(out, ((0, 0), (1, 0)), mode='constant')
359  changed = True
360  if changed:
361  out *= (np.mean(kernel) / np.mean(out)) # need to re-scale to same mean for FFT
362  return out
363 
364  @staticmethod
365  def _fixEvenKernel(kernel):
366  """Take a kernel with even dimensions and make them odd, centered correctly.
367 
368  Parameters
369  ----------
370  kernel : `numpy.array`
371  a numpy.array
372 
373  Returns
374  -------
375  out : `numpy.array`
376  a fixed kernel numpy.array
377  """
378  # Make sure the peak (close to a delta-function) is in the center!
379  maxloc = np.unravel_index(np.argmax(kernel), kernel.shape)
380  out = np.roll(kernel, kernel.shape[0]//2 - maxloc[0], axis=0)
381  out = np.roll(out, out.shape[1]//2 - maxloc[1], axis=1)
382  # Make sure it is odd-dimensioned by trimming it.
383  if (out.shape[0] % 2) == 0:
384  maxloc = np.unravel_index(np.argmax(out), out.shape)
385  if out.shape[0] - maxloc[0] > maxloc[0]:
386  out = out[:-1, :]
387  else:
388  out = out[1:, :]
389  if out.shape[1] - maxloc[1] > maxloc[1]:
390  out = out[:, :-1]
391  else:
392  out = out[:, 1:]
393  return out
394 
395  @staticmethod
396  def _doConvolve(exposure, kernel):
397  """Convolve an Exposure with a decorrelation convolution kernel.
398 
399  Parameters
400  ----------
401  exposure : `lsst.afw.image.Exposure`
402  Input exposure to be convolved.
403  kernel : `numpy.array`
404  Input 2-d numpy.array to convolve the image with
405 
406  Returns
407  -------
408  out : `lsst.afw.image.Exposure`
409  a new Exposure with the convolved pixels and the (possibly
410  re-centered) kernel.
411 
412  Notes
413  -----
414  We re-center the kernel if necessary and return the possibly re-centered kernel
415  """
416  kernelImg = afwImage.ImageD(kernel.shape[0], kernel.shape[1])
417  kernelImg.getArray()[:, :] = kernel
418  kern = afwMath.FixedKernel(kernelImg)
419  maxloc = np.unravel_index(np.argmax(kernel), kernel.shape)
420  kern.setCtrX(maxloc[0])
421  kern.setCtrY(maxloc[1])
422  outExp = exposure.clone() # Do this to keep WCS, PSF, masks, etc.
423  convCntrl = afwMath.ConvolutionControl(False, True, 0)
424  afwMath.convolve(outExp.getMaskedImage(), exposure.getMaskedImage(), kern, convCntrl)
425 
426  return outExp, kern
427 
428 
430  """Task to be used as an ImageMapper for performing
431  A&L decorrelation on subimages on a grid across a A&L difference image.
432 
433  This task subclasses DecorrelateALKernelTask in order to implement
434  all of that task's configuration parameters, as well as its `run` method.
435  """
436 
437  ConfigClass = DecorrelateALKernelConfig
438  _DefaultName = 'ip_diffim_decorrelateALKernelMapper'
439 
440  def __init__(self, *args, **kwargs):
441  DecorrelateALKernelTask.__init__(self, *args, **kwargs)
442 
443  def run(self, subExposure, expandedSubExposure, fullBBox,
444  template, science, alTaskResult=None, psfMatchingKernel=None,
445  preConvKernel=None, **kwargs):
446  """Perform decorrelation operation on `subExposure`, using
447  `expandedSubExposure` to allow for invalid edge pixels arising from
448  convolutions.
449 
450  This method performs A&L decorrelation on `subExposure` using
451  local measures for image variances and PSF. `subExposure` is a
452  sub-exposure of the non-decorrelated A&L diffim. It also
453  requires the corresponding sub-exposures of the template
454  (`template`) and science (`science`) exposures.
455 
456  Parameters
457  ----------
458  subExposure : `lsst.afw.image.Exposure`
459  the sub-exposure of the diffim
460  expandedSubExposure : `lsst.afw.image.Exposure`
461  the expanded sub-exposure upon which to operate
462  fullBBox : `lsst.afw.geom.BoundingBox`
463  the bounding box of the original exposure
464  template : `lsst.afw.image.Exposure`
465  the corresponding sub-exposure of the template exposure
466  science : `lsst.afw.image.Exposure`
467  the corresponding sub-exposure of the science exposure
468  alTaskResult : `lsst.pipe.base.Struct`
469  the result of A&L image differencing on `science` and
470  `template`, importantly containing the resulting
471  `psfMatchingKernel`. Can be `None`, only if
472  `psfMatchingKernel` is not `None`.
473  psfMatchingKernel : Alternative parameter for passing the
474  A&L `psfMatchingKernel` directly.
475  preConvKernel : If not None, then pre-filtering was applied
476  to science exposure, and this is the pre-convolution
477  kernel.
478  kwargs :
479  additional keyword arguments propagated from
480  `ImageMapReduceTask.run`.
481 
482  Returns
483  -------
484  A `pipeBase.Struct` containing:
485 
486  - ``subExposure`` : the result of the `subExposure` processing.
487  - ``decorrelationKernel`` : the decorrelation kernel, currently
488  not used.
489 
490  Notes
491  -----
492  This `run` method accepts parameters identical to those of
493  `ImageMapper.run`, since it is called from the
494  `ImageMapperTask`. See that class for more information.
495  """
496  templateExposure = template # input template
497  scienceExposure = science # input science image
498  if alTaskResult is None and psfMatchingKernel is None:
499  raise RuntimeError('Both alTaskResult and psfMatchingKernel cannot be None')
500  psfMatchingKernel = alTaskResult.psfMatchingKernel if alTaskResult is not None else psfMatchingKernel
501 
502  # subExp and expandedSubExp are subimages of the (un-decorrelated) diffim!
503  # So here we compute corresponding subimages of templateExposure and scienceExposure
504  subExp2 = scienceExposure.Factory(scienceExposure, expandedSubExposure.getBBox())
505  subExp1 = templateExposure.Factory(templateExposure, expandedSubExposure.getBBox())
506 
507  # Prevent too much log INFO verbosity from DecorrelateALKernelTask.run
508  logLevel = self.log.getLevel()
509  self.log.setLevel(lsst.log.WARN)
510  res = DecorrelateALKernelTask.run(self, subExp2, subExp1, expandedSubExposure,
511  psfMatchingKernel, preConvKernel)
512  self.log.setLevel(logLevel) # reset the log level
513 
514  diffim = res.correctedExposure.Factory(res.correctedExposure, subExposure.getBBox())
515  out = pipeBase.Struct(subExposure=diffim, decorrelationKernel=res.correctionKernel)
516  return out
517 
518 
520  """Configuration parameters for the ImageMapReduceTask to direct it to use
521  DecorrelateALKernelMapper as its mapper for A&L decorrelation.
522  """
523  mapper = pexConfig.ConfigurableField(
524  doc='A&L decorrelation task to run on each sub-image',
525  target=DecorrelateALKernelMapper
526  )
527 
528 
529 class DecorrelateALKernelSpatialConfig(pexConfig.Config):
530  """Configuration parameters for the DecorrelateALKernelSpatialTask.
531  """
532  decorrelateConfig = pexConfig.ConfigField(
533  dtype=DecorrelateALKernelConfig,
534  doc='DecorrelateALKernel config to use when running on complete exposure (non spatially-varying)',
535  )
536 
537  decorrelateMapReduceConfig = pexConfig.ConfigField(
538  dtype=DecorrelateALKernelMapReduceConfig,
539  doc='DecorrelateALKernelMapReduce config to use when running on each sub-image (spatially-varying)',
540  )
541 
542  ignoreMaskPlanes = pexConfig.ListField(
543  dtype=str,
544  doc="""Mask planes to ignore for sigma-clipped statistics""",
545  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
546  )
547 
548  def setDefaults(self):
549  self.decorrelateMapReduceConfig.gridStepX = self.decorrelateMapReduceConfig.gridStepY = 40
550  self.decorrelateMapReduceConfig.cellSizeX = self.decorrelateMapReduceConfig.cellSizeY = 41
551  self.decorrelateMapReduceConfig.borderSizeX = self.decorrelateMapReduceConfig.borderSizeY = 8
552  self.decorrelateMapReduceConfig.reducer.reduceOperation = 'average'
553 
554 
555 class DecorrelateALKernelSpatialTask(pipeBase.Task):
556  """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference
557 
558  Notes
559  -----
560 
561  Pipe-task that removes the neighboring-pixel covariance in an
562  image difference that are added when the template image is
563  convolved with the Alard-Lupton PSF matching kernel.
564 
565  This task is a simple wrapper around @ref DecorrelateALKernelTask,
566  which takes a `spatiallyVarying` parameter in its `run` method. If
567  it is `False`, then it simply calls the `run` method of @ref
568  DecorrelateALKernelTask. If it is True, then it uses the @ref
569  ImageMapReduceTask framework to break the exposures into
570  subExposures on a grid, and performs the `run` method of @ref
571  DecorrelateALKernelTask on each subExposure. This enables it to
572  account for spatially-varying PSFs and noise in the exposures when
573  performing the decorrelation.
574 
575  This task has no standalone example, however it is applied as a
576  subtask of pipe.tasks.imageDifference.ImageDifferenceTask.
577  There is also an example of its use in `tests/testImageDecorrelation.py`.
578  """
579  ConfigClass = DecorrelateALKernelSpatialConfig
580  _DefaultName = "ip_diffim_decorrelateALKernelSpatial"
581 
582  def __init__(self, *args, **kwargs):
583  """Create the image decorrelation Task
584 
585  Parameters
586  ----------
587  args :
588  arguments to be passed to
589  `lsst.pipe.base.task.Task.__init__`
590  kwargs :
591  additional keyword arguments to be passed to
592  `lsst.pipe.base.task.Task.__init__`
593  """
594  pipeBase.Task.__init__(self, *args, **kwargs)
595 
597  self.statsControl.setNumSigmaClip(3.)
598  self.statsControl.setNumIter(3)
599  self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes))
600 
601  def computeVarianceMean(self, exposure):
602  """Compute the mean of the variance plane of `exposure`.
603  """
604  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
605  exposure.getMaskedImage().getMask(),
606  afwMath.MEANCLIP, self.statsControl)
607  var = statObj.getValue(afwMath.MEANCLIP)
608  return var
609 
610  def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel,
611  spatiallyVarying=True, preConvKernel=None):
612  """Perform decorrelation of an image difference exposure.
613 
614  Decorrelates the diffim due to the convolution of the
615  templateExposure with the A&L psfMatchingKernel. If
616  `spatiallyVarying` is True, it utilizes the spatially varying
617  matching kernel via the `imageMapReduce` framework to perform
618  spatially-varying decorrelation on a grid of subExposures.
619 
620  Parameters
621  ----------
622  scienceExposure : `lsst.afw.image.Exposure`
623  the science Exposure used for PSF matching
624  templateExposure : `lsst.afw.image.Exposure`
625  the template Exposure used for PSF matching
626  subtractedExposure : `lsst.afw.image.Exposure`
627  the subtracted Exposure produced by `ip_diffim.ImagePsfMatchTask.subtractExposures()`
628  psfMatchingKernel :
629  an (optionally spatially-varying) PSF matching kernel produced
630  by `ip_diffim.ImagePsfMatchTask.subtractExposures()`
631  spatiallyVarying : `bool`
632  if True, perform the spatially-varying operation
633  preConvKernel : `lsst.meas.algorithms.Psf`
634  if not none, the scienceExposure has been pre-filtered with this kernel. (Currently
635  this option is experimental.)
636 
637  Returns
638  -------
639  results : `lsst.pipe.base.Struct`
640  a structure containing:
641 
642  - ``correctedExposure`` : the decorrelated diffim
643 
644  """
645  self.log.info('Running A&L decorrelation: spatiallyVarying=%r' % spatiallyVarying)
646 
647  svar = self.computeVarianceMean(scienceExposure)
648  tvar = self.computeVarianceMean(templateExposure)
649  if np.isnan(svar) or np.isnan(tvar): # Should not happen unless entire image has been masked.
650  # Double check that one of the exposures is all NaNs
651  if (np.all(np.isnan(scienceExposure.getMaskedImage().getImage().getArray())) or
652  np.all(np.isnan(templateExposure.getMaskedImage().getImage().getArray()))):
653  self.log.warn('Template or science image is entirely NaNs: skipping decorrelation.')
654  if np.isnan(svar):
655  svar = 1e-9
656  if np.isnan(tvar):
657  tvar = 1e-9
658 
659  var = self.computeVarianceMean(subtractedExposure)
660 
661  if spatiallyVarying:
662  self.log.info("Variance (science, template): (%f, %f)", svar, tvar)
663  self.log.info("Variance (uncorrected diffim): %f", var)
664  config = self.config.decorrelateMapReduceConfig
665  task = ImageMapReduceTask(config=config)
666  results = task.run(subtractedExposure, science=scienceExposure,
667  template=templateExposure, psfMatchingKernel=psfMatchingKernel,
668  preConvKernel=preConvKernel, forceEvenSized=True)
669  results.correctedExposure = results.exposure
670 
671  # Make sure masks of input image are propagated to diffim
672  def gm(exp):
673  return exp.getMaskedImage().getMask()
674  gm(results.correctedExposure)[:, :] = gm(subtractedExposure)
675 
676  var = self.computeVarianceMean(results.correctedExposure)
677  self.log.info("Variance (corrected diffim): %f", var)
678 
679  else:
680  config = self.config.decorrelateConfig
681  task = DecorrelateALKernelTask(config=config)
682  results = task.run(scienceExposure, templateExposure,
683  subtractedExposure, psfMatchingKernel, preConvKernel=preConvKernel)
684 
685  return results
def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel, spatiallyVarying=True, preConvKernel=None)
def run(self, exposure, templateExposure, subtractedExposure, psfMatchingKernel, preConvKernel=None, xcen=None, ycen=None, svar=None, tvar=None)
Statistics makeStatistics(lsst::afw::math::MaskedVector< EntryT > const &mv, std::vector< WeightPixel > const &vweights, int const flags, StatisticsControl const &sctrl=StatisticsControl())
def computeCorrectedDiffimPsf(kappa, psf, svar=0.04, tvar=0.04)
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, bool doNormalize, bool doCopyEdge=false)
def run(self, subExposure, expandedSubExposure, fullBBox, template, science, alTaskResult=None, psfMatchingKernel=None, preConvKernel=None, kwargs)