lsst.ip.diffim  22.0.1-17-ge5a99e88+50b11f7f3e
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.math as afwMath
27 import lsst.geom as geom
28 import lsst.log
29 import lsst.meas.algorithms as measAlg
30 import lsst.pex.config as pexConfig
31 import lsst.pipe.base as pipeBase
32 
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  """Configuration parameters for the DecorrelateALKernelTask
44  """
45 
46  ignoreMaskPlanes = pexConfig.ListField(
47  dtype=str,
48  doc="""Mask planes to ignore for sigma-clipped statistics""",
49  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
50  )
51 
52 
53 class DecorrelateALKernelTask(pipeBase.Task):
54  """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference
55 
56  Notes
57  -----
58 
59  Pipe-task that removes the neighboring-pixel covariance in an
60  image difference that are added when the template image is
61  convolved with the Alard-Lupton PSF matching kernel.
62 
63  The image differencing pipeline task @link
64  ip.diffim.psfMatch.PsfMatchTask PSFMatchTask@endlink and @link
65  ip.diffim.psfMatch.PsfMatchConfigAL PSFMatchConfigAL@endlink uses
66  the Alard and Lupton (1998) method for matching the PSFs of the
67  template and science exposures prior to subtraction. The
68  Alard-Lupton method identifies a matching kernel, which is then
69  (typically) convolved with the template image to perform PSF
70  matching. This convolution has the effect of adding covariance
71  between neighboring pixels in the template image, which is then
72  added to the image difference by subtraction.
73 
74  The pixel covariance may be corrected by whitening the noise of
75  the image difference. This task performs such a decorrelation by
76  computing a decorrelation kernel (based upon the A&L matching
77  kernel and variances in the template and science images) and
78  convolving the image difference with it. This process is described
79  in detail in [DMTN-021](http://dmtn-021.lsst.io).
80 
81  This task has no standalone example, however it is applied as a
82  subtask of pipe.tasks.imageDifference.ImageDifferenceTask.
83  """
84  ConfigClass = DecorrelateALKernelConfig
85  _DefaultName = "ip_diffim_decorrelateALKernel"
86 
87  def __init__(self, *args, **kwargs):
88  """Create the image decorrelation Task
89 
90  Parameters
91  ----------
92  args :
93  arguments to be passed to ``lsst.pipe.base.task.Task.__init__``
94  kwargs :
95  keyword arguments to be passed to ``lsst.pipe.base.task.Task.__init__``
96  """
97  pipeBase.Task.__init__(self, *args, **kwargs)
98 
100  self.statsControlstatsControl.setNumSigmaClip(3.)
101  self.statsControlstatsControl.setNumIter(3)
102  self.statsControlstatsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes))
103 
104  def computeVarianceMean(self, exposure):
105  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
106  exposure.getMaskedImage().getMask(),
107  afwMath.MEANCLIP, self.statsControlstatsControl)
108  var = statObj.getValue(afwMath.MEANCLIP)
109  return var
110 
111  @pipeBase.timeMethod
112  def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel,
113  preConvKernel=None, xcen=None, ycen=None, svar=None, tvar=None, templateMatched=True):
114  """Perform decorrelation of an image difference exposure.
115 
116  Decorrelates the diffim due to the convolution of the templateExposure with the
117  A&L PSF matching kernel. Currently can accept a spatially varying matching kernel but in
118  this case it simply uses a static kernel from the center of the exposure. The decorrelation
119  is described in [DMTN-021, Equation 1](http://dmtn-021.lsst.io/#equation-1), where
120  `exposure` is I_1; templateExposure is I_2; `subtractedExposure` is D(k);
121  `psfMatchingKernel` is kappa; and svar and tvar are their respective
122  variances (see below).
123 
124  Parameters
125  ----------
126  scienceExposure : `lsst.afw.image.Exposure`
127  The original science exposure (before `preConvKernel` applied).
128  templateExposure : `lsst.afw.image.Exposure`
129  The original template exposure warped into the science exposure dimensions.
130  subtractedExposure : `lsst.afw.iamge.Exposure`
131  the subtracted exposure produced by
132  `ip_diffim.ImagePsfMatchTask.subtractExposures()`. The `subtractedExposure` must
133  inherit its PSF from `exposure`, see notes below.
134  psfMatchingKernel : `lsst.afw.detection.Psf`
135  An (optionally spatially-varying) PSF matching kernel produced
136  by `ip_diffim.ImagePsfMatchTask.subtractExposures()`.
137  preConvKernel : `lsst.afw.math.Kernel`, optional
138  if not None, then the `scienceExposure` was pre-convolved with this kernel.
139  Allowed only if ``templateMatched==True``.
140  xcen : `float`, optional
141  X-pixel coordinate to use for computing constant matching kernel to use
142  If `None` (default), then use the center of the image.
143  ycen : `float`, optional
144  Y-pixel coordinate to use for computing constant matching kernel to use
145  If `None` (default), then use the center of the image.
146  svar : `float`, optional
147  Image variance for science image
148  If `None` (default) then compute the variance over the entire input science image.
149  tvar : `float`, optional
150  Image variance for template image
151  If `None` (default) then compute the variance over the entire input template image.
152  templateMatched : `bool`, optional
153  If True, the template exposure was matched (convolved) to the science exposure.
154  See also notes below.
155 
156  Returns
157  -------
158  result : `lsst.pipe.base.Struct`
159  - ``correctedExposure`` : the decorrelated diffim
160 
161  Notes
162  -----
163  The `subtractedExposure` is NOT updated. The returned `correctedExposure` has an updated but
164  spatially fixed PSF. It is calculated as the center of image PSF corrected by the center of
165  image matching kernel.
166 
167  If ``templateMatched==True``, the templateExposure was matched (convolved)
168  to the ``scienceExposure`` by ``psfMatchingKernel``. Otherwise the ``scienceExposure``
169  was matched (convolved) by ``psfMatchingKernel``.
170 
171  This task discards the variance plane of ``subtractedExposure`` and re-computes
172  it from the variance planes of ``scienceExposure`` and ``templateExposure``.
173  The image plane of ``subtractedExposure`` must be at the photometric level
174  set by the AL PSF matching in `ImagePsfMatchTask.subtractExposures`.
175  The assumptions about the photometric level are controlled by the
176  `templateMatched` option in this task.
177 
178  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 
181  We are also using a constant accross-the-image measure of sigma (sqrt(variance)) to compute
182  the decorrelation kernel.
183 
184  TODO DM-23857 As part of the spatially varying correction implementation
185  consider whether returning a Struct is still necessary.
186  """
187  if preConvKernel is not None and not templateMatched:
188  raise ValueError("Pre-convolution and the matching of the "
189  "science exposure is not supported.")
190 
191  spatialKernel = psfMatchingKernel
192  kimg = afwImage.ImageD(spatialKernel.getDimensions())
193  bbox = subtractedExposure.getBBox()
194  if xcen is None:
195  xcen = (bbox.getBeginX() + bbox.getEndX()) / 2.
196  if ycen is None:
197  ycen = (bbox.getBeginY() + bbox.getEndY()) / 2.
198  self.log.info("Using matching kernel computed at (%d, %d)", xcen, ycen)
199  spatialKernel.computeImage(kimg, False, xcen, ycen)
200 
201  if svar is None:
202  svar = self.computeVarianceMeancomputeVarianceMean(scienceExposure)
203  if tvar is None:
204  tvar = self.computeVarianceMeancomputeVarianceMean(templateExposure)
205  self.log.info("Original variance plane means. Science:%.5e, warped template:%.5e)",
206  svar, tvar)
207 
208  if templateMatched:
209  # Regular subtraction, we convolved the template
210  self.log.info("Decorrelation after template image convolution")
211  expVar = svar
212  matchedVar = tvar
213  exposure = scienceExposure
214  matchedExposure = templateExposure
215  else:
216  # We convolved the science image
217  self.log.info("Decorrelation after science image convolution")
218  expVar = tvar
219  matchedVar = svar
220  exposure = templateExposure
221  matchedExposure = scienceExposure
222 
223  # Should not happen unless entire image has been masked, which could happen
224  # if this is a small subimage of the main exposure. In this case, just return a full NaN
225  # exposure
226  if np.isnan(expVar) or np.isnan(matchedVar):
227  # Double check that one of the exposures is all NaNs
228  if (np.all(np.isnan(exposure.image.array))
229  or np.all(np.isnan(matchedExposure.image.array))):
230  self.log.warning('Template or science image is entirely NaNs: skipping decorrelation.')
231  outExposure = subtractedExposure.clone()
232  return pipeBase.Struct(correctedExposure=outExposure, )
233 
234  # The maximal correction value converges to sqrt(matchedVar/expVar).
235  # Correction divergence warning if the correction exceeds 4 orders of magnitude.
236  mOverExpVar = matchedVar/expVar
237  if mOverExpVar > 1e8:
238  self.log.warning("Diverging correction: matched image variance is "
239  " much larger than the unconvolved one's"
240  ", matchedVar/expVar:%.2e", mOverExpVar)
241 
242  oldVarMean = self.computeVarianceMeancomputeVarianceMean(subtractedExposure)
243  self.log.info("Variance plane mean of uncorrected diffim: %f", oldVarMean)
244 
245  if preConvKernel is not None:
246  self.log.info('Using a pre-convolution kernel as part of decorrelation correction.')
247  kimg2 = afwImage.ImageD(preConvKernel.getDimensions())
248  preConvKernel.computeImage(kimg2, False)
249  pckArr = kimg2.array
250 
251  kArr = kimg.array
252  diffExpArr = subtractedExposure.image.array
253  psfImg = subtractedExposure.getPsf().computeKernelImage(geom.Point2D(xcen, ycen))
254  psfDim = psfImg.getDimensions()
255  psfArr = psfImg.array
256 
257  # Determine the common shape
258  kSum = np.sum(kArr)
259  kSumSq = kSum*kSum
260  self.log.debug("Matching kernel sum: %.3e", kSum)
261  preSum = 1.
262  if preConvKernel is None:
263  self.computeCommonShapecomputeCommonShape(kArr.shape, psfArr.shape, diffExpArr.shape)
264  corrft = self.computeCorrectioncomputeCorrection(kArr, expVar, matchedVar)
265  else:
266  preSum = np.sum(pckArr)
267  self.log.debug("Pre-convolution kernel sum: %.3e", preSum)
268  self.computeCommonShapecomputeCommonShape(pckArr.shape, kArr.shape,
269  psfArr.shape, diffExpArr.shape)
270  corrft = self.computeCorrectioncomputeCorrection(kArr, expVar, matchedVar, preConvArr=pckArr)
271 
272  diffExpArr = self.computeCorrectedImagecomputeCorrectedImage(corrft, diffExpArr)
273  corrPsfArr = self.computeCorrectedDiffimPsfcomputeCorrectedDiffimPsf(corrft, psfArr)
274 
275  psfcI = afwImage.ImageD(psfDim)
276  psfcI.array = corrPsfArr
277  psfcK = afwMath.FixedKernel(psfcI)
278  psfNew = measAlg.KernelPsf(psfcK)
279 
280  correctedExposure = subtractedExposure.clone()
281  correctedExposure.image.array[...] = diffExpArr # Allow for numpy type casting
282  # The subtracted exposure variance plane is already correlated, we cannot propagate
283  # it through another convolution; instead we need to use the uncorrelated originals
284  # The whitening should scale it to expVar + matchedVar on average
285  varImg = correctedExposure.variance.array
286  # Allow for numpy type casting
287  varImg[...] = preSum*preSum*exposure.variance.array + kSumSq*matchedExposure.variance.array
288  if not templateMatched:
289  # ImagePsfMatch.subtractExposures re-scales the difference in
290  # the science image convolution mode
291  varImg /= kSumSq
292  correctedExposure.setPsf(psfNew)
293 
294  newVarMean = self.computeVarianceMeancomputeVarianceMean(correctedExposure)
295  self.log.info("Variance plane mean of corrected diffim: %.5e", newVarMean)
296 
297  # TODO DM-23857 As part of the spatially varying correction implementation
298  # consider whether returning a Struct is still necessary.
299  return pipeBase.Struct(correctedExposure=correctedExposure, )
300 
301  def computeCommonShape(self, *shapes):
302  """Calculate the common shape for FFT operations. Set `self.freqSpaceShape`
303  internally.
304 
305  Parameters
306  ----------
307  shapes : one or more `tuple` of `int`
308  Shapes of the arrays. All must have the same dimensionality.
309  At least one shape must be provided.
310 
311  Returns
312  -------
313  None.
314 
315  Notes
316  -----
317  For each dimension, gets the smallest even number greater than or equal to
318  `N1+N2-1` where `N1` and `N2` are the two largest values.
319  In case of only one shape given, rounds up to even each dimension value.
320  """
321  S = np.array(shapes, dtype=int)
322  if len(shapes) > 2:
323  S.sort(axis=0)
324  S = S[-2:]
325  if len(shapes) > 1:
326  commonShape = np.sum(S, axis=0) - 1
327  else:
328  commonShape = S[0]
329  commonShape[commonShape % 2 != 0] += 1
330  self.freqSpaceShapefreqSpaceShape = tuple(commonShape)
331  self.log.info("Common frequency space shape %s", self.freqSpaceShapefreqSpaceShape)
332 
333  @staticmethod
334  def padCenterOriginArray(A, newShape: tuple, useInverse=False):
335  """Zero pad an image where the origin is at the center and replace the
336  origin to the corner as required by the periodic input of FFT. Implement also
337  the inverse operation, crop the padding and re-center data.
338 
339  Parameters
340  ----------
341  A : `numpy.ndarray`
342  An array to copy from.
343  newShape : `tuple` of `int`
344  The dimensions of the resulting array. For padding, the resulting array
345  must be larger than A in each dimension. For the inverse operation this
346  must be the original, before padding size of the array.
347  useInverse : bool, optional
348  Selector of forward, add padding, operation (False)
349  or its inverse, crop padding, operation (True).
350 
351  Returns
352  -------
353  R : `numpy.ndarray`
354  The padded or unpadded array with shape of `newShape` and the same dtype as A.
355 
356  Notes
357  -----
358  For odd dimensions, the splitting is rounded to
359  put the center pixel into the new corner origin (0,0). This is to be consistent
360  e.g. for a dirac delta kernel that is originally located at the center pixel.
361  """
362 
363  # The forward and inverse operations should round odd dimension halves at the opposite
364  # sides to get the pixels back to their original positions.
365  if not useInverse:
366  # Forward operation: First and second halves with respect to the axes of A.
367  firstHalves = [x//2 for x in A.shape]
368  secondHalves = [x-y for x, y in zip(A.shape, firstHalves)]
369  else:
370  # Inverse operation: Opposite rounding
371  secondHalves = [x//2 for x in newShape]
372  firstHalves = [x-y for x, y in zip(newShape, secondHalves)]
373 
374  R = np.zeros_like(A, shape=newShape)
375  R[-firstHalves[0]:, -firstHalves[1]:] = A[:firstHalves[0], :firstHalves[1]]
376  R[:secondHalves[0], -firstHalves[1]:] = A[-secondHalves[0]:, :firstHalves[1]]
377  R[:secondHalves[0], :secondHalves[1]] = A[-secondHalves[0]:, -secondHalves[1]:]
378  R[-firstHalves[0]:, :secondHalves[1]] = A[:firstHalves[0], -secondHalves[1]:]
379  return R
380 
381  def computeCorrection(self, kappa, svar, tvar, preConvArr=None):
382  """Compute the Lupton decorrelation post-convolution kernel for decorrelating an
383  image difference, based on the PSF-matching kernel.
384 
385  Parameters
386  ----------
387  kappa : `numpy.ndarray`
388  A matching kernel 2-d numpy.array derived from Alard & Lupton PSF matching.
389  svar : `float`
390  Average variance of science image used for PSF matching.
391  tvar : `float`
392  Average variance of the template (matched) image used for PSF matching.
393  preConvArr : `numpy.ndarray`, optional
394  If not None, then pre-filtering was applied
395  to science exposure, and this is the pre-convolution kernel.
396 
397  Returns
398  -------
399  corrft : `numpy.ndarray` of `float`
400  The frequency space representation of the correction. The array is real (dtype float).
401  Shape is `self.freqSpaceShape`.
402 
403  Notes
404  -----
405  The maximum correction factor converges to `sqrt(tvar/svar)` towards high frequencies.
406  This should be a plausible value.
407  """
408  kSum = np.sum(kappa)
409  kappa = self.padCenterOriginArraypadCenterOriginArray(kappa, self.freqSpaceShapefreqSpaceShape)
410  kft = np.fft.fft2(kappa)
411  kftAbsSq = np.real(np.conj(kft) * kft)
412  # If there is no pre-convolution kernel, use placeholder scalars
413  if preConvArr is None:
414  preSum = 1.
415  preAbsSq = 1.
416  else:
417  preSum = np.sum(preConvArr)
418  preConvArr = self.padCenterOriginArraypadCenterOriginArray(preConvArr, self.freqSpaceShapefreqSpaceShape)
419  preK = np.fft.fft2(preConvArr)
420  preAbsSq = np.real(np.conj(preK)*preK)
421 
422  denom = svar * preAbsSq + tvar * kftAbsSq
423  # Division by zero protection, though we don't expect to hit it
424  # (rather we'll have numerical noise)
425  tiny = np.finfo(kftAbsSq.dtype).tiny * 1000.
426  flt = denom < tiny
427  sumFlt = np.sum(flt)
428  if sumFlt > 0:
429  self.log.warning("Avoid zero division. Skip decorrelation "
430  "at %f divergent frequencies.", sumFlt)
431  denom[flt] = 1.
432  kft = np.sqrt((svar * preSum*preSum + tvar * kSum*kSum) / denom)
433  # Don't do any correction at these frequencies
434  # the difference image should be close to zero anyway, so can't be decorrelated
435  if sumFlt > 0:
436  kft[flt] = 1.
437  return kft
438 
439  def computeCorrectedDiffimPsf(self, corrft, psfOld):
440  """Compute the (decorrelated) difference image's new PSF.
441 
442  Parameters
443  ----------
444  corrft : `numpy.ndarray`
445  The frequency space representation of the correction calculated by
446  `computeCorrection`. Shape must be `self.freqSpaceShape`.
447  psfOld : `numpy.ndarray`
448  The psf of the difference image to be corrected.
449 
450  Returns
451  -------
452  psfNew : `numpy.ndarray`
453  The corrected psf, same shape as `psfOld`, sum normed to 1.
454 
455  Notes
456  -----
457  There is no algorithmic guarantee that the corrected psf can
458  meaningfully fit to the same size as the original one.
459  """
460  psfShape = psfOld.shape
461  psfNew = self.padCenterOriginArraypadCenterOriginArray(psfOld, self.freqSpaceShapefreqSpaceShape)
462  psfNew = np.fft.fft2(psfNew)
463  psfNew *= corrft
464  psfNew = np.fft.ifft2(psfNew)
465  psfNew = psfNew.real
466  psfNew = self.padCenterOriginArraypadCenterOriginArray(psfNew, psfShape, useInverse=True)
467  psfNew = psfNew/psfNew.sum()
468  return psfNew
469 
470  def computeCorrectedImage(self, corrft, imgOld):
471  """Compute the decorrelated difference image.
472 
473  Parameters
474  ----------
475  corrft : `numpy.ndarray`
476  The frequency space representation of the correction calculated by
477  `computeCorrection`. Shape must be `self.freqSpaceShape`.
478  imgOld : `numpy.ndarray`
479  The difference image to be corrected.
480 
481  Returns
482  -------
483  imgNew : `numpy.ndarray`
484  The corrected image, same size as the input.
485  """
486  expShape = imgOld.shape
487  imgNew = np.copy(imgOld)
488  filtInf = np.isinf(imgNew)
489  filtNan = np.isnan(imgNew)
490  imgNew[filtInf] = np.nan
491  imgNew[filtInf | filtNan] = np.nanmean(imgNew)
492  imgNew = self.padCenterOriginArraypadCenterOriginArray(imgNew, self.freqSpaceShapefreqSpaceShape)
493  imgNew = np.fft.fft2(imgNew)
494  imgNew *= corrft
495  imgNew = np.fft.ifft2(imgNew)
496  imgNew = imgNew.real
497  imgNew = self.padCenterOriginArraypadCenterOriginArray(imgNew, expShape, useInverse=True)
498  imgNew[filtNan] = np.nan
499  imgNew[filtInf] = np.inf
500  return imgNew
501 
502 
504  """Task to be used as an ImageMapper for performing
505  A&L decorrelation on subimages on a grid across a A&L difference image.
506 
507  This task subclasses DecorrelateALKernelTask in order to implement
508  all of that task's configuration parameters, as well as its `run` method.
509  """
510 
511  ConfigClass = DecorrelateALKernelConfig
512  _DefaultName = 'ip_diffim_decorrelateALKernelMapper'
513 
514  def __init__(self, *args, **kwargs):
515  DecorrelateALKernelTask.__init__(self, *args, **kwargs)
516 
517  def run(self, subExposure, expandedSubExposure, fullBBox,
518  template, science, alTaskResult=None, psfMatchingKernel=None,
519  preConvKernel=None, **kwargs):
520  """Perform decorrelation operation on `subExposure`, using
521  `expandedSubExposure` to allow for invalid edge pixels arising from
522  convolutions.
523 
524  This method performs A&L decorrelation on `subExposure` using
525  local measures for image variances and PSF. `subExposure` is a
526  sub-exposure of the non-decorrelated A&L diffim. It also
527  requires the corresponding sub-exposures of the template
528  (`template`) and science (`science`) exposures.
529 
530  Parameters
531  ----------
532  subExposure : `lsst.afw.image.Exposure`
533  the sub-exposure of the diffim
534  expandedSubExposure : `lsst.afw.image.Exposure`
535  the expanded sub-exposure upon which to operate
536  fullBBox : `lsst.geom.Box2I`
537  the bounding box of the original exposure
538  template : `lsst.afw.image.Exposure`
539  the corresponding sub-exposure of the template exposure
540  science : `lsst.afw.image.Exposure`
541  the corresponding sub-exposure of the science exposure
542  alTaskResult : `lsst.pipe.base.Struct`
543  the result of A&L image differencing on `science` and
544  `template`, importantly containing the resulting
545  `psfMatchingKernel`. Can be `None`, only if
546  `psfMatchingKernel` is not `None`.
547  psfMatchingKernel : Alternative parameter for passing the
548  A&L `psfMatchingKernel` directly.
549  preConvKernel : If not None, then pre-filtering was applied
550  to science exposure, and this is the pre-convolution
551  kernel.
552  kwargs :
553  additional keyword arguments propagated from
554  `ImageMapReduceTask.run`.
555 
556  Returns
557  -------
558  A `pipeBase.Struct` containing:
559 
560  - ``subExposure`` : the result of the `subExposure` processing.
561  - ``decorrelationKernel`` : the decorrelation kernel, currently
562  not used.
563 
564  Notes
565  -----
566  This `run` method accepts parameters identical to those of
567  `ImageMapper.run`, since it is called from the
568  `ImageMapperTask`. See that class for more information.
569  """
570  templateExposure = template # input template
571  scienceExposure = science # input science image
572  if alTaskResult is None and psfMatchingKernel is None:
573  raise RuntimeError('Both alTaskResult and psfMatchingKernel cannot be None')
574  psfMatchingKernel = alTaskResult.psfMatchingKernel if alTaskResult is not None else psfMatchingKernel
575 
576  # subExp and expandedSubExp are subimages of the (un-decorrelated) diffim!
577  # So here we compute corresponding subimages of templateExposure and scienceExposure
578  subExp2 = scienceExposure.Factory(scienceExposure, expandedSubExposure.getBBox())
579  subExp1 = templateExposure.Factory(templateExposure, expandedSubExposure.getBBox())
580 
581  # Prevent too much log INFO verbosity from DecorrelateALKernelTask.run
582  logLevel = self.log.getLevel()
583  self.log.setLevel(lsst.log.WARN)
584  res = DecorrelateALKernelTask.run(self, subExp2, subExp1, expandedSubExposure,
585  psfMatchingKernel, preConvKernel)
586  self.log.setLevel(logLevel) # reset the log level
587 
588  diffim = res.correctedExposure.Factory(res.correctedExposure, subExposure.getBBox())
589  out = pipeBase.Struct(subExposure=diffim, )
590  return out
591 
592 
594  """Configuration parameters for the ImageMapReduceTask to direct it to use
595  DecorrelateALKernelMapper as its mapper for A&L decorrelation.
596  """
597  mapper = pexConfig.ConfigurableField(
598  doc='A&L decorrelation task to run on each sub-image',
599  target=DecorrelateALKernelMapper
600  )
601 
602 
603 class DecorrelateALKernelSpatialConfig(pexConfig.Config):
604  """Configuration parameters for the DecorrelateALKernelSpatialTask.
605  """
606  decorrelateConfig = pexConfig.ConfigField(
607  dtype=DecorrelateALKernelConfig,
608  doc='DecorrelateALKernel config to use when running on complete exposure (non spatially-varying)',
609  )
610 
611  decorrelateMapReduceConfig = pexConfig.ConfigField(
612  dtype=DecorrelateALKernelMapReduceConfig,
613  doc='DecorrelateALKernelMapReduce config to use when running on each sub-image (spatially-varying)',
614  )
615 
616  ignoreMaskPlanes = pexConfig.ListField(
617  dtype=str,
618  doc="""Mask planes to ignore for sigma-clipped statistics""",
619  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
620  )
621 
622  def setDefaults(self):
623  self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.gridStepX = self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.gridStepY = 40
624  self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.cellSizeX = self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.cellSizeY = 41
625  self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.borderSizeX = self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.borderSizeY = 8
626  self.decorrelateMapReduceConfigdecorrelateMapReduceConfig.reducer.reduceOperation = 'average'
627 
628 
629 class DecorrelateALKernelSpatialTask(pipeBase.Task):
630  """Decorrelate the effect of convolution by Alard-Lupton matching kernel in image difference
631 
632  Notes
633  -----
634 
635  Pipe-task that removes the neighboring-pixel covariance in an
636  image difference that are added when the template image is
637  convolved with the Alard-Lupton PSF matching kernel.
638 
639  This task is a simple wrapper around @ref DecorrelateALKernelTask,
640  which takes a `spatiallyVarying` parameter in its `run` method. If
641  it is `False`, then it simply calls the `run` method of @ref
642  DecorrelateALKernelTask. If it is True, then it uses the @ref
643  ImageMapReduceTask framework to break the exposures into
644  subExposures on a grid, and performs the `run` method of @ref
645  DecorrelateALKernelTask on each subExposure. This enables it to
646  account for spatially-varying PSFs and noise in the exposures when
647  performing the decorrelation.
648 
649  This task has no standalone example, however it is applied as a
650  subtask of pipe.tasks.imageDifference.ImageDifferenceTask.
651  There is also an example of its use in `tests/testImageDecorrelation.py`.
652  """
653  ConfigClass = DecorrelateALKernelSpatialConfig
654  _DefaultName = "ip_diffim_decorrelateALKernelSpatial"
655 
656  def __init__(self, *args, **kwargs):
657  """Create the image decorrelation Task
658 
659  Parameters
660  ----------
661  args :
662  arguments to be passed to
663  `lsst.pipe.base.task.Task.__init__`
664  kwargs :
665  additional keyword arguments to be passed to
666  `lsst.pipe.base.task.Task.__init__`
667  """
668  pipeBase.Task.__init__(self, *args, **kwargs)
669 
671  self.statsControlstatsControl.setNumSigmaClip(3.)
672  self.statsControlstatsControl.setNumIter(3)
673  self.statsControlstatsControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.ignoreMaskPlanes))
674 
675  def computeVarianceMean(self, exposure):
676  """Compute the mean of the variance plane of `exposure`.
677  """
678  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
679  exposure.getMaskedImage().getMask(),
680  afwMath.MEANCLIP, self.statsControlstatsControl)
681  var = statObj.getValue(afwMath.MEANCLIP)
682  return var
683 
684  def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel,
685  spatiallyVarying=True, preConvKernel=None, templateMatched=True):
686  """Perform decorrelation of an image difference exposure.
687 
688  Decorrelates the diffim due to the convolution of the
689  templateExposure with the A&L psfMatchingKernel. If
690  `spatiallyVarying` is True, it utilizes the spatially varying
691  matching kernel via the `imageMapReduce` framework to perform
692  spatially-varying decorrelation on a grid of subExposures.
693 
694  Parameters
695  ----------
696  scienceExposure : `lsst.afw.image.Exposure`
697  the science Exposure used for PSF matching
698  templateExposure : `lsst.afw.image.Exposure`
699  the template Exposure used for PSF matching
700  subtractedExposure : `lsst.afw.image.Exposure`
701  the subtracted Exposure produced by `ip_diffim.ImagePsfMatchTask.subtractExposures()`
702  psfMatchingKernel :
703  an (optionally spatially-varying) PSF matching kernel produced
704  by `ip_diffim.ImagePsfMatchTask.subtractExposures()`
705  spatiallyVarying : `bool`
706  if True, perform the spatially-varying operation
707  preConvKernel : `lsst.meas.algorithms.Psf`
708  if not none, the scienceExposure has been pre-filtered with this kernel. (Currently
709  this option is experimental.)
710  templateMatched : `bool`, optional
711  If True, the template exposure was matched (convolved) to the science exposure.
712 
713  Returns
714  -------
715  results : `lsst.pipe.base.Struct`
716  a structure containing:
717 
718  - ``correctedExposure`` : the decorrelated diffim
719 
720  """
721  self.log.info('Running A&L decorrelation: spatiallyVarying=%r', spatiallyVarying)
722 
723  svar = self.computeVarianceMeancomputeVarianceMean(scienceExposure)
724  tvar = self.computeVarianceMeancomputeVarianceMean(templateExposure)
725  if np.isnan(svar) or np.isnan(tvar): # Should not happen unless entire image has been masked.
726  # Double check that one of the exposures is all NaNs
727  if (np.all(np.isnan(scienceExposure.image.array))
728  or np.all(np.isnan(templateExposure.image.array))):
729  self.log.warning('Template or science image is entirely NaNs: skipping decorrelation.')
730  if np.isnan(svar):
731  svar = 1e-9
732  if np.isnan(tvar):
733  tvar = 1e-9
734 
735  var = self.computeVarianceMeancomputeVarianceMean(subtractedExposure)
736 
737  if spatiallyVarying:
738  self.log.info("Variance (science, template): (%f, %f)", svar, tvar)
739  self.log.info("Variance (uncorrected diffim): %f", var)
740  config = self.config.decorrelateMapReduceConfig
741  task = ImageMapReduceTask(config=config)
742  results = task.run(subtractedExposure, science=scienceExposure,
743  template=templateExposure, psfMatchingKernel=psfMatchingKernel,
744  preConvKernel=preConvKernel, forceEvenSized=True,
745  templateMatched=templateMatched)
746  results.correctedExposure = results.exposure
747 
748  # Make sure masks of input image are propagated to diffim
749  def gm(exp):
750  return exp.getMaskedImage().getMask()
751  gm(results.correctedExposure)[:, :] = gm(subtractedExposure)
752 
753  var = self.computeVarianceMeancomputeVarianceMean(results.correctedExposure)
754  self.log.info("Variance (corrected diffim): %f", var)
755 
756  else:
757  config = self.config.decorrelateConfig
758  task = DecorrelateALKernelTask(config=config)
759  results = task.run(scienceExposure, templateExposure,
760  subtractedExposure, psfMatchingKernel, preConvKernel=preConvKernel,
761  templateMatched=templateMatched)
762 
763  return results
def run(self, subExposure, expandedSubExposure, fullBBox, template, science, alTaskResult=None, psfMatchingKernel=None, preConvKernel=None, **kwargs)
def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel, spatiallyVarying=True, preConvKernel=None, templateMatched=True)
def computeCorrection(self, kappa, svar, tvar, preConvArr=None)
def run(self, scienceExposure, templateExposure, subtractedExposure, psfMatchingKernel, preConvKernel=None, xcen=None, ycen=None, svar=None, tvar=None, templateMatched=True)
def padCenterOriginArray(A, tuple newShape, useInverse=False)
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())