lsst.ip.diffim  13.0-28-gf4bc96c+11
zogy.py
Go to the documentation of this file.
1 from __future__ import absolute_import, division, print_function
2 #
3 # LSST Data Management System
4 # Copyright 2016 AURA/LSST.
5 #
6 # This product includes software developed by the
7 # LSST Project (http://www.lsst.org/).
8 #
9 # This program is free software: you can redistribute it and/or modify
10 # it under the terms of the GNU General Public License as published by
11 # the Free Software Foundation, either version 3 of the License, or
12 # (at your option) any later version.
13 #
14 # This program is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 # GNU General Public License for more details.
18 #
19 # You should have received a copy of the LSST License Statement and
20 # the GNU General Public License along with this program. If not,
21 # see <https://www.lsstcorp.org/LegalNotices/>.
22 #
23 
24 import numpy as np
25 
26 import lsst.afw.image as afwImage
27 import lsst.afw.geom as afwGeom
28 import lsst.meas.algorithms as measAlg
29 import lsst.afw.math as afwMath
30 import lsst.pex.config as pexConfig
31 import lsst.pipe.base as pipeBase
32 
33 from .imageMapReduce import (ImageMapReduceConfig, ImageMapper,
34  ImageMapReduceTask)
35 from .imagePsfMatch import (ImagePsfMatchTask, ImagePsfMatchConfig,
36  subtractAlgorithmRegistry)
37 
38 __all__ = ["ZogyTask", "ZogyConfig",
39  "ZogyMapper", "ZogyMapReduceConfig",
40  "ZogyImagePsfMatchConfig", "ZogyImagePsfMatchTask"]
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 `ZogyMapper` 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 MIN_KERNEL = 1.0e-4
123 
124 
125 class ZogyTask(pipeBase.Task):
126  """Task to perform ZOGY proper image subtraction. See module-level documentation for
127  additional details.
128 
129  In all methods, im1 is R (reference, or template) and im2 is N (new, or science).
130  """
131  ConfigClass = ZogyConfig
132  _DefaultName = "ip_diffim_Zogy"
133 
134  def __init__(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None,
135  psf1=None, psf2=None, *args, **kwargs):
136  """Create the ZOGY task.
137 
138  Parameters
139  ----------
140  templateExposure : `lsst.afw.image.Exposure`
141  Template exposure ("Reference image" in ZOGY (2016)).
142  scienceExposure : `lsst.afw.image.Exposure`
143  Science exposure ("New image" in ZOGY (2016)). Must have already been
144  registered and photmetrically matched to template.
145  sig1 : `float`
146  (Optional) sqrt(variance) of `templateExposure`. If `None`, it is
147  computed from the sqrt(mean) of the `templateExposure` variance image.
148  sig2 : `float`
149  (Optional) sqrt(variance) of `scienceExposure`. If `None`, it is
150  computed from the sqrt(mean) of the `scienceExposure` variance image.
151  psf1 : 2D `numpy.array`
152  (Optional) 2D array containing the PSF image for the template. If
153  `None`, it is extracted from the PSF taken at the center of `templateExposure`.
154  psf2 : 2D `numpy.array`
155  (Optional) 2D array containing the PSF image for the science img. If
156  `None`, it is extracted from the PSF taken at the center of `scienceExposure`.
157  args :
158  additional arguments to be passed to
159  `lsst.pipe.base.task.Task.__init__`
160  kwargs :
161  additional keyword arguments to be passed to
162  `lsst.pipe.base.task.Task.__init__`
163  """
164  pipeBase.Task.__init__(self, *args, **kwargs)
165  self.template = self.science = None
166  self.setup(templateExposure=templateExposure, scienceExposure=scienceExposure,
167  sig1=sig1, sig2=sig2, psf1=psf1, psf2=psf2, *args, **kwargs)
168 
169  def setup(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None,
170  psf1=None, psf2=None, correctBackground=False, *args, **kwargs):
171  """Set up the ZOGY task.
172 
173  Parameters
174  ----------
175  templateExposure : `lsst.afw.image.Exposure`
176  Template exposure ("Reference image" in ZOGY (2016)).
177  scienceExposure : `lsst.afw.image.Exposure`
178  Science exposure ("New image" in ZOGY (2016)). Must have already been
179  registered and photmetrically matched to template.
180  sig1 : `float`
181  (Optional) sqrt(variance) of `templateExposure`. If `None`, it is
182  computed from the sqrt(mean) of the `templateExposure` variance image.
183  sig2 : `float`
184  (Optional) sqrt(variance) of `scienceExposure`. If `None`, it is
185  computed from the sqrt(mean) of the `scienceExposure` variance image.
186  psf1 : 2D `numpy.array`
187  (Optional) 2D array containing the PSF image for the template. If
188  `None`, it is extracted from the PSF taken at the center of `templateExposure`.
189  psf2 : 2D `numpy.array`
190  (Optional) 2D array containing the PSF image for the science img. If
191  `None`, it is extracted from the PSF taken at the center of `scienceExposure`.
192  correctBackground : `bool`
193  (Optional) subtract sigma-clipped mean of exposures. Zogy doesn't correct
194  nonzero backgrounds (unlike AL) so subtract them here.
195  args :
196  additional arguments to be passed to
197  `lsst.pipe.base.task.Task.__init__`
198  kwargs :
199  additional keyword arguments to be passed to
200  `lsst.pipe.base.task.Task.__init__`
201  """
202  if self.template is None and templateExposure is None:
203  return
204  if self.science is None and scienceExposure is None:
205  return
206 
207  self.template = templateExposure
208  self.science = scienceExposure
209 
210  self.statsControl = afwMath.StatisticsControl()
211  self.statsControl.setNumSigmaClip(3.)
212  self.statsControl.setNumIter(3)
213  self.statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(
214  self.config.ignoreMaskPlanes))
215 
216  self.im1 = self.template.getMaskedImage().getImage().getArray()
217  self.im2 = self.science.getMaskedImage().getImage().getArray()
218  self.im1_var = self.template.getMaskedImage().getVariance().getArray()
219  self.im2_var = self.science.getMaskedImage().getVariance().getArray()
220 
221  def selectPsf(psf, exposure):
222  if psf is not None:
223  return psf
224  else:
225  bbox1 = self.template.getBBox()
226  xcen = (bbox1.getBeginX() + bbox1.getEndX()) / 2.
227  ycen = (bbox1.getBeginY() + bbox1.getEndY()) / 2.
228  return exposure.getPsf().computeKernelImage(afwGeom.Point2D(xcen, ycen)).getArray()
229 
230  self.im1_psf = selectPsf(psf1, self.template)
231  self.im2_psf = selectPsf(psf2, self.science)
232 
233  # Make sure PSFs are the same size. Messy, but should work for all cases.
234  psf1 = self.im1_psf
235  psf2 = self.im2_psf
236  pShape1 = psf1.shape
237  pShape2 = psf2.shape
238  if (pShape1[0] < pShape2[0]):
239  psf1 = np.pad(psf1, ((0, pShape2[0] - pShape1[0]), (0, 0)), mode='constant', constant_values=0.)
240  elif (pShape2[0] < pShape1[0]):
241  psf2 = np.pad(psf2, ((0, pShape1[0] - pShape2[0]), (0, 0)), mode='constant', constant_values=0.)
242  if (pShape1[1] < pShape2[1]):
243  psf1 = np.pad(psf1, ((0, 0), (0, pShape2[1] - pShape1[1])), mode='constant', constant_values=0.)
244  elif (pShape2[1] < pShape1[1]):
245  psf2 = np.pad(psf2, ((0, 0), (0, pShape1[1] - pShape2[1])), mode='constant', constant_values=0.)
246 
247  # PSFs' centers may be offset relative to each other; now fix that!
248  maxLoc1 = np.unravel_index(np.argmax(psf1), psf1.shape)
249  maxLoc2 = np.unravel_index(np.argmax(psf2), psf2.shape)
250  # *Very* rarely happens but if they're off by >1 pixel, do it more than once.
251  while (maxLoc1[0] != maxLoc2[0]) or (maxLoc1[1] != maxLoc2[1]):
252  if maxLoc1[0] > maxLoc2[0]:
253  psf2[1:, :] = psf2[:-1, :]
254  elif maxLoc1[0] < maxLoc2[0]:
255  psf1[1:, :] = psf1[:-1, :]
256  if maxLoc1[1] > maxLoc2[1]:
257  psf2[:, 1:] = psf2[:, :-1]
258  elif maxLoc1[1] < maxLoc2[1]:
259  psf1[:, 1:] = psf1[:, :-1]
260  maxLoc1 = np.unravel_index(np.argmax(psf1), psf1.shape)
261  maxLoc2 = np.unravel_index(np.argmax(psf2), psf2.shape)
262 
263  # Make sure there are no div-by-zeros
264  psf1[psf1 < MIN_KERNEL] = MIN_KERNEL
265  psf2[psf2 < MIN_KERNEL] = MIN_KERNEL
266 
267  self.im1_psf = psf1
268  self.im2_psf = psf2
269 
270  self.sig1 = np.sqrt(self._computeVarianceMean(self.template)) if sig1 is None else sig1
271  self.sig2 = np.sqrt(self._computeVarianceMean(self.science)) if sig2 is None else sig2
272  # if sig1 or sig2 are NaN, then the entire region being Zogy-ed is masked.
273  # Don't worry about it - the result will be masked but avoid warning messages.
274  if np.isnan(self.sig1) or self.sig1 == 0:
275  self.sig1 = 1.
276  if np.isnan(self.sig2) or self.sig2 == 0:
277  self.sig2 = 1.
278 
279  # Zogy doesn't correct nonzero backgrounds (unlike AL) so subtract them here.
280  if correctBackground:
281  def _subtractImageMean(exposure):
282  """Compute the sigma-clipped mean of the image of `exposure`."""
283  mi = exposure.getMaskedImage()
284  statObj = afwMath.makeStatistics(mi.getImage(), mi.getMask(),
285  afwMath.MEANCLIP, self.statsControl)
286  mean = statObj.getValue(afwMath.MEANCLIP)
287  if not np.isnan(mean):
288  mi -= mean
289 
290  _subtractImageMean(self.template)
291  _subtractImageMean(self.science)
292 
293  self.Fr = self.config.templateFluxScaling # default is 1
294  self.Fn = self.config.scienceFluxScaling # default is 1
295  self.padSize = self.config.padSize # default is 7
296 
297  def _computeVarianceMean(self, exposure):
298  """Compute the sigma-clipped mean of the variance image of `exposure`.
299  """
300  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
301  exposure.getMaskedImage().getMask(),
302  afwMath.MEANCLIP, self.statsControl)
303  var = statObj.getValue(afwMath.MEANCLIP)
304  return var
305 
306  @staticmethod
307  def _padPsfToSize(psf, size):
308  """Zero-pad `psf` to the dimensions given by `size`.
309 
310  Parameters
311  ----------
312  psf : 2D `numpy.array`
313  Input psf to be padded
314  size : `list`
315  Two element list containing the dimensions to pad the `psf` to
316 
317  Returns
318  -------
319  psf : 2D `numpy.array`
320  The padded copy of the input `psf`.
321  """
322  newArr = np.zeros(size)
323  offset = [size[0]//2 - psf.shape[0]//2 - 1, size[1]//2 - psf.shape[1]//2 - 1]
324  tmp = newArr[offset[0]:(psf.shape[0] + offset[0]), offset[1]:(psf.shape[1] + offset[1])]
325  tmp[:, :] = psf
326  return newArr
327 
328  def computePrereqs(self, psf1=None, psf2=None, padSize=0):
329  """Compute standard ZOGY quantities used by (nearly) all methods.
330 
331  Many of the ZOGY calculations require similar quantities, including
332  FFTs of the PSFs, and the "denominator" term (e.g. in eq. 13 of
333  ZOGY manuscript (2016). This function consolidates many of those
334  operations.
335 
336  Parameters
337  ----------
338  psf1 : 2D `numpy.array`
339  (Optional) Input psf of template, override if already padded
340  psf2 : 2D `numpy.array`
341  (Optional) Input psf of science image, override if already padded
342 
343  Returns
344  -------
345  A lsst.pipe.base.Struct containing:
346  - Pr : 2D `numpy.array`, the (possibly zero-padded) template PSF
347  - Pn : 2D `numpy.array`, the (possibly zero-padded) science PSF
348  - Pr_hat : 2D `numpy.array`, the FFT of `Pr`
349  - Pn_hat : 2D `numpy.array`, the FFT of `Pn`
350  - denom : 2D `numpy.array`, the denominator of equation (13) in ZOGY (2016) manuscript
351  - Fd : `float`, the relative flux scaling factor between science and template
352  """
353  psf1 = self.im1_psf if psf1 is None else psf1
354  psf2 = self.im2_psf if psf2 is None else psf2
355  padSize = self.padSize if padSize is None else padSize
356  Pr, Pn = psf1, psf2
357  if padSize > 0:
358  Pr = ZogyTask._padPsfToSize(psf1, (psf1.shape[0] + padSize, psf1.shape[1] + padSize))
359  Pn = ZogyTask._padPsfToSize(psf2, (psf2.shape[0] + padSize, psf2.shape[1] + padSize))
360  # Make sure there are no div-by-zeros
361  psf1[np.abs(psf1) <= MIN_KERNEL] = MIN_KERNEL
362  psf2[np.abs(psf2) <= MIN_KERNEL] = MIN_KERNEL
363 
364  sigR, sigN = self.sig1, self.sig2
365  Pr_hat = np.fft.fft2(Pr)
366  Pr_hat2 = np.conj(Pr_hat) * Pr_hat
367  Pn_hat = np.fft.fft2(Pn)
368  Pn_hat2 = np.conj(Pn_hat) * Pn_hat
369  denom = np.sqrt((sigN**2 * self.Fr**2 * Pr_hat2) + (sigR**2 * self.Fn**2 * Pn_hat2))
370  Fd = self.Fr * self.Fn / np.sqrt(sigN**2 * self.Fr**2 + sigR**2 * self.Fn**2)
371 
372  res = pipeBase.Struct(
373  Pr=Pr, Pn=Pn, Pr_hat=Pr_hat, Pn_hat=Pn_hat, denom=denom, Fd=Fd
374  )
375  return res
376 
377  # In all functions, im1 is R (reference, or template) and im2 is N (new, or science)
378  def computeDiffimFourierSpace(self, debug=False, returnMatchedTemplate=False, **kwargs):
379  """Compute ZOGY diffim `D` as proscribed in ZOGY (2016) manuscript
380 
381  Compute the ZOGY eqn. (13):
382  $$
383  \widehat{D} = \frac{Fr\widehat{Pr}\widehat{N} -
384  F_n\widehat{Pn}\widehat{R}}{\sqrt{\sigma_n^2 Fr^2
385  |\widehat{Pr}|^2 + \sigma_r^2 F_n^2 |\widehat{Pn}|^2}}
386  $$
387  where $D$ is the optimal difference image, $R$ and $N$ are the
388  reference and "new" image, respectively, $Pr$ and $P_n$ are their
389  PSFs, $Fr$ and $Fn$ are their flux-based zero-points (which we
390  will set to one here), $\sigma_r^2$ and $\sigma_n^2$ are their
391  variance, and $\widehat{D}$ denotes the FT of $D$.
392 
393  Returns
394  -------
395  A `lsst.pipe.base.Struct` containing:
396  - D : 2D `numpy.array`, the proper image difference
397  - D_var : 2D `numpy.array`, the variance image for `D`
398  """
399  # Do all in fourier space (needs image-sized PSFs)
400  psf1 = ZogyTask._padPsfToSize(self.im1_psf, self.im1.shape)
401  psf2 = ZogyTask._padPsfToSize(self.im2_psf, self.im2.shape)
402 
403  preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs
404 
405  def _filterKernel(K, trim_amount):
406  # Filter the wings of Kn, Kr, set to zero
407  ps = trim_amount
408  K[:ps, :] = K[-ps:, :] = 0
409  K[:, :ps] = K[:, -ps:] = 0
410  return K
411 
412  Kr_hat = self.Fr * preqs.Pr_hat / preqs.denom
413  Kn_hat = self.Fn * preqs.Pn_hat / preqs.denom
414  if debug and self.config.doTrimKernels: # default False
415  # Suggestion from Barak to trim Kr and Kn to remove artifacts
416  # Here we just filter them (in image space) to keep them the same size
417  ps = (Kn_hat.shape[1] - 80)//2
418  Kn = _filterKernel(np.fft.ifft2(Kn_hat), ps)
419  Kn_hat = np.fft.fft2(Kn)
420  Kr = _filterKernel(np.fft.ifft2(Kr_hat), ps)
421  Kr_hat = np.fft.fft2(Kr)
422 
423  def processImages(im1, im2, doAdd=False):
424  # Some masked regions are NaN or infinite!, and FFTs no likey.
425  im1[np.isinf(im1)] = np.nan
426  im1[np.isnan(im1)] = np.nanmean(im1)
427  im2[np.isinf(im2)] = np.nan
428  im2[np.isnan(im2)] = np.nanmean(im2)
429 
430  R_hat = np.fft.fft2(im1)
431  N_hat = np.fft.fft2(im2)
432 
433  D_hat = Kr_hat * N_hat
434  D_hat_R = Kn_hat * R_hat
435  if not doAdd:
436  D_hat -= D_hat_R
437  else:
438  D_hat += D_hat_R
439 
440  D = np.fft.ifft2(D_hat)
441  D = np.fft.ifftshift(D.real) / preqs.Fd
442 
443  R = None
444  if returnMatchedTemplate:
445  R = np.fft.ifft2(D_hat_R)
446  R = np.fft.ifftshift(R.real) / preqs.Fd
447 
448  return D, R
449 
450  # First do the image
451  D, R = processImages(self.im1, self.im2, doAdd=False)
452  # Do the exact same thing to the var images, except add them
453  D_var, R_var = processImages(self.im1_var, self.im2_var, doAdd=True)
454 
455  return pipeBase.Struct(D=D, D_var=D_var, R=R, R_var=R_var)
456 
457  def _doConvolve(self, exposure, kernel, recenterKernel=False):
458  """! Convolve an Exposure with a decorrelation convolution kernel.
459 
460  Parameters
461  ----------
462  exposure : `lsst.afw.image.Exposure` to be convolved.
463  kernel : 2D `numpy.array` to convolve the image with
464 
465  Returns
466  -------
467  A new `lsst.afw.image.Exposure` with the convolved pixels and the (possibly
468  re-centered) kernel.
469 
470  Notes
471  -----
472  - We optionally re-center the kernel if necessary and return the possibly
473  re-centered kernel
474  """
475  kernelImg = afwImage.ImageD(kernel.shape[1], kernel.shape[0])
476  kernelImg.getArray()[:, :] = kernel
477  kern = afwMath.FixedKernel(kernelImg)
478  if recenterKernel:
479  maxloc = np.unravel_index(np.argmax(kernel), kernel.shape)
480  kern.setCtrX(maxloc[0])
481  kern.setCtrY(maxloc[1])
482  outExp = exposure.clone() # Do this to keep WCS, PSF, masks, etc.
483  convCntrl = afwMath.ConvolutionControl(doNormalize=False, doCopyEdge=False,
484  maxInterpolationDistance=0)
485  try:
486  afwMath.convolve(outExp.getMaskedImage(), exposure.getMaskedImage(), kern, convCntrl)
487  except:
488  # Allow exposure to actually be an image/maskedImage
489  afwMath.convolve(outExp, exposure, kern, convCntrl)
490 
491  return outExp, kern
492 
493  def computeDiffimImageSpace(self, padSize=None, debug=False, **kwargs):
494  """Compute ZOGY diffim `D` using image-space convlutions
495 
496  This method is still being debugged as it results in artifacts
497  when the PSFs are noisy (see module-level docstring). Thus
498  there are several options still enabled by the `debug` flag,
499  which are disabled by defult.
500 
501  Parameters
502  ----------
503  padSize : `int`, the amount to pad the PSFs by
504  debug : `bool`, flag to enable debugging tests and options
505 
506  Returns
507  -------
508  D : `lsst.afw.Exposure`
509  the proper image difference, including correct variance,
510  masks, and PSF
511  """
512  preqs = self.computePrereqs(padSize=padSize)
513 
514  delta = 0.
515  if debug:
516  delta = 1. # Regularize the ratio, a possible option to remove artifacts
517  Kr_hat = (preqs.Pr_hat + delta) / (preqs.denom + delta)
518  Kn_hat = (preqs.Pn_hat + delta) / (preqs.denom + delta)
519  Kr = np.fft.ifft2(Kr_hat).real
520  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
521  Kn = np.fft.ifft2(Kn_hat).real
522  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
523 
524  def _trimKernel(self, K, trim_amount):
525  # Trim out the wings of Kn, Kr (see notebook #15)
526  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
527  ps = trim_amount
528  K = K[ps:-ps, ps:-ps]
529  return K
530 
531  padSize = self.padSize if padSize is None else padSize
532  # Enabling this block (debug=True) makes it slightly faster, but ~25% worse artifacts:
533  if debug and self.config.doTrimKernels: # default False
534  # Filtering also makes it slightly faster (zeros are ignored in convolution)
535  # but a bit worse. Filter the wings of Kn, Kr (see notebook #15)
536  Kn = _trimKernel(Kn, padSize)
537  Kr = _trimKernel(Kr, padSize)
538 
539  # Note these are reverse-labelled, this is CORRECT!
540  exp1, _ = self._doConvolve(self.template, Kn)
541  exp2, _ = self._doConvolve(self.science, Kr)
542  D = exp2
543  tmp = D.getMaskedImage()
544  tmp -= exp1.getMaskedImage()
545  tmp /= preqs.Fd
546  return pipeBase.Struct(D=D, R=exp1)
547 
548  def _setNewPsf(self, exposure, psfArr):
549  """Utility method to set an exposure's PSF when provided as a 2-d numpy.array
550  """
551  bbox = exposure.getBBox()
552  center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.)
553  center = afwGeom.Point2D(center[0], center[1])
554  psfI = afwImage.ImageD(psfArr.shape[1], psfArr.shape[0])
555  psfI.getArray()[:, :] = psfArr
556  psfK = afwMath.FixedKernel(psfI)
557  psfNew = measAlg.KernelPsf(psfK, center)
558  exposure.setPsf(psfNew)
559  return exposure
560 
561  def computeDiffim(self, inImageSpace=None, padSize=None,
562  returnMatchedTemplate=False, **kwargs):
563  """Wrapper method to compute ZOGY proper diffim
564 
565  This method should be used as the public interface for
566  computing the ZOGY diffim.
567 
568  Parameters
569  ----------
570  inImageSpace : `bool`
571  Override config `inImageSpace` parameter
572  padSize : `int`
573  Override config `padSize` parameter
574  returnMatchedTemplate : bool
575  Include the PSF-matched template in the results Struct
576  **kwargs : `dict`
577  additional keyword arguments to be passed to
578  `computeDiffimFourierSpace` or `computeDiffimImageSpace`.
579 
580  Returns
581  -------
582  An lsst.pipe.base.Struct containing:
583  - D : `lsst.afw.Exposure`
584  the proper image difference, including correct variance,
585  masks, and PSF
586  - R : `lsst.afw.Exposure`
587  If `returnMatchedTemplate` is True, the PSF-matched template
588  exposure
589  """
590  R = None
591  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
592  if inImageSpace:
593  padSize = self.padSize if padSize is None else padSize
594  res = self.computeDiffimImageSpace(padSize=padSize, **kwargs)
595  D = res.D
596  if returnMatchedTemplate:
597  R = res.R
598  else:
599  res = self.computeDiffimFourierSpace(**kwargs)
600  D = self.science.clone()
601  D.getMaskedImage().getImage().getArray()[:, :] = res.D
602  D.getMaskedImage().getVariance().getArray()[:, :] = res.D_var
603  if returnMatchedTemplate:
604  R = self.science.clone()
605  R.getMaskedImage().getImage().getArray()[:, :] = res.R
606  R.getMaskedImage().getVariance().getArray()[:, :] = res.R_var
607 
608  psf = self.computeDiffimPsf()
609  D = self._setNewPsf(D, psf)
610  return pipeBase.Struct(D=D, R=R)
611 
612  def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None):
613  """Compute the ZOGY diffim PSF (ZOGY manuscript eq. 14)
614 
615  Parameters
616  ----------
617  padSize : `int`
618  Override config `padSize` parameter
619  keepFourier : `bool`
620  Return the FFT of the diffim PSF (do not inverse-FFT it)
621  psf1 : 2D `numpy.array`
622  (Optional) Input psf of template, override if already padded
623  psf2 : 2D `numpy.array`
624  (Optional) Input psf of science image, override if already padded
625 
626  Returns
627  -------
628  Pd : 2D `numpy.array`, the diffim PSF (or FFT of PSF if `keepFourier=True`)
629  """
630  preqs = self.computePrereqs(psf1=psf1, psf2=psf2, padSize=padSize)
631 
632  Pd_hat_numerator = (self.Fr * self.Fn * preqs.Pr_hat * preqs.Pn_hat)
633  Pd_hat = Pd_hat_numerator / (preqs.Fd * preqs.denom)
634 
635  if keepFourier:
636  return Pd_hat
637 
638  Pd = np.fft.ifft2(Pd_hat)
639  Pd = np.fft.ifftshift(Pd).real
640 
641  return Pd
642 
643  def _computeVarAstGradients(self, xVarAst=0., yVarAst=0., inImageSpace=False,
644  R_hat=None, Kr_hat=None, Kr=None,
645  N_hat=None, Kn_hat=None, Kn=None):
646  """Compute the astrometric noise correction terms
647 
648  Compute the correction for estimated astrometric noise as
649  proscribed in ZOGY (2016), section 3.3. All convolutions
650  performed either in real (image) or Fourier space.
651 
652  Parameters
653  ----------
654  xVarAst, yVarAst : `float`
655  estimated astrometric noise (variance of astrometric registration errors)
656  inImageSpace : `bool`
657  Perform all convolutions in real (image) space rather than Fourier space
658  R_hat : 2-D `numpy.array`
659  (Optional) FFT of template image, only required if `inImageSpace=False`
660  Kr_hat : 2-D `numpy.array`
661  FFT of Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=False`
662  Kr : 2-D `numpy.array`
663  Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=True`.
664  Kr is associated with the template (reference).
665  N_hat : 2-D `numpy.array`
666  FFT of science image, only required if `inImageSpace=False`
667  Kn_hat : 2-D `numpy.array`
668  FFT of Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=False`
669  Kn : 2-D `numpy.array`
670  Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=True`.
671  Kn is associated with the science (new) image.
672 
673  Returns
674  -------
675  VastSR, VastSN : 2-D `numpy.array`s containing the values in eqs. 30 and 32 of
676  ZOGY (2016).
677  """
678  VastSR = VastSN = 0.
679  if xVarAst + yVarAst > 0: # Do the astrometric variance correction
680  if inImageSpace:
681  S_R, _ = self._doConvolve(self.template, Kr)
682  S_R = S_R.getMaskedImage().getImage().getArray()
683  else:
684  S_R = np.fft.ifft2(R_hat * Kr_hat)
685  gradRx, gradRy = np.gradient(S_R)
686  VastSR = xVarAst * gradRx**2. + yVarAst * gradRy**2.
687 
688  if inImageSpace:
689  S_N, _ = self._doConvolve(self.science, Kn)
690  S_N = S_N.getMaskedImage().getImage().getArray()
691  else:
692  S_N = np.fft.ifft2(N_hat * Kn_hat)
693  gradNx, gradNy = np.gradient(S_N)
694  VastSN = xVarAst * gradNx**2. + yVarAst * gradNy**2.
695 
696  return VastSR, VastSN
697 
698  def computeScorrFourierSpace(self, xVarAst=0., yVarAst=0., **kwargs):
699  """Compute corrected likelihood image, optimal for source detection
700 
701  Compute ZOGY S_corr image. This image can be thresholded for
702  detection without optimal filtering, and the variance image is
703  corrected to account for astrometric noise (errors in
704  astrometric registration whether systematic or due to effects
705  such as DCR). The calculations here are all performed in
706  Fourier space, as proscribed in ZOGY (2016).
707 
708  Parameters
709  ----------
710  xVarAst, yVarAst : `float`
711  estimated astrometric noise (variance of astrometric registration errors)
712 
713  Returns
714  -------
715  A lsst.pipe.base.Struct containing:
716  - S : `numpy.array`, the likelihood image S (eq. 12 of ZOGY (2016))
717  - S_var : the corrected variance image (denominator of eq. 25 of ZOGY (2016))
718  - Dpsf : the PSF of the diffim D, likely never to be used.
719  """
720  # Some masked regions are NaN or infinite!, and FFTs no likey.
721  def fix_nans(im):
722  """Replace any NaNs or Infs with the mean of the image."""
723  isbad = ~np.isfinite(im)
724  if np.any(isbad):
725  im[isbad] = np.nan
726  im[isbad] = np.nanmean(im)
727  return im
728 
729  self.im1 = fix_nans(self.im1)
730  self.im2 = fix_nans(self.im2)
731  self.im1_var = fix_nans(self.im1_var)
732  self.im2_var = fix_nans(self.im2_var)
733 
734  # Do all in fourier space (needs image-sized PSFs)
735  psf1 = ZogyTask._padPsfToSize(self.im1_psf, self.im1.shape)
736  psf2 = ZogyTask._padPsfToSize(self.im2_psf, self.im2.shape)
737 
738  preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs
739 
740  # Compute D_hat here (don't need D then, for speed)
741  R_hat = np.fft.fft2(self.im1)
742  N_hat = np.fft.fft2(self.im2)
743  D_hat = self.Fr * preqs.Pr_hat * N_hat - self.Fn * preqs.Pn_hat * R_hat
744  D_hat /= preqs.denom
745 
746  Pd_hat = self.computeDiffimPsf(padSize=0, keepFourier=True, psf1=psf1, psf2=psf2)
747  Pd_bar = np.conj(Pd_hat)
748  S = np.fft.ifft2(D_hat * Pd_bar)
749 
750  # Adjust the variance planes of the two images to contribute to the final detection
751  # (eq's 26-29).
752  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat
753  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
754  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat
755  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
756 
757  Kr_hat2 = np.fft.fft2(np.fft.ifft2(Kr_hat)**2.)
758  Kn_hat2 = np.fft.fft2(np.fft.ifft2(Kn_hat)**2.)
759  var1c_hat = Kr_hat2 * np.fft.fft2(self.im1_var)
760  var2c_hat = Kn_hat2 * np.fft.fft2(self.im2_var)
761 
762  # Do the astrometric variance correction
763  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=False,
764  R_hat=R_hat, Kr_hat=Kr_hat,
765  N_hat=N_hat, Kn_hat=Kn_hat)
766 
767  S_var = np.sqrt(np.fft.ifftshift(np.fft.ifft2(var1c_hat + var2c_hat)) + fGradR + fGradN)
768  S_var *= preqs.Fd
769 
770  S = np.fft.ifftshift(np.fft.ifft2(Kn_hat * N_hat - Kr_hat * R_hat))
771  S *= preqs.Fd
772 
773  Pd = self.computeDiffimPsf(padSize=0)
774  return pipeBase.Struct(S=S.real, S_var=S_var.real, Dpsf=Pd)
775 
776  def computeScorrImageSpace(self, xVarAst=0., yVarAst=0., padSize=None, **kwargs):
777  """Compute corrected likelihood image, optimal for source detection
778 
779  Compute ZOGY S_corr image. This image can be thresholded for
780  detection without optimal filtering, and the variance image is
781  corrected to account for astrometric noise (errors in
782  astrometric registration whether systematic or due to effects
783  such as DCR). The calculations here are all performed in
784  Real (image) space.
785 
786  Parameters
787  ----------
788  xVarAst, yVarAst : `float`
789  estimated astrometric noise (variance of astrometric registration errors)
790 
791  Returns
792  -------
793  a tuple containing:
794  - S : `lsst.afw.image.Exposure`, the likelihood exposure S (eq. 12 of ZOGY (2016)),
795  including corrected variance, masks, and PSF
796  - D : `lsst.afw.image.Exposure`, the proper image difference, including correct
797  variance, masks, and PSF
798  """
799  # Do convolutions in image space
800  preqs = self.computePrereqs(padSize=0)
801 
802  padSize = self.padSize if padSize is None else padSize
803  D = self.computeDiffimImageSpace(padSize=padSize).D
804  Pd = self.computeDiffimPsf()
805  D = self._setNewPsf(D, Pd)
806  Pd_bar = np.fliplr(np.flipud(Pd))
807  S, _ = self._doConvolve(D, Pd_bar)
808  tmp = S.getMaskedImage()
809  tmp *= preqs.Fd
810 
811  # Adjust the variance planes of the two images to contribute to the final detection
812  # (eq's 26-29).
813  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat
814  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
815  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat
816  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
817 
818  Kr = np.fft.ifft2(Kr_hat).real
819  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
820  Kn = np.fft.ifft2(Kn_hat).real
821  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
822  var1c, _ = self._doConvolve(self.template.getMaskedImage().getVariance(), Kr**2.)
823  var2c, _ = self._doConvolve(self.science.getMaskedImage().getVariance(), Kn**2.)
824 
825  # Do the astrometric variance correction
826  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=True,
827  Kr=Kr, Kn=Kn)
828 
829  Smi = S.getMaskedImage()
830  Smi *= preqs.Fd
831  S_var = np.sqrt(var1c.getArray() + var2c.getArray() + fGradR + fGradN)
832  S.getMaskedImage().getVariance().getArray()[:, :] = S_var
833  S = self._setNewPsf(S, Pd)
834 
835  # also return diffim since it was calculated and might be desired
836  return pipeBase.Struct(S=S, D=D)
837 
838  def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, **kwargs):
839  """Wrapper method to compute ZOGY corrected likelihood image, optimal for
840  source detection
841 
842  This method should be used as the public interface for
843  computing the ZOGY S_corr.
844 
845  Parameters
846  ----------
847  xVarAst, yVarAst : float
848  estimated astrometric noise (variance of astrometric registration errors)
849  inImageSpace : bool
850  Override config `inImageSpace` parameter
851  padSize : int
852  Override config `padSize` parameter
853 
854  Returns
855  -------
856  S : lsst.afw.image.Exposure, the likelihood exposure S (eq. 12 of ZOGY (2016)),
857  including corrected variance, masks, and PSF
858  """
859  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
860  if inImageSpace:
861  res = self.computeScorrImageSpace(xVarAst=xVarAst, yVarAst=yVarAst, padSize=padSize)
862  S = res.S
863  else:
864  res = self.computeScorrFourierSpace(xVarAst=xVarAst, yVarAst=yVarAst)
865 
866  S = self.science.clone()
867  S.getMaskedImage().getImage().getArray()[:, :] = res.S
868  S.getMaskedImage().getVariance().getArray()[:, :] = res.S_var
869  S = self._setNewPsf(S, res.Dpsf)
870 
871  return pipeBase.Struct(S=S)
872 
873 
875  """Task to be used as an ImageMapper for performing
876  ZOGY image subtraction on a grid of subimages.
877  """
878  ConfigClass = ZogyConfig
879  _DefaultName = 'ip_diffim_ZogyMapper'
880 
881  def __init__(self, *args, **kwargs):
882  ImageMapper.__init__(self, *args, **kwargs)
883 
884  def run(self, subExposure, expandedSubExposure, fullBBox, template,
885  **kwargs):
886  """Perform ZOGY proper image subtraction on sub-images
887 
888  This method performs ZOGY proper image subtraction on
889  `subExposure` using local measures for image variances and
890  PSF. `subExposure` is a sub-exposure of the science image. It
891  also requires the corresponding sub-exposures of the template
892  (`template`). The operations are actually performed on
893  `expandedSubExposure` to allow for invalid edge pixels arising
894  from convolutions, which are then removed.
895 
896  Parameters
897  ----------
898  subExposure : lsst.afw.image.Exposure
899  the sub-exposure of the diffim
900  expandedSubExposure : lsst.afw.image.Exposure
901  the expanded sub-exposure upon which to operate
902  fullBBox : lsst.afw.geom.BoundingBox
903  the bounding box of the original exposure
904  template : lsst.afw.image.Exposure
905  the template exposure, from which a corresponding sub-exposure
906  is extracted
907  kwargs :
908  additional keyword arguments propagated from
909  `ImageMapReduceTask.run`. These include:
910  - doScorr : bool
911  Compute and return the corrected likelihood image S_corr
912  rather than the proper image difference
913  - inImageSpace : bool
914  Perform all convolutions in real (image) space rather than
915  in Fourier space. This option currently leads to artifacts
916  when using real (measured and noisy) PSFs, thus it is set
917  to `False` by default.
918  These kwargs may also include arguments to be propagated to
919  `ZogyTask.computeDiffim` and `ZogyTask.computeScorr`.
920 
921  Returns
922  -------
923  A `lsst.pipe.base.Struct` containing the result of the
924  `subExposure` processing, labelled 'subExposure'. In this case
925  it is either the subExposure of the proper image difference D,
926  or (if `doScorr==True`) the corrected likelihood exposure `S`.
927 
928  Notes
929  -----
930  This `run` method accepts parameters identical to those of
931  `ImageMapper.run`, since it is called from the
932  `ImageMapperTask`. See that class for more information.
933  """
934  bbox = subExposure.getBBox()
935  center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.)
936  center = afwGeom.Point2D(center[0], center[1])
937 
938  imageSpace = kwargs.pop('inImageSpace', False)
939  doScorr = kwargs.pop('doScorr', False)
940  sigmas = kwargs.pop('sigmas', None)
941  padSize = kwargs.pop('padSize', 7)
942 
943  # Psf and image for science img (index 2)
944  subExp2 = expandedSubExposure
945 
946  # Psf and image for template img (index 1)
947  subExp1 = template.Factory(template, expandedSubExposure.getBBox())
948 
949  if sigmas is None:
950  sig1 = sig2 = None
951  else:
952  # for testing, can use the input sigma (e.g., global value for entire exposure)
953  sig1, sig2 = sigmas[0], sigmas[1]
954 
955  def _makePsfSquare(psf):
956  # Sometimes CoaddPsf does this. Make it square.
957  if psf.shape[0] < psf.shape[1]:
958  psf = np.pad(psf, ((0, psf.shape[1] - psf.shape[0]), (0, 0)), mode='constant',
959  constant_values=0.)
960  elif psf.shape[0] > psf.shape[1]:
961  psf = np.pad(psf, ((0, 0), (0, psf.shape[0] - psf.shape[1])), mode='constant',
962  constant_values=0.)
963  return psf
964 
965  psf2 = subExp2.getPsf().computeKernelImage(center).getArray()
966  psf2 = _makePsfSquare(psf2)
967 
968  psf1 = template.getPsf().computeKernelImage(center).getArray()
969  psf1 = _makePsfSquare(psf1)
970 
971  # from diffimTests.diffimTests ...
972  if subExp1.getDimensions()[0] < psf1.shape[0] or subExp1.getDimensions()[1] < psf1.shape[1]:
973  return pipeBase.Struct(subExposure=subExposure)
974 
975  def _filterPsf(psf):
976  """Filter a noisy Psf to remove artifacts. Subject of future research."""
977  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
978  if psf.shape[0] == 41: # its from a measured psf
979  psf = psf.copy()
980  psf[psf < 0] = 0
981  psf[0:10, :] = psf[:, 0:10] = psf[31:41, :] = psf[:, 31:41] = 0
982  psf /= psf.sum()
983 
984  return psf
985 
986  psf1b = psf2b = None
987  if self.config.doFilterPsfs: # default True
988  # Note this *really* helps for measured psfs.
989  psf1b = _filterPsf(psf1)
990  psf2b = _filterPsf(psf2)
991 
992  config = ZogyConfig()
993  if imageSpace is True:
994  config.inImageSpace = imageSpace
995  config.padSize = padSize # Don't need padding if doing all in fourier space
996  task = ZogyTask(templateExposure=subExp1, scienceExposure=subExp2,
997  sig1=sig1, sig2=sig2, psf1=psf1b, psf2=psf2b, config=config)
998 
999  if not doScorr:
1000  res = task.computeDiffim(**kwargs)
1001  D = res.D
1002  else:
1003  res = task.computeScorr(**kwargs)
1004  D = res.S
1005 
1006  outExp = D.Factory(D, subExposure.getBBox())
1007  out = pipeBase.Struct(subExposure=outExp)
1008  return out
1009 
1010 
1012  """Config to be passed to ImageMapReduceTask
1013 
1014  This config targets the imageMapper to use the ZogyMapper.
1015  """
1016  mapper = pexConfig.ConfigurableField(
1017  doc='Zogy task to run on each sub-image',
1018  target=ZogyMapper
1019  )
1020 
1021 
1023  """Config for the ZogyImagePsfMatchTask"""
1024 
1025  zogyConfig = pexConfig.ConfigField(
1026  dtype=ZogyConfig,
1027  doc='ZogyTask config to use when running on complete exposure (non spatially-varying)',
1028  )
1029 
1030  zogyMapReduceConfig = pexConfig.ConfigField(
1031  dtype=ZogyMapReduceConfig,
1032  doc='ZogyMapReduce config to use when running Zogy on each sub-image (spatially-varying)',
1033  )
1034 
1035  def setDefaults(self):
1036  self.zogyMapReduceConfig.gridStepX = self.zogyMapReduceConfig.gridStepY = 40
1037  self.zogyMapReduceConfig.cellSizeX = self.zogyMapReduceConfig.cellSizeY = 41
1038  self.zogyMapReduceConfig.borderSizeX = self.zogyMapReduceConfig.borderSizeY = 8
1039  self.zogyMapReduceConfig.reducer.reduceOperation = 'average'
1040  self.zogyConfig.inImageSpace = False
1041 
1042 
1044  """Task to perform Zogy PSF matching and image subtraction.
1045 
1046  This class inherits from ImagePsfMatchTask to contain the _warper
1047  subtask and related methods.
1048  """
1049 
1050  ConfigClass = ZogyImagePsfMatchConfig
1051 
1052  def __init__(self, *args, **kwargs):
1053  ImagePsfMatchTask.__init__(self, *args, **kwargs)
1054 
1055  def _computeImageMean(self, exposure):
1056  """Compute the sigma-clipped mean of the pixels image of `exposure`.
1057  """
1058  statsControl = afwMath.StatisticsControl()
1059  statsControl.setNumSigmaClip(3.)
1060  statsControl.setNumIter(3)
1061  ignoreMaskPlanes = ("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
1062  statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(ignoreMaskPlanes))
1063  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getImage(),
1064  exposure.getMaskedImage().getMask(),
1065  afwMath.MEANCLIP | afwMath.MEDIAN, statsControl)
1066  mn = statObj.getValue(afwMath.MEANCLIP)
1067  med = statObj.getValue(afwMath.MEDIAN)
1068  return mn, med
1069 
1070  def subtractExposures(self, templateExposure, scienceExposure,
1071  doWarping=True, spatiallyVarying=True, inImageSpace=False,
1072  doPreConvolve=False):
1073  """Register, PSF-match, and subtract two Exposures using the ZOGY algorithm.
1074 
1075  Do the following, in order:
1076  - Warp templateExposure to match scienceExposure, if their WCSs do not already match
1077  - Compute subtracted exposure ZOGY image subtraction algorithm on the two exposures
1078 
1079  Parameters
1080  ----------
1081  templateExposure : `lsst.afw.image.Exposure`
1082  exposure to PSF-match to scienceExposure. The exposure's mean value is subtracted
1083  in-place.
1084  scienceExposure : `lsst.afw.image.Exposure`
1085  reference Exposure. The exposure's mean value is subtracted in-place.
1086  doWarping : `bool`
1087  what to do if templateExposure's and scienceExposure's WCSs do not match:
1088  - if True then warp templateExposure to match scienceExposure
1089  - if False then raise an Exception
1090  spatiallyVarying : bool
1091  If True, perform the operation over a grid of patches across the two exposures
1092  inImageSpace : `bool`
1093  If True, perform the Zogy convolutions in image space rather than in frequency space.
1094  doPreConvolve : `bool`
1095  ***Currently not implemented.*** If True assume we are to compute the match filter-convolved
1096  exposure which can be thresholded for detection. In the case of Zogy this would mean
1097  we compute the Scorr image.
1098 
1099  Returns
1100  -------
1101  A `lsst.pipe.base.Struct` containing these fields:
1102  - subtractedExposure: subtracted Exposure
1103  - warpedExposure: templateExposure after warping to match scienceExposure (if doWarping true)
1104  """
1105 
1106  mn1 = self._computeImageMean(templateExposure)
1107  mn2 = self._computeImageMean(scienceExposure)
1108  self.log.info("Exposure means=%f, %f; median=%f, %f:" % (mn1[0], mn2[0], mn1[1], mn2[1]))
1109  if not np.isnan(mn1[0]) and np.abs(mn1[0]) > 1:
1110  mi = templateExposure.getMaskedImage()
1111  mi -= mn1[0]
1112  if not np.isnan(mn2[0]) and np.abs(mn2[0]) > 1:
1113  mi = scienceExposure.getMaskedImage()
1114  mi -= mn2[0]
1115 
1116  self.log.info('Running Zogy algorithm: spatiallyVarying=%r' % spatiallyVarying)
1117 
1118  if not self._validateWcs(templateExposure, scienceExposure):
1119  if doWarping:
1120  self.log.info("Astrometrically registering template to science image")
1121  # Also warp the PSF
1122  xyTransform = afwImage.XYTransformFromWcsPair(scienceExposure.getWcs(),
1123  templateExposure.getWcs())
1124  psfWarped = measAlg.WarpedPsf(templateExposure.getPsf(), xyTransform)
1125  templateExposure = self._warper.warpExposure(scienceExposure.getWcs(),
1126  templateExposure,
1127  destBBox=scienceExposure.getBBox())
1128 
1129  templateExposure.setPsf(psfWarped)
1130  else:
1131  self.log.error("ERROR: Input images not registered")
1132  raise RuntimeError("Input images not registered")
1133 
1134  def gm(exp):
1135  return exp.getMaskedImage().getMask()
1136 
1137  def ga(exp):
1138  return exp.getMaskedImage().getImage().getArray()
1139 
1140  if self.config.zogyConfig.inImageSpace:
1141  inImageSpace = True # Override
1142  self.log.info('Running Zogy algorithm: inImageSpace=%r' % inImageSpace)
1143  if spatiallyVarying:
1144  config = self.config.zogyMapReduceConfig
1145  task = ImageMapReduceTask(config=config)
1146  results = task.run(scienceExposure, template=templateExposure, inImageSpace=inImageSpace,
1147  doScorr=doPreConvolve, forceEvenSized=False)
1148  results.D = results.exposure
1149  # The CoaddPsf, when used for detection does not utilize its spatially-varying
1150  # properties; it simply computes the PSF at its getAveragePosition().
1151  # TODO: we need to get it to return the matchedExposure (convolved template)
1152  # too, for dipole fitting; but the imageMapReduce task might need to be engineered
1153  # for this purpose.
1154  else:
1155  config = self.config.zogyConfig
1156  task = ZogyTask(scienceExposure=scienceExposure, templateExposure=templateExposure,
1157  config=config)
1158  if not doPreConvolve:
1159  results = task.computeDiffim(inImageSpace=inImageSpace)
1160  results.matchedExposure = results.R
1161  else:
1162  results = task.computeScorr(inImageSpace=inImageSpace)
1163  results.D = results.S
1164 
1165  # Make sure masks of input images are propagated to diffim
1166  mask = results.D.getMaskedImage().getMask()
1167  mask |= scienceExposure.getMaskedImage().getMask()
1168  mask |= templateExposure.getMaskedImage().getMask()
1169  results.D.getMaskedImage().getMask()[:, :] = mask
1170  badBitsNan = mask.addMaskPlane('UNMASKEDNAN')
1171  resultsArr = results.D.getMaskedImage().getMask().getArray()
1172  resultsArr[np.isnan(resultsArr)] |= badBitsNan
1173  resultsArr[np.isnan(scienceExposure.getMaskedImage().getImage().getArray())] |= badBitsNan
1174  resultsArr[np.isnan(templateExposure.getMaskedImage().getImage().getArray())] |= badBitsNan
1175 
1176  results.subtractedExposure = results.D
1177  results.warpedExposure = templateExposure
1178  return results
1179 
1180  def subtractMaskedImages(self, templateExposure, scienceExposure,
1181  doWarping=True, spatiallyVarying=True, inImageSpace=False,
1182  doPreConvolve=False):
1183  raise NotImplementedError
1184 
1185 
1186 subtractAlgorithmRegistry.register('zogy', ZogyImagePsfMatchTask)
def _computeImageMean(self, exposure)
Definition: zogy.py:1055
def _computeVarAstGradients(self, xVarAst=0., yVarAst=0., inImageSpace=False, R_hat=None, Kr_hat=None, Kr=None, N_hat=None, Kn_hat=None, Kn=None)
Definition: zogy.py:645
Configuration for image-to-image Psf matching.
def computeDiffim(self, inImageSpace=None, padSize=None, returnMatchedTemplate=False, kwargs)
Definition: zogy.py:562
def computeDiffimFourierSpace(self, debug=False, returnMatchedTemplate=False, kwargs)
Definition: zogy.py:378
def __init__(self, args, kwargs)
Definition: zogy.py:881
Psf-match two MaskedImages or Exposures using the sources in the images.
def computePrereqs(self, psf1=None, psf2=None, padSize=0)
Definition: zogy.py:328
def _padPsfToSize(psf, size)
Definition: zogy.py:307
def run(self, subExposure, expandedSubExposure, fullBBox, template, kwargs)
Definition: zogy.py:885
def subtractExposures(self, templateExposure, scienceExposure, doWarping=True, spatiallyVarying=True, inImageSpace=False, doPreConvolve=False)
Definition: zogy.py:1072
def setup(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None, psf1=None, psf2=None, correctBackground=False, args, kwargs)
Definition: zogy.py:170
def __init__(self, templateExposure=None, scienceExposure=None, sig1=None, sig2=None, psf1=None, psf2=None, args, kwargs)
Definition: zogy.py:135
def computeScorrFourierSpace(self, xVarAst=0., yVarAst=0., kwargs)
Definition: zogy.py:698
def _computeVarianceMean(self, exposure)
Definition: zogy.py:297
def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None)
Definition: zogy.py:612
def _doConvolve(self, exposure, kernel, recenterKernel=False)
Convolve an Exposure with a decorrelation convolution kernel.
Definition: zogy.py:457
def _setNewPsf(self, exposure, psfArr)
Definition: zogy.py:548
def subtractMaskedImages(self, templateExposure, scienceExposure, doWarping=True, spatiallyVarying=True, inImageSpace=False, doPreConvolve=False)
Definition: zogy.py:1182
def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, kwargs)
Definition: zogy.py:838
def computeDiffimImageSpace(self, padSize=None, debug=False, kwargs)
Definition: zogy.py:493
def _validateWcs(self, templateExposure, scienceExposure)
Return True if the WCS of the two Exposures have the same origin and extent.
def computeScorrImageSpace(self, xVarAst=0., yVarAst=0., padSize=None, kwargs)
Definition: zogy.py:776
def __init__(self, args, kwargs)
Definition: zogy.py:1052