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