lsst.ip.diffim  13.0-22-g3839dbb+8
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Groups Pages
zogy.py
Go to the documentation of this file.
1 from __future__ import absolute_import, division, print_function
2 from future import standard_library
3 standard_library.install_aliases()
4 #
5 # LSST Data Management System
6 # Copyright 2016 AURA/LSST.
7 #
8 # This product includes software developed by the
9 # LSST Project (http://www.lsst.org/).
10 #
11 # This program is free software: you can redistribute it and/or modify
12 # it under the terms of the GNU General Public License as published by
13 # the Free Software Foundation, either version 3 of the License, or
14 # (at your option) any later version.
15 #
16 # This program is distributed in the hope that it will be useful,
17 # but WITHOUT ANY WARRANTY; without even the implied warranty of
18 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
19 # GNU General Public License for more details.
20 #
21 # You should have received a copy of the LSST License Statement and
22 # the GNU General Public License along with this program. If not,
23 # see <https://www.lsstcorp.org/LegalNotices/>.
24 #
25 
26 import numpy as np
27 
28 import lsst.afw.image as afwImage
29 import lsst.afw.geom as afwGeom
30 import lsst.meas.algorithms as measAlg
31 import lsst.afw.math as afwMath
32 import lsst.pex.config as pexConfig
33 import lsst.pipe.base as pipeBase
34 import lsst.log
35 
36 from .imageMapReduce import (ImageMapReduceConfig, ImageMapperSubtask,
37  ImageMapperSubtaskConfig)
38 
39 __all__ = ["ZogyTask", "ZogyConfig",
40  "ZogyMapperSubtask", "ZogyMapReduceConfig"]
41 
42 
43 """Tasks for performing the "Proper image subtraction" algorithm of
44 Zackay, et al. (2016), hereafter simply referred to as 'ZOGY (2016)'.
45 
46 `ZogyTask` contains methods to perform the basic estimation of the
47 ZOGY diffim `D`, its updated PSF, and the variance-normalized
48 likelihood image `S_corr`. We have implemented ZOGY using the
49 proscribed methodology, computing all convolutions in Fourier space,
50 and also variants in which the convolutions are performed in real
51 (image) space. The former is faster and results in fewer artifacts
52 when the PSFs are noisy (i.e., measured, for example, via
53 `PsfEx`). The latter is presumed to be preferred as it can account for
54 masks correctly with fewer "ringing" artifacts from edge effects or
55 saturated stars, but noisy PSFs result in their own smaller
56 artifacts. Removal of these artifacts is a subject of continuing
57 research. Currently, we "pad" the PSFs when performing the
58 subtractions in real space, which reduces, but does not entirely
59 eliminate these artifacts.
60 
61 All methods in `ZogyTask` assume template and science images are
62 already accurately photometrically and astrometrically registered.
63 
64 `ZogyMapperSubtask` is a wrapper which runs `ZogyTask` in the
65 `ImageMapReduce` framework, computing of ZOGY diffim's on small,
66 overlapping sub-images, thereby enabling complete ZOGY diffim's which
67 account for spatially-varying noise and PSFs across the two input
68 exposures. An example of the use of this task is in the `testZogy.py`
69 unit test.
70 """
71 
72 
73 class ZogyConfig(pexConfig.Config):
74  """Configuration parameters for the ZogyTask
75  """
76  inImageSpace = pexConfig.Field(
77  dtype=bool,
78  default=False,
79  doc="""Perform all convolutions in real (image) space rather than Fourier space.
80  Currently if True, this results in artifacts when using real (noisy) PSFs."""
81  )
82 
83  padSize = pexConfig.Field(
84  dtype=int,
85  default=7,
86  doc="""Number of pixels to pad PSFs to avoid artifacts (when inImageSpace is True)"""
87  )
88 
89  templateFluxScaling = pexConfig.Field(
90  dtype=float,
91  default=1.,
92  doc="""Template flux scaling factor (Fr in ZOGY paper)"""
93  )
94 
95  scienceFluxScaling = pexConfig.Field(
96  dtype=float,
97  default=1.,
98  doc="""Science flux scaling factor (Fn in ZOGY paper)"""
99  )
100 
101  doTrimKernels = pexConfig.Field(
102  dtype=bool,
103  default=False,
104  doc="""Trim kernels for image-space ZOGY. Speeds up convolutions and shrinks artifacts.
105  Subject of future research."""
106  )
107 
108  doFilterPsfs = pexConfig.Field(
109  dtype=bool,
110  default=True,
111  doc="""Filter PSFs for image-space ZOGY. Aids in reducing artifacts.
112  Subject of future research."""
113  )
114 
115  ignoreMaskPlanes = pexConfig.ListField(
116  dtype=str,
117  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE"),
118  doc="""Mask planes to ignore for statistics"""
119  )
120 
121 
122 class ZogyTask(pipeBase.Task):
123  """Task to perform ZOGY proper image subtraction. See module-level documentation for
124  additional details.
125 
126  In all methods, im1 is R (reference, or template) and im2 is N (new, or science).
127  """
128  ConfigClass = ZogyConfig
129  _DefaultName = "ip_diffim_Zogy"
130 
131  def __init__(self, templateExposure, scienceExposure, sig1=None, sig2=None,
132  psf1=None, psf2=None, *args, **kwargs):
133  """Create the ZOGY task.
134 
135  Parameters
136  ----------
137  templateExposure : lsst.afw.image.Exposure
138  Template exposure ("Reference image" in ZOGY (2016)).
139  scienceExposure : lsst.afw.image.Exposure
140  Science exposure ("New image" in ZOGY (2016)). Must have already been
141  registered and photmetrically matched to template.
142  sig1 : float
143  (Optional) sqrt(variance) of `templateExposure`. If `None`, it is
144  computed from the sqrt(mean) of the `templateExposure` variance image.
145  sig2 : float
146  (Optional) sqrt(variance) of `scienceExposure`. If `None`, it is
147  computed from the sqrt(mean) of the `scienceExposure` variance image.
148  psf1 : 2D numpy.array
149  (Optional) 2D array containing the PSF image for the template. If
150  `None`, it is extracted from the PSF taken at the center of `templateExposure`.
151  psf2 : 2D numpy.array
152  (Optional) 2D array containing the PSF image for the science img. If
153  `None`, it is extracted from the PSF taken at the center of `scienceExposure`.
154  args :
155  additional arguments to be passed to
156  `lsst.pipe.base.task.Task.__init__`
157  kwargs :
158  additional keyword arguments to be passed to
159  `lsst.pipe.base.task.Task.__init__`
160  """
161  pipeBase.Task.__init__(self, *args, **kwargs)
162 
163  self.template = templateExposure
164  self.science = scienceExposure
165 
166  self.statsControl = afwMath.StatisticsControl()
167  self.statsControl.setNumSigmaClip(3.)
168  self.statsControl.setNumIter(3)
169  self.statsControl.setAndMask(afwImage.Mask\
170  .getPlaneBitMask(self.config.ignoreMaskPlanes))
171 
172  self.im1 = self.template.getMaskedImage().getImage().getArray()
173  self.im2 = self.science.getMaskedImage().getImage().getArray()
174  self.im1_var = self.template.getMaskedImage().getVariance().getArray()
175  self.im2_var = self.science.getMaskedImage().getVariance().getArray()
176 
177  def selectPsf(psf, exposure):
178  if psf is not None:
179  return psf
180  else:
181  bbox1 = self.template.getBBox()
182  xcen = (bbox1.getBeginX() + bbox1.getEndX()) / 2.
183  ycen = (bbox1.getBeginY() + bbox1.getEndY()) / 2.
184  return exposure.getPsf().computeKernelImage(afwGeom.Point2D(xcen, ycen)).getArray()
185 
186  self.im1_psf = selectPsf(psf1, self.template)
187  self.im2_psf = selectPsf(psf2, self.science)
188 
189  # Make sure PSFs are the same size. Assume they're square...
190  if self.im1_psf.shape[0] < self.im2_psf.shape[0]:
191  self.im1_psf = np.pad(self.im1_psf, (((self.im2_psf.shape[0] - self.im1_psf.shape[0])//2,
192  (self.im2_psf.shape[1] - self.im1_psf.shape[1])//2)),
193  mode='constant', constant_values=0)
194  elif self.im2_psf.shape[0] < self.im1_psf.shape[0]:
195  self.im2_psf = np.pad(self.im2_psf, (((self.im1_psf.shape[0] - self.im2_psf.shape[0])//2,
196  (self.im1_psf.shape[1] - self.im2_psf.shape[1])//2)),
197  mode='constant', constant_values=0)
198 
199  self.sig1 = np.sqrt(self._computeVarianceMean(self.template)) if sig1 is None else sig1
200  self.sig2 = np.sqrt(self._computeVarianceMean(self.science)) if sig2 is None else sig2
201 
202  self.Fr = self.config.templateFluxScaling # default is 1
203  self.Fn = self.config.scienceFluxScaling # default is 1
204  self.padSize = self.config.padSize # default is 7
205 
206  def _computeVarianceMean(self, exposure):
207  """Compute the sigma-clipped mean of the variance image of `exposure`.
208  """
209  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
210  exposure.getMaskedImage().getMask(),
211  afwMath.MEANCLIP, self.statsControl)
212  var = statObj.getValue(afwMath.MEANCLIP)
213  return var
214 
215  @staticmethod
216  def _padPsfToSize(psf, size):
217  """Zero-pad `psf` to the dimensions given by `size`.
218 
219  Parameters
220  ----------
221  psf : 2D numpy.array
222  Input psf to be padded
223  size : list
224  Two element list containing the dimensions to pad the `psf` to
225 
226  Returns
227  -------
228  psf : 2D numpy.array
229  The padded copy of the input `psf`.
230  """
231  padSize0 = size[0] # im.shape[0]//2 - psf.shape[0]//2
232  padSize1 = size[1] # im.shape[1]//2 - psf.shape[1]//2
233  # Hastily assume the psf is odd-sized...
234  if padSize0 > 0 or padSize1 > 0:
235  if padSize0 < 0:
236  padSize0 = 0
237  if padSize1 < 0:
238  padSize1 = 0
239  psf = np.pad(psf, ((padSize0, padSize0-1), (padSize1, padSize1-1)), mode='constant',
240  constant_values=0)
241  return psf
242 
243  @staticmethod
244  def _padPsfToImageSize(psf, im):
245  """Zero-pad `psf` to same dimensions as im.
246 
247  Parameters
248  ----------
249  psf : 2D numpy.array
250  Input psf to be padded
251  im : lsst.afw.Image, MaskedImage or Exposure
252  Dimensions of this image are used to set the PSF pad dimensions
253 
254  Returns
255  -------
256  psf : 2D numpy.array
257  The padded copy of the input `psf`.
258  """
259  return ZogyTask._padPsfToSize(psf, (im.shape[0]//2 - psf.shape[0]//2,
260  im.shape[1]//2 - psf.shape[1]//2))
261 
262  def computePrereqs(self, psf1=None, psf2=None, padSize=0):
263  """Compute standard ZOGY quantities used by (nearly) all methods.
264 
265  Many of the ZOGY calculations require similar quantities, including
266  FFTs of the PSFs, and the "denominator" term (e.g. in eq. 13 of
267  ZOGY manuscript (2016). This function consolidates many of those
268  operations.
269 
270  Parameters
271  ----------
272  psf1 : 2D numpy.array
273  (Optional) Input psf of template, override if already padded
274  psf2 : 2D numpy.array
275  (Optional) Input psf of science image, override if already padded
276 
277  Returns
278  -------
279  A lsst.pipe.base.Struct containing:
280  - Pr : 2D numpy.array, the (possibly zero-padded) template PSF
281  - Pn : 2D numpy.array, the (possibly zero-padded) science PSF
282  - Pr_hat : 2D numpy.array, the FFT of `Pr`
283  - Pn_hat : 2D numpy.array, the FFT of `Pn`
284  - denom : 2D numpy.array, the denominator of equation (13) in ZOGY (2016) manuscript
285  - Fd : float, the relative flux scaling factor between science and template
286  """
287  psf1 = self.im1_psf if psf1 is None else psf1
288  psf2 = self.im2_psf if psf2 is None else psf2
289  padSize = self.padSize if padSize is None else padSize
290  Pr, Pn = psf1, psf2
291  if padSize > 0:
292  Pr = ZogyTask._padPsfToSize(psf1, (padSize, padSize))
293  Pn = ZogyTask._padPsfToSize(psf2, (padSize, padSize))
294 
295  sigR, sigN = self.sig1, self.sig2
296  Pr_hat = np.fft.fft2(Pr)
297  Pr_hat2 = np.conj(Pr_hat) * Pr_hat # np.abs(Pr_hat)**2)
298  Pn_hat = np.fft.fft2(Pn)
299  Pn_hat2 = np.conj(Pn_hat) * Pn_hat # np.abs(Pn_hat)**2))
300  denom = np.sqrt((sigN**2 * self.Fr**2 * Pr_hat2) + (sigR**2 * self.Fn**2 * Pn_hat2))
301  Fd = self.Fr*self.Fn / np.sqrt(sigN**2 * self.Fr**2 + sigR**2 * self.Fn**2)
302 
303  res = pipeBase.Struct(
304  Pr=Pr, Pn=Pn, Pr_hat=Pr_hat, Pn_hat=Pn_hat, denom=denom, Fd=Fd
305  )
306  return res
307 
308  # In all functions, im1 is R (reference, or template) and im2 is N (new, or science)
309  def computeDiffimFourierSpace(self, debug=False, **kwargs):
310  """Compute ZOGY diffim `D` as proscribed in ZOGY (2016) manuscript
311 
312  Compute the ZOGY eqn. (13):
313  $$
314  \widehat{D} = \frac{Fr\widehat{Pr}\widehat{N} -
315  F_n\widehat{Pn}\widehat{R}}{\sqrt{\sigma_n^2 Fr^2
316  |\widehat{Pr}|^2 + \sigma_r^2 F_n^2 |\widehat{Pn}|^2}}
317  $$
318  where $D$ is the optimal difference image, $R$ and $N$ are the
319  reference and "new" image, respectively, $Pr$ and $P_n$ are their
320  PSFs, $Fr$ and $Fn$ are their flux-based zero-points (which we
321  will set to one here), $\sigma_r^2$ and $\sigma_n^2$ are their
322  variance, and $\widehat{D}$ denotes the FT of $D$.
323 
324  Returns
325  -------
326  A lsst.pipe.base.Struct containing:
327  - D : 2D numpy.array, the proper image difference
328  - D_var : 2D numpy.array, the variance image for `D`
329  """
330  # Do all in fourier space (needs image-sized PSFs)
331  psf1 = ZogyTask._padPsfToImageSize(self.im1_psf, self.im1)
332  psf2 = ZogyTask._padPsfToImageSize(self.im2_psf, self.im2)
333 
334  preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs
335 
336  def _filterKernel(K, trim_amount):
337  # Filter the wings of Kn, Kr, set to zero
338  ps = trim_amount
339  K[:ps, :] = K[-ps:, :] = 0
340  K[:, :ps] = K[:, -ps:] = 0
341  return K
342 
343  Kr_hat = self.Fr * preqs.Pr_hat / preqs.denom
344  Kn_hat = self.Fn * preqs.Pn_hat / preqs.denom
345  if debug and self.config.doTrimKernels: # default False
346  # Suggestion from Barak to trim Kr and Kn to remove artifacts
347  # Here we just filter them (in image space) to keep them the same size
348  ps = (Kn_hat.shape[1] - 80)//2
349  Kn = _filterKernel(np.fft.ifft2(Kn_hat), ps)
350  Kn_hat = np.fft.fft2(Kn)
351  Kr = _filterKernel(np.fft.ifft2(Kr_hat), ps)
352  Kr_hat = np.fft.fft2(Kr)
353 
354  def processImages(im1, im2, doAdd=False):
355  # Some masked regions are NaN or infinite!, and FFTs no likey.
356  im1[np.isinf(im1)] = np.nan
357  im1[np.isnan(im1)] = np.nanmean(im1)
358  im2[np.isinf(im2)] = np.nan
359  im2[np.isnan(im2)] = np.nanmean(im2)
360 
361  R_hat = np.fft.fft2(im1)
362  N_hat = np.fft.fft2(im2)
363 
364  D_hat = Kr_hat * N_hat
365  if not doAdd:
366  D_hat -= Kn_hat * R_hat
367  else:
368  D_hat += Kn_hat * R_hat
369 
370  D = np.fft.ifft2(D_hat)
371  D = np.fft.ifftshift(D.real) / preqs.Fd
372  return D
373 
374  # First do the image
375  D = processImages(self.im1, self.im2, doAdd=False)
376  # Do the exact same thing to the var images, except add them
377  D_var = processImages(self.im1_var, self.im2_var, doAdd=True)
378 
379  return pipeBase.Struct(D=D, D_var=D_var)
380 
381  def _doConvolve(self, exposure, kernel, recenterKernel=False):
382  """! Convolve an Exposure with a decorrelation convolution kernel.
383 
384  Parameters
385  ----------
386  exposure : lsst.afw.image.Exposure to be convolved.
387  kernel : 2D numpy.array to convolve the image with
388 
389  Returns
390  -------
391  A new lsst.afw.image.Exposure with the convolved pixels and the (possibly
392  re-centered) kernel.
393 
394  Notes
395  -----
396  - We optionally re-center the kernel if necessary and return the possibly
397  re-centered kernel
398  """
399  kernelImg = afwImage.ImageD(kernel.shape[0], kernel.shape[1])
400  kernelImg.getArray()[:, :] = kernel
401  kern = afwMath.FixedKernel(kernelImg)
402  if recenterKernel:
403  maxloc = np.unravel_index(np.argmax(kernel), kernel.shape)
404  kern.setCtrX(maxloc[0])
405  kern.setCtrY(maxloc[1])
406  outExp = exposure.clone() # Do this to keep WCS, PSF, masks, etc.
407  convCntrl = afwMath.ConvolutionControl(doNormalize=False, doCopyEdge=False,
408  maxInterpolationDistance=0)
409  try:
410  afwMath.convolve(outExp.getMaskedImage(), exposure.getMaskedImage(), kern, convCntrl)
411  except:
412  # Allow exposure to actually be an image/maskedImage
413  afwMath.convolve(outExp, exposure, kern, convCntrl)
414 
415  return outExp, kern
416 
417  def computeDiffimImageSpace(self, padSize=None, debug=False, **kwargs):
418  """Compute ZOGY diffim `D` using image-space convlutions
419 
420  This method is still being debugged as it results in artifacts
421  when the PSFs are noisy (see module-level docstring). Thus
422  there are several options still enabled by the `debug` flag,
423  which are disabled by defult.
424 
425  Parameters
426  ----------
427  padSize : int, the amount to pad the PSFs by
428  debug : bool, flag to enable debugging tests and options
429 
430  Returns
431  -------
432  D : lsst.afw.Exposure
433  the proper image difference, including correct variance,
434  masks, and PSF
435  """
436  preqs = self.computePrereqs(padSize=padSize)
437 
438  delta = 0.
439  if debug:
440  delta = 1. # Regularize the ratio, a possible option to remove artifacts
441  Kr_hat = (preqs.Pr_hat + delta) / (preqs.denom + delta)
442  Kn_hat = (preqs.Pn_hat + delta) / (preqs.denom + delta)
443  Kr = np.fft.ifft2(Kr_hat).real
444  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
445  Kn = np.fft.ifft2(Kn_hat).real
446  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
447 
448  def _trimKernel(self, K, trim_amount):
449  # Trim out the wings of Kn, Kr (see notebook #15)
450  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
451  ps = trim_amount
452  K = K[ps:-ps, ps:-ps]
453  return K
454 
455  padSize = self.padSize if padSize is None else padSize
456  # Enabling this block (debug=True) makes it slightly faster, but ~25% worse artifacts:
457  if debug and self.config.doTrimKernels: # default False
458  # Filtering also makes it slightly faster (zeros are ignored in convolution)
459  # but a bit worse. Filter the wings of Kn, Kr (see notebook #15)
460  Kn = _trimKernel(Kn, padSize)
461  Kr = _trimKernel(Kr, padSize)
462 
463  # Note these are reverse-labelled, this is CORRECT!
464  exp1, _ = self._doConvolve(self.template, Kn)
465  exp2, _ = self._doConvolve(self.science, Kr)
466  D = exp2
467  tmp = D.getMaskedImage()
468  tmp -= exp1.getMaskedImage()
469  tmp /= preqs.Fd
470  return D
471 
472  def _setNewPsf(self, exposure, psfArr):
473  """Utility method to set an exposure's PSF when provided as a 2-d numpy.array
474  """
475  psfI = afwImage.ImageD(psfArr.shape[0], psfArr.shape[1])
476  psfI.getArray()[:, :] = psfArr
477  psfK = afwMath.FixedKernel(psfI)
478  psfNew = measAlg.KernelPsf(psfK)
479  exposure.setPsf(psfNew)
480  return exposure
481 
482  def computeDiffim(self, inImageSpace=None, padSize=None, **kwargs):
483  """Wrapper method to compute ZOGY proper diffim
484 
485  This method should be used as the public interface for
486  computing the ZOGY diffim.
487 
488  Parameters
489  ----------
490  inImageSpace : bool
491  Override config `inImageSpace` parameter
492  padSize : int
493  Override config `padSize` parameter
494  **kwargs : dict
495  additional keyword arguments to be passed to
496  `computeDiffimFourierSpace` or `computeDiffimImageSpace`.
497 
498  Returns
499  -------
500  D : lsst.afw.Exposure
501  the proper image difference, including correct variance,
502  masks, and PSF
503  """
504  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
505  if inImageSpace:
506  padSize = self.padSize if padSize is None else padSize
507  D = self.computeDiffimImageSpace(padSize=padSize, **kwargs)
508  else:
509  res = self.computeDiffimFourierSpace(**kwargs)
510  D = self.science.clone()
511  D.getMaskedImage().getImage().getArray()[:, :] = res.D
512  D.getMaskedImage().getVariance().getArray()[:, :] = res.D_var
513 
514  psf = self.computeDiffimPsf()
515  D = self._setNewPsf(D, psf)
516  return D
517 
518  def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None):
519  """Compute the ZOGY diffim PSF (ZOGY manuscript eq. 14)
520 
521  Parameters
522  ----------
523  padSize : int
524  Override config `padSize` parameter
525  keepFourier : bool
526  Return the FFT of the diffim PSF (do not inverse-FFT it)
527  psf1 : 2D numpy.array
528  (Optional) Input psf of template, override if already padded
529  psf2 : 2D numpy.array
530  (Optional) Input psf of science image, override if already padded
531 
532  Returns
533  -------
534  Pd : 2D numpy.array, the diffim PSF (or FFT of PSF if `keepFourier=True`)
535  """
536  preqs = self.computePrereqs(psf1=psf1, psf2=psf2, padSize=padSize)
537 
538  Pd_hat_numerator = (self.Fr * self.Fn * preqs.Pr_hat * preqs.Pn_hat)
539  Pd_hat = Pd_hat_numerator / (preqs.Fd * preqs.denom)
540 
541  if keepFourier:
542  return Pd_hat
543 
544  Pd = np.fft.ifft2(Pd_hat)
545  Pd = np.fft.ifftshift(Pd).real
546 
547  return Pd
548 
549  def _computeVarAstGradients(self, xVarAst=0., yVarAst=0., inImageSpace=False,
550  R_hat=None, Kr_hat=None, Kr=None,
551  N_hat=None, Kn_hat=None, Kn=None):
552  """Compute the astrometric noise correction terms
553 
554  Compute the correction for estimated astrometric noise as
555  proscribed in ZOGY (2016), section 3.3. All convolutions
556  performed either in real (image) or Fourier space.
557 
558  Parameters
559  ----------
560  xVarAst, yVarAst : float
561  estimated astrometric noise (variance of astrometric registration errors)
562  inImageSpace : bool
563  Perform all convolutions in real (image) space rather than Fourier space
564  R_hat : 2-D numpy.array
565  (Optional) FFT of template image, only required if `inImageSpace=False`
566  Kr_hat : 2-D numpy.array
567  FFT of Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=False`
568  Kr : 2-D numpy.array
569  Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=True`.
570  Kr is associated with the template (reference).
571  N_hat : 2-D numpy.array
572  FFT of science image, only required if `inImageSpace=False`
573  Kn_hat : 2-D numpy.array
574  FFT of Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=False`
575  Kn : 2-D numpy.array
576  Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=True`.
577  Kn is associated with the science (new) image.
578 
579  Returns
580  -------
581  VastSR, VastSN : 2-D numpy.arrays containing the values in eqs. 30 and 32 of
582  ZOGY (2016).
583  """
584  VastSR = VastSN = 0.
585  if xVarAst + yVarAst > 0: # Do the astrometric variance correction
586  if inImageSpace:
587  S_R, _ = self._doConvolve(self.template, Kr)
588  S_R = S_R.getMaskedImage().getImage().getArray()
589  else:
590  S_R = np.fft.ifft2(R_hat * Kr_hat)
591  gradRx, gradRy = np.gradient(S_R)
592  VastSR = xVarAst * gradRx**2. + yVarAst * gradRy**2.
593 
594  if inImageSpace:
595  S_N, _ = self._doConvolve(self.science, Kn)
596  S_N = S_N.getMaskedImage().getImage().getArray()
597  else:
598  S_N = np.fft.ifft2(N_hat * Kn_hat)
599  gradNx, gradNy = np.gradient(S_N)
600  VastSN = xVarAst * gradNx**2. + yVarAst * gradNy**2.
601 
602  return VastSR, VastSN
603 
604  def computeScorrFourierSpace(self, xVarAst=0., yVarAst=0., **kwargs):
605  """Compute corrected likelihood image, optimal for source detection
606 
607  Compute ZOGY S_corr image. This image can be thresholded for
608  detection without optimal filtering, and the variance image is
609  corrected to account for astrometric noise (errors in
610  astrometric registration whether systematic or due to effects
611  such as DCR). The calculations here are all performed in
612  Fourier space, as proscribed in ZOGY (2016).
613 
614  Parameters
615  ----------
616  xVarAst, yVarAst : float
617  estimated astrometric noise (variance of astrometric registration errors)
618 
619  Returns
620  -------
621  A lsst.pipe.base.Struct containing:
622  - S : numpy.array, the likelihood image S (eq. 12 of ZOGY (2016))
623  - S_var : the corrected variance image (denominator of eq. 25 of ZOGY (2016))
624  - Dpsf : the PSF of the diffim D, likely never to be used.
625  """
626  # Do all in fourier space (needs image-sized PSFs)
627  psf1 = ZogyTask._padPsfToImageSize(self.im1_psf, self.im1)
628  psf2 = ZogyTask._padPsfToImageSize(self.im2_psf, self.im2)
629 
630  preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs
631 
632  # Compute D_hat here (don't need D then, for speed)
633  R_hat = np.fft.fft2(self.im1)
634  N_hat = np.fft.fft2(self.im2)
635  D_hat = self.Fr * preqs.Pr_hat * N_hat - self.Fn * preqs.Pn_hat * R_hat
636  D_hat /= preqs.denom
637 
638  Pd_hat = self.computeDiffimPsf(padSize=0, keepFourier=True, psf1=psf1, psf2=psf2)
639  Pd_bar = np.conj(Pd_hat)
640  S = np.fft.ifft2(D_hat * Pd_bar)
641 
642  # Adjust the variance planes of the two images to contribute to the final detection
643  # (eq's 26-29).
644  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat # np.abs(preqs.Pn_hat)**2.
645  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
646  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat # np.abs(preqs.Pr_hat)**2.
647  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
648 
649  Kr_hat2 = np.fft.fft2(np.fft.ifft2(Kr_hat)**2.)
650  Kn_hat2 = np.fft.fft2(np.fft.ifft2(Kn_hat)**2.)
651  var1c_hat = Kr_hat2 * np.fft.fft2(self.im1_var)
652  var2c_hat = Kn_hat2 * np.fft.fft2(self.im2_var)
653 
654  # Do the astrometric variance correction
655  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=False,
656  R_hat=R_hat, Kr_hat=Kr_hat,
657  N_hat=N_hat, Kn_hat=Kn_hat)
658 
659  S_var = np.sqrt(np.fft.ifftshift(np.fft.ifft2(var1c_hat + var2c_hat)) + fGradR + fGradN)
660  S_var *= preqs.Fd
661 
662  S = np.fft.ifftshift(np.fft.ifft2(Kn_hat * N_hat - Kr_hat * R_hat))
663  S *= preqs.Fd
664 
665  Pd = self.computeDiffimPsf(padSize=0)
666  return pipeBase.Struct(S=S.real, S_var=S_var.real, Dpsf=Pd)
667 
668  def computeScorrImageSpace(self, xVarAst=0., yVarAst=0., padSize=None, **kwargs):
669  """Compute corrected likelihood image, optimal for source detection
670 
671  Compute ZOGY S_corr image. This image can be thresholded for
672  detection without optimal filtering, and the variance image is
673  corrected to account for astrometric noise (errors in
674  astrometric registration whether systematic or due to effects
675  such as DCR). The calculations here are all performed in
676  Real (image) space.
677 
678  Parameters
679  ----------
680  xVarAst, yVarAst : float
681  estimated astrometric noise (variance of astrometric registration errors)
682 
683  Returns
684  -------
685  a tuple containing:
686  - S : lsst.afw.image.Exposure, the likelihood exposure S (eq. 12 of ZOGY (2016)),
687  including corrected variance, masks, and PSF
688  - D : lsst.afw.Exposure, the proper image difference, including correct
689  variance, masks, and PSF
690  """
691  # Do convolutions in image space
692  preqs = self.computePrereqs(padSize=0)
693 
694  padSize = self.padSize if padSize is None else padSize
695  D = self.computeDiffimImageSpace(padSize=padSize)
696  Pd = self.computeDiffimPsf()
697  D = self._setNewPsf(D, Pd)
698  Pd_bar = np.fliplr(np.flipud(Pd))
699  S, _ = self._doConvolve(D, Pd_bar)
700  tmp = S.getMaskedImage()
701  tmp *= preqs.Fd
702 
703  # Adjust the variance planes of the two images to contribute to the final detection
704  # (eq's 26-29).
705  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat # np.abs(Pn_hat)**2. # np.abs(Pn_hat*np.conj(Pn_hat))
706  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
707  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat # np.abs(Pr_hat)**2. # np.abs(Pr_hat*np.conj(Pr_hat))
708  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
709 
710  Kr = np.fft.ifft2(Kr_hat).real
711  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
712  Kn = np.fft.ifft2(Kn_hat).real
713  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
714  var1c, _ = self._doConvolve(self.template.getMaskedImage().getVariance(), Kr**2.)
715  var2c, _ = self._doConvolve(self.science.getMaskedImage().getVariance(), Kn**2.)
716 
717  # Do the astrometric variance correction
718  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=True,
719  Kr=Kr, Kn=Kn)
720 
721  Smi = S.getMaskedImage()
722  Smi *= preqs.Fd
723  S_var = np.sqrt(var1c.getArray() + var2c.getArray() + fGradR + fGradN)
724  S.getMaskedImage().getVariance().getArray()[:, :] = S_var
725  S = self._setNewPsf(S, Pd)
726 
727  return S, D # also return diffim since it was calculated and might be desired
728 
729  def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, **kwargs):
730  """Wrapper method to compute ZOGY corrected likelihood image, optimal for
731  source detection
732 
733  This method should be used as the public interface for
734  computing the ZOGY S_corr.
735 
736  Parameters
737  ----------
738  xVarAst, yVarAst : float
739  estimated astrometric noise (variance of astrometric registration errors)
740  inImageSpace : bool
741  Override config `inImageSpace` parameter
742  padSize : int
743  Override config `padSize` parameter
744 
745  Returns
746  -------
747  S : lsst.afw.image.Exposure, the likelihood exposure S (eq. 12 of ZOGY (2016)),
748  including corrected variance, masks, and PSF
749  """
750  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
751  if inImageSpace:
752  S, _ = self.computeScorrImageSpace(xVarAst=xVarAst, yVarAst=yVarAst, padSize=padSize)
753  else:
754  res = self.computeScorrFourierSpace(xVarAst=xVarAst, yVarAst=yVarAst)
755 
756  S = self.science.clone()
757  S.getMaskedImage().getImage().getArray()[:, :] = res.S
758  S.getMaskedImage().getVariance().getArray()[:, :] = res.S_var
759  S = self._setNewPsf(S, res.Dpsf)
760 
761  return S
762 
763 
764 class ZogyMapperSubtask(ZogyTask, ImageMapperSubtask):
765  """Task to be used as an ImageMapperSubtask for performing
766  ZOGY image subtraction on a grid of subimages.
767  """
768  ConfigClass = ZogyConfig
769  _DefaultName = 'ip_diffim_ZogyMapper'
770 
771  def __init__(self, *args, **kwargs):
772  ImageMapperSubtask.__init__(self, *args, **kwargs)
773 
774  def run(self, subExposure, expandedSubExposure, fullBBox, template,
775  **kwargs):
776  """Perform ZOGY proper image subtraction on sub-images
777 
778  This method performs ZOGY proper image subtraction on
779  `subExposure` using local measures for image variances and
780  PSF. `subExposure` is a sub-exposure of the science image. It
781  also requires the corresponding sub-exposures of the template
782  (`template`). The operations are actually performed on
783  `expandedSubExposure` to allow for invalid edge pixels arising
784  from convolutions, which are then removed.
785 
786  Parameters
787  ----------
788  subExposure : lsst.afw.image.Exposure
789  the sub-exposure of the diffim
790  expandedSubExposure : lsst.afw.image.Exposure
791  the expanded sub-exposure upon which to operate
792  fullBBox : lsst.afw.geom.BoundingBox
793  the bounding box of the original exposure
794  template : lsst.afw.image.Exposure
795  the template exposure, from which a corresponding sub-exposure
796  is extracted
797  kwargs :
798  additional keyword arguments propagated from
799  `ImageMapReduceTask.run`. These include:
800  - doScorr : bool
801  Compute and return the corrected likelihood image S_corr
802  rather than the proper image difference
803  - inImageSpace : bool
804  Perform all convolutions in real (image) space rather than
805  in Fourier space. This option currently leads to artifacts
806  when using real (measured and noisy) PSFs, thus it is set
807  to `False` by default.
808  These kwargs may also include arguments to be propagated to
809  `ZogyTask.computeDiffim` and `ZogyTask.computeScorr`.
810 
811  Returns
812  -------
813  A `lsst.pipe.base.Struct` containing the result of the
814  `subExposure` processing, labelled 'subExposure'. In this case
815  it is either the subExposure of the proper image difference D,
816  or (if `doScorr==True`) the corrected likelihood exposure `S`.
817 
818  Notes
819  -----
820  This `run` method accepts parameters identical to those of
821  `ImageMapperSubtask.run`, since it is called from the
822  `ImageMapperTask`. See that class for more information.
823  """
824  bbox = subExposure.getBBox()
825  center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.)
826  center = afwGeom.Point2D(center[0], center[1])
827 
828  imageSpace = kwargs.pop('inImageSpace', False)
829  doScorr = kwargs.pop('doScorr', False)
830  sigmas = kwargs.pop('sigmas', None)
831  padSize = kwargs.pop('padSize', 7)
832 
833  # Psf and image for science img (index 2)
834  subExp2 = expandedSubExposure
835 
836  # Psf and image for template img (index 1)
837  subExp1 = template.Factory(template, expandedSubExposure.getBBox())
838 
839  if sigmas is None:
840  sig1 = sig2 = None
841  else:
842  # for testing, can use the input sigma (e.g., global value for entire exposure)
843  sig1, sig2 = sigmas[0], sigmas[1]
844 
845  def _makePsfSquare(psf):
846  if psf.shape[0] < psf.shape[1]:
847  # Sometimes CoaddPsf does this. Make it square.
848  psf = np.pad(psf, ((1, 1), (0, 0)), mode='constant')
849  elif psf.shape[0] > psf.shape[1]:
850  psf = np.pad(psf, ((0, 0), (1, 1)), mode='constant')
851  return psf
852 
853  psf2 = subExp2.getPsf().computeKernelImage(center).getArray()
854  psf2 = _makePsfSquare(psf2)
855 
856  psf1 = template.getPsf().computeKernelImage(center).getArray()
857  psf1 = _makePsfSquare(psf1)
858 
859  # from diffimTests.diffimTests ...
860  if subExp1.getDimensions()[0] < psf1.shape[0] or subExp1.getDimensions()[1] < psf1.shape[1]:
861  return pipeBase.Struct(subExposure=subExposure)
862 
863  def _filterPsf(psf):
864  """Filter a noisy Psf to remove artifacts. Subject of future research."""
865  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
866  if psf.shape[0] == 41: # its from a measured psf
867  psf = psf.copy()
868  psf[psf < 0] = 0
869  psf[0:10, :] = psf[:, 0:10] = psf[31:41, :] = psf[:, 31:41] = 0
870  psf /= psf.sum()
871 
872  return psf
873 
874  psf1b = psf2b = None
875  if self.config.doFilterPsfs: # default True
876  # Note this *really* helps for measured psfs.
877  psf1b = _filterPsf(psf1)
878  psf2b = _filterPsf(psf2)
879 
880  config = ZogyConfig()
881  if imageSpace is True:
882  config.inImageSpace = imageSpace
883  config.padSize = padSize # Don't need padding if doing all in fourier space
884  task = ZogyTask(templateExposure=subExp1, scienceExposure=subExp2,
885  sig1=sig1, sig2=sig2, psf1=psf1b, psf2=psf2b, config=config)
886 
887  if not doScorr:
888  D = task.computeDiffim(**kwargs)
889  else:
890  D = task.computeScorr(**kwargs)
891 
892  outExp = D.Factory(D, subExposure.getBBox())
893  out = pipeBase.Struct(subExposure=outExp)
894  return out
895 
896 
897 class ZogyMapReduceConfig(ImageMapReduceConfig):
898  mapperSubtask = pexConfig.ConfigurableField(
899  doc='Zogy subtask to run on each sub-image',
900  target=ZogyMapperSubtask
901  )
def _doConvolve
Convolve an Exposure with a decorrelation convolution kernel.
Definition: zogy.py:381