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