Coverage for python/lsst/ip/diffim/imageDecorrelation.py : 16%

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