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