lsst.ip.diffim  15.0-2-g5218728+1
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 
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 AttributeError:
488  # Allow exposure to actually be an image/maskedImage
489  # (getMaskedImage will throw AttributeError in that case)
490  afwMath.convolve(outExp, exposure, kern, convCntrl)
491 
492  return outExp, kern
493 
494  def computeDiffimImageSpace(self, padSize=None, debug=False, **kwargs):
495  """Compute ZOGY diffim `D` using image-space convlutions
496 
497  This method is still being debugged as it results in artifacts
498  when the PSFs are noisy (see module-level docstring). Thus
499  there are several options still enabled by the `debug` flag,
500  which are disabled by defult.
501 
502  Parameters
503  ----------
504  padSize : `int`, the amount to pad the PSFs by
505  debug : `bool`, flag to enable debugging tests and options
506 
507  Returns
508  -------
509  D : `lsst.afw.Exposure`
510  the proper image difference, including correct variance,
511  masks, and PSF
512  """
513  preqs = self.computePrereqs(padSize=padSize)
514 
515  delta = 0.
516  if debug:
517  delta = 1. # Regularize the ratio, a possible option to remove artifacts
518  Kr_hat = (preqs.Pr_hat + delta) / (preqs.denom + delta)
519  Kn_hat = (preqs.Pn_hat + delta) / (preqs.denom + delta)
520  Kr = np.fft.ifft2(Kr_hat).real
521  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
522  Kn = np.fft.ifft2(Kn_hat).real
523  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
524 
525  def _trimKernel(self, K, trim_amount):
526  # Trim out the wings of Kn, Kr (see notebook #15)
527  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
528  ps = trim_amount
529  K = K[ps:-ps, ps:-ps]
530  return K
531 
532  padSize = self.padSize if padSize is None else padSize
533  # Enabling this block (debug=True) makes it slightly faster, but ~25% worse artifacts:
534  if debug and self.config.doTrimKernels: # default False
535  # Filtering also makes it slightly faster (zeros are ignored in convolution)
536  # but a bit worse. Filter the wings of Kn, Kr (see notebook #15)
537  Kn = _trimKernel(Kn, padSize)
538  Kr = _trimKernel(Kr, padSize)
539 
540  # Note these are reverse-labelled, this is CORRECT!
541  exp1, _ = self._doConvolve(self.template, Kn)
542  exp2, _ = self._doConvolve(self.science, Kr)
543  D = exp2
544  tmp = D.getMaskedImage()
545  tmp -= exp1.getMaskedImage()
546  tmp /= preqs.Fd
547  return pipeBase.Struct(D=D, R=exp1)
548 
549  def _setNewPsf(self, exposure, psfArr):
550  """Utility method to set an exposure's PSF when provided as a 2-d numpy.array
551  """
552  bbox = exposure.getBBox()
553  center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.)
554  center = afwGeom.Point2D(center[0], center[1])
555  psfI = afwImage.ImageD(psfArr.shape[1], psfArr.shape[0])
556  psfI.getArray()[:, :] = psfArr
557  psfK = afwMath.FixedKernel(psfI)
558  psfNew = measAlg.KernelPsf(psfK, center)
559  exposure.setPsf(psfNew)
560  return exposure
561 
562  def computeDiffim(self, inImageSpace=None, padSize=None,
563  returnMatchedTemplate=False, **kwargs):
564  """Wrapper method to compute ZOGY proper diffim
565 
566  This method should be used as the public interface for
567  computing the ZOGY diffim.
568 
569  Parameters
570  ----------
571  inImageSpace : `bool`
572  Override config `inImageSpace` parameter
573  padSize : `int`
574  Override config `padSize` parameter
575  returnMatchedTemplate : bool
576  Include the PSF-matched template in the results Struct
577  **kwargs : `dict`
578  additional keyword arguments to be passed to
579  `computeDiffimFourierSpace` or `computeDiffimImageSpace`.
580 
581  Returns
582  -------
583  An lsst.pipe.base.Struct containing:
584  - D : `lsst.afw.Exposure`
585  the proper image difference, including correct variance,
586  masks, and PSF
587  - R : `lsst.afw.Exposure`
588  If `returnMatchedTemplate` is True, the PSF-matched template
589  exposure
590  """
591  R = None
592  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
593  if inImageSpace:
594  padSize = self.padSize if padSize is None else padSize
595  res = self.computeDiffimImageSpace(padSize=padSize, **kwargs)
596  D = res.D
597  if returnMatchedTemplate:
598  R = res.R
599  else:
600  res = self.computeDiffimFourierSpace(**kwargs)
601  D = self.science.clone()
602  D.getMaskedImage().getImage().getArray()[:, :] = res.D
603  D.getMaskedImage().getVariance().getArray()[:, :] = res.D_var
604  if returnMatchedTemplate:
605  R = self.science.clone()
606  R.getMaskedImage().getImage().getArray()[:, :] = res.R
607  R.getMaskedImage().getVariance().getArray()[:, :] = res.R_var
608 
609  psf = self.computeDiffimPsf()
610  D = self._setNewPsf(D, psf)
611  return pipeBase.Struct(D=D, R=R)
612 
613  def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None):
614  """Compute the ZOGY diffim PSF (ZOGY manuscript eq. 14)
615 
616  Parameters
617  ----------
618  padSize : `int`
619  Override config `padSize` parameter
620  keepFourier : `bool`
621  Return the FFT of the diffim PSF (do not inverse-FFT it)
622  psf1 : 2D `numpy.array`
623  (Optional) Input psf of template, override if already padded
624  psf2 : 2D `numpy.array`
625  (Optional) Input psf of science image, override if already padded
626 
627  Returns
628  -------
629  Pd : 2D `numpy.array`, the diffim PSF (or FFT of PSF if `keepFourier=True`)
630  """
631  preqs = self.computePrereqs(psf1=psf1, psf2=psf2, padSize=padSize)
632 
633  Pd_hat_numerator = (self.Fr * self.Fn * preqs.Pr_hat * preqs.Pn_hat)
634  Pd_hat = Pd_hat_numerator / (preqs.Fd * preqs.denom)
635 
636  if keepFourier:
637  return Pd_hat
638 
639  Pd = np.fft.ifft2(Pd_hat)
640  Pd = np.fft.ifftshift(Pd).real
641 
642  return Pd
643 
644  def _computeVarAstGradients(self, xVarAst=0., yVarAst=0., inImageSpace=False,
645  R_hat=None, Kr_hat=None, Kr=None,
646  N_hat=None, Kn_hat=None, Kn=None):
647  """Compute the astrometric noise correction terms
648 
649  Compute the correction for estimated astrometric noise as
650  proscribed in ZOGY (2016), section 3.3. All convolutions
651  performed either in real (image) or Fourier space.
652 
653  Parameters
654  ----------
655  xVarAst, yVarAst : `float`
656  estimated astrometric noise (variance of astrometric registration errors)
657  inImageSpace : `bool`
658  Perform all convolutions in real (image) space rather than Fourier space
659  R_hat : 2-D `numpy.array`
660  (Optional) FFT of template image, only required if `inImageSpace=False`
661  Kr_hat : 2-D `numpy.array`
662  FFT of Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=False`
663  Kr : 2-D `numpy.array`
664  Kr kernel (eq. 28 of ZOGY (2016)), only required if `inImageSpace=True`.
665  Kr is associated with the template (reference).
666  N_hat : 2-D `numpy.array`
667  FFT of science image, only required if `inImageSpace=False`
668  Kn_hat : 2-D `numpy.array`
669  FFT of Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=False`
670  Kn : 2-D `numpy.array`
671  Kn kernel (eq. 29 of ZOGY (2016)), only required if `inImageSpace=True`.
672  Kn is associated with the science (new) image.
673 
674  Returns
675  -------
676  VastSR, VastSN : 2-D `numpy.array`s containing the values in eqs. 30 and 32 of
677  ZOGY (2016).
678  """
679  VastSR = VastSN = 0.
680  if xVarAst + yVarAst > 0: # Do the astrometric variance correction
681  if inImageSpace:
682  S_R, _ = self._doConvolve(self.template, Kr)
683  S_R = S_R.getMaskedImage().getImage().getArray()
684  else:
685  S_R = np.fft.ifft2(R_hat * Kr_hat)
686  gradRx, gradRy = np.gradient(S_R)
687  VastSR = xVarAst * gradRx**2. + yVarAst * gradRy**2.
688 
689  if inImageSpace:
690  S_N, _ = self._doConvolve(self.science, Kn)
691  S_N = S_N.getMaskedImage().getImage().getArray()
692  else:
693  S_N = np.fft.ifft2(N_hat * Kn_hat)
694  gradNx, gradNy = np.gradient(S_N)
695  VastSN = xVarAst * gradNx**2. + yVarAst * gradNy**2.
696 
697  return VastSR, VastSN
698 
699  def computeScorrFourierSpace(self, xVarAst=0., yVarAst=0., **kwargs):
700  """Compute corrected likelihood image, optimal for source detection
701 
702  Compute ZOGY S_corr image. This image can be thresholded for
703  detection without optimal filtering, and the variance image is
704  corrected to account for astrometric noise (errors in
705  astrometric registration whether systematic or due to effects
706  such as DCR). The calculations here are all performed in
707  Fourier space, as proscribed in ZOGY (2016).
708 
709  Parameters
710  ----------
711  xVarAst, yVarAst : `float`
712  estimated astrometric noise (variance of astrometric registration errors)
713 
714  Returns
715  -------
716  A lsst.pipe.base.Struct containing:
717  - S : `numpy.array`, the likelihood image S (eq. 12 of ZOGY (2016))
718  - S_var : the corrected variance image (denominator of eq. 25 of ZOGY (2016))
719  - Dpsf : the PSF of the diffim D, likely never to be used.
720  """
721  # Some masked regions are NaN or infinite!, and FFTs no likey.
722  def fix_nans(im):
723  """Replace any NaNs or Infs with the mean of the image."""
724  isbad = ~np.isfinite(im)
725  if np.any(isbad):
726  im[isbad] = np.nan
727  im[isbad] = np.nanmean(im)
728  return im
729 
730  self.im1 = fix_nans(self.im1)
731  self.im2 = fix_nans(self.im2)
732  self.im1_var = fix_nans(self.im1_var)
733  self.im2_var = fix_nans(self.im2_var)
734 
735  # Do all in fourier space (needs image-sized PSFs)
736  psf1 = ZogyTask._padPsfToSize(self.im1_psf, self.im1.shape)
737  psf2 = ZogyTask._padPsfToSize(self.im2_psf, self.im2.shape)
738 
739  preqs = self.computePrereqs(psf1, psf2, padSize=0) # already padded the PSFs
740 
741  # Compute D_hat here (don't need D then, for speed)
742  R_hat = np.fft.fft2(self.im1)
743  N_hat = np.fft.fft2(self.im2)
744  D_hat = self.Fr * preqs.Pr_hat * N_hat - self.Fn * preqs.Pn_hat * R_hat
745  D_hat /= preqs.denom
746 
747  Pd_hat = self.computeDiffimPsf(padSize=0, keepFourier=True, psf1=psf1, psf2=psf2)
748  Pd_bar = np.conj(Pd_hat)
749  S = np.fft.ifft2(D_hat * Pd_bar)
750 
751  # Adjust the variance planes of the two images to contribute to the final detection
752  # (eq's 26-29).
753  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat
754  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
755  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat
756  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
757 
758  Kr_hat2 = np.fft.fft2(np.fft.ifft2(Kr_hat)**2.)
759  Kn_hat2 = np.fft.fft2(np.fft.ifft2(Kn_hat)**2.)
760  var1c_hat = Kr_hat2 * np.fft.fft2(self.im1_var)
761  var2c_hat = Kn_hat2 * np.fft.fft2(self.im2_var)
762 
763  # Do the astrometric variance correction
764  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=False,
765  R_hat=R_hat, Kr_hat=Kr_hat,
766  N_hat=N_hat, Kn_hat=Kn_hat)
767 
768  S_var = np.sqrt(np.fft.ifftshift(np.fft.ifft2(var1c_hat + var2c_hat)) + fGradR + fGradN)
769  S_var *= preqs.Fd
770 
771  S = np.fft.ifftshift(np.fft.ifft2(Kn_hat * N_hat - Kr_hat * R_hat))
772  S *= preqs.Fd
773 
774  Pd = self.computeDiffimPsf(padSize=0)
775  return pipeBase.Struct(S=S.real, S_var=S_var.real, Dpsf=Pd)
776 
777  def computeScorrImageSpace(self, xVarAst=0., yVarAst=0., padSize=None, **kwargs):
778  """Compute corrected likelihood image, optimal for source detection
779 
780  Compute ZOGY S_corr image. This image can be thresholded for
781  detection without optimal filtering, and the variance image is
782  corrected to account for astrometric noise (errors in
783  astrometric registration whether systematic or due to effects
784  such as DCR). The calculations here are all performed in
785  Real (image) space.
786 
787  Parameters
788  ----------
789  xVarAst, yVarAst : `float`
790  estimated astrometric noise (variance of astrometric registration errors)
791 
792  Returns
793  -------
794  a tuple containing:
795  - S : `lsst.afw.image.Exposure`, the likelihood exposure S (eq. 12 of ZOGY (2016)),
796  including corrected variance, masks, and PSF
797  - D : `lsst.afw.image.Exposure`, the proper image difference, including correct
798  variance, masks, and PSF
799  """
800  # Do convolutions in image space
801  preqs = self.computePrereqs(padSize=0)
802 
803  padSize = self.padSize if padSize is None else padSize
804  D = self.computeDiffimImageSpace(padSize=padSize).D
805  Pd = self.computeDiffimPsf()
806  D = self._setNewPsf(D, Pd)
807  Pd_bar = np.fliplr(np.flipud(Pd))
808  S, _ = self._doConvolve(D, Pd_bar)
809  tmp = S.getMaskedImage()
810  tmp *= preqs.Fd
811 
812  # Adjust the variance planes of the two images to contribute to the final detection
813  # (eq's 26-29).
814  Pn_hat2 = np.conj(preqs.Pn_hat) * preqs.Pn_hat
815  Kr_hat = self.Fr * self.Fn**2. * np.conj(preqs.Pr_hat) * Pn_hat2 / preqs.denom**2.
816  Pr_hat2 = np.conj(preqs.Pr_hat) * preqs.Pr_hat
817  Kn_hat = self.Fn * self.Fr**2. * np.conj(preqs.Pn_hat) * Pr_hat2 / preqs.denom**2.
818 
819  Kr = np.fft.ifft2(Kr_hat).real
820  Kr = np.roll(np.roll(Kr, -1, 0), -1, 1)
821  Kn = np.fft.ifft2(Kn_hat).real
822  Kn = np.roll(np.roll(Kn, -1, 0), -1, 1)
823  var1c, _ = self._doConvolve(self.template.getMaskedImage().getVariance(), Kr**2.)
824  var2c, _ = self._doConvolve(self.science.getMaskedImage().getVariance(), Kn**2.)
825 
826  # Do the astrometric variance correction
827  fGradR, fGradN = self._computeVarAstGradients(xVarAst, yVarAst, inImageSpace=True,
828  Kr=Kr, Kn=Kn)
829 
830  Smi = S.getMaskedImage()
831  Smi *= preqs.Fd
832  S_var = np.sqrt(var1c.getArray() + var2c.getArray() + fGradR + fGradN)
833  S.getMaskedImage().getVariance().getArray()[:, :] = S_var
834  S = self._setNewPsf(S, Pd)
835 
836  # also return diffim since it was calculated and might be desired
837  return pipeBase.Struct(S=S, D=D)
838 
839  def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, **kwargs):
840  """Wrapper method to compute ZOGY corrected likelihood image, optimal for
841  source detection
842 
843  This method should be used as the public interface for
844  computing the ZOGY S_corr.
845 
846  Parameters
847  ----------
848  xVarAst, yVarAst : float
849  estimated astrometric noise (variance of astrometric registration errors)
850  inImageSpace : bool
851  Override config `inImageSpace` parameter
852  padSize : int
853  Override config `padSize` parameter
854 
855  Returns
856  -------
857  S : lsst.afw.image.Exposure, the likelihood exposure S (eq. 12 of ZOGY (2016)),
858  including corrected variance, masks, and PSF
859  """
860  inImageSpace = self.config.inImageSpace if inImageSpace is None else inImageSpace
861  if inImageSpace:
862  res = self.computeScorrImageSpace(xVarAst=xVarAst, yVarAst=yVarAst, padSize=padSize)
863  S = res.S
864  else:
865  res = self.computeScorrFourierSpace(xVarAst=xVarAst, yVarAst=yVarAst)
866 
867  S = self.science.clone()
868  S.getMaskedImage().getImage().getArray()[:, :] = res.S
869  S.getMaskedImage().getVariance().getArray()[:, :] = res.S_var
870  S = self._setNewPsf(S, res.Dpsf)
871 
872  return pipeBase.Struct(S=S)
873 
874 
876  """Task to be used as an ImageMapper for performing
877  ZOGY image subtraction on a grid of subimages.
878  """
879  ConfigClass = ZogyConfig
880  _DefaultName = 'ip_diffim_ZogyMapper'
881 
882  def __init__(self, *args, **kwargs):
883  ImageMapper.__init__(self, *args, **kwargs)
884 
885  def run(self, subExposure, expandedSubExposure, fullBBox, template,
886  **kwargs):
887  """Perform ZOGY proper image subtraction on sub-images
888 
889  This method performs ZOGY proper image subtraction on
890  `subExposure` using local measures for image variances and
891  PSF. `subExposure` is a sub-exposure of the science image. It
892  also requires the corresponding sub-exposures of the template
893  (`template`). The operations are actually performed on
894  `expandedSubExposure` to allow for invalid edge pixels arising
895  from convolutions, which are then removed.
896 
897  Parameters
898  ----------
899  subExposure : lsst.afw.image.Exposure
900  the sub-exposure of the diffim
901  expandedSubExposure : lsst.afw.image.Exposure
902  the expanded sub-exposure upon which to operate
903  fullBBox : lsst.afw.geom.BoundingBox
904  the bounding box of the original exposure
905  template : lsst.afw.image.Exposure
906  the template exposure, from which a corresponding sub-exposure
907  is extracted
908  kwargs :
909  additional keyword arguments propagated from
910  `ImageMapReduceTask.run`. These include:
911  - doScorr : bool
912  Compute and return the corrected likelihood image S_corr
913  rather than the proper image difference
914  - inImageSpace : bool
915  Perform all convolutions in real (image) space rather than
916  in Fourier space. This option currently leads to artifacts
917  when using real (measured and noisy) PSFs, thus it is set
918  to `False` by default.
919  These kwargs may also include arguments to be propagated to
920  `ZogyTask.computeDiffim` and `ZogyTask.computeScorr`.
921 
922  Returns
923  -------
924  A `lsst.pipe.base.Struct` containing the result of the
925  `subExposure` processing, labelled 'subExposure'. In this case
926  it is either the subExposure of the proper image difference D,
927  or (if `doScorr==True`) the corrected likelihood exposure `S`.
928 
929  Notes
930  -----
931  This `run` method accepts parameters identical to those of
932  `ImageMapper.run`, since it is called from the
933  `ImageMapperTask`. See that class for more information.
934  """
935  bbox = subExposure.getBBox()
936  center = ((bbox.getBeginX() + bbox.getEndX()) // 2., (bbox.getBeginY() + bbox.getEndY()) // 2.)
937  center = afwGeom.Point2D(center[0], center[1])
938 
939  imageSpace = kwargs.pop('inImageSpace', False)
940  doScorr = kwargs.pop('doScorr', False)
941  sigmas = kwargs.pop('sigmas', None)
942  padSize = kwargs.pop('padSize', 7)
943 
944  # Psf and image for science img (index 2)
945  subExp2 = expandedSubExposure
946 
947  # Psf and image for template img (index 1)
948  subExp1 = template.Factory(template, expandedSubExposure.getBBox())
949 
950  if sigmas is None:
951  sig1 = sig2 = None
952  else:
953  # for testing, can use the input sigma (e.g., global value for entire exposure)
954  sig1, sig2 = sigmas[0], sigmas[1]
955 
956  def _makePsfSquare(psf):
957  # Sometimes CoaddPsf does this. Make it square.
958  if psf.shape[0] < psf.shape[1]:
959  psf = np.pad(psf, ((0, psf.shape[1] - psf.shape[0]), (0, 0)), mode='constant',
960  constant_values=0.)
961  elif psf.shape[0] > psf.shape[1]:
962  psf = np.pad(psf, ((0, 0), (0, psf.shape[0] - psf.shape[1])), mode='constant',
963  constant_values=0.)
964  return psf
965 
966  psf2 = subExp2.getPsf().computeKernelImage(center).getArray()
967  psf2 = _makePsfSquare(psf2)
968 
969  psf1 = template.getPsf().computeKernelImage(center).getArray()
970  psf1 = _makePsfSquare(psf1)
971 
972  # from diffimTests.diffimTests ...
973  if subExp1.getDimensions()[0] < psf1.shape[0] or subExp1.getDimensions()[1] < psf1.shape[1]:
974  return pipeBase.Struct(subExposure=subExposure)
975 
976  def _filterPsf(psf):
977  """Filter a noisy Psf to remove artifacts. Subject of future research."""
978  # only necessary if it's from a measured psf and PsfEx seems to always make PSFs of size 41x41
979  if psf.shape[0] == 41: # its from a measured psf
980  psf = psf.copy()
981  psf[psf < 0] = 0
982  psf[0:10, :] = psf[:, 0:10] = psf[31:41, :] = psf[:, 31:41] = 0
983  psf /= psf.sum()
984 
985  return psf
986 
987  psf1b = psf2b = None
988  if self.config.doFilterPsfs: # default True
989  # Note this *really* helps for measured psfs.
990  psf1b = _filterPsf(psf1)
991  psf2b = _filterPsf(psf2)
992 
993  config = ZogyConfig()
994  if imageSpace is True:
995  config.inImageSpace = imageSpace
996  config.padSize = padSize # Don't need padding if doing all in fourier space
997  task = ZogyTask(templateExposure=subExp1, scienceExposure=subExp2,
998  sig1=sig1, sig2=sig2, psf1=psf1b, psf2=psf2b, config=config)
999 
1000  if not doScorr:
1001  res = task.computeDiffim(**kwargs)
1002  D = res.D
1003  else:
1004  res = task.computeScorr(**kwargs)
1005  D = res.S
1006 
1007  outExp = D.Factory(D, subExposure.getBBox())
1008  out = pipeBase.Struct(subExposure=outExp)
1009  return out
1010 
1011 
1013  """Config to be passed to ImageMapReduceTask
1014 
1015  This config targets the imageMapper to use the ZogyMapper.
1016  """
1017  mapper = pexConfig.ConfigurableField(
1018  doc='Zogy task to run on each sub-image',
1019  target=ZogyMapper
1020  )
1021 
1022 
1024  """Config for the ZogyImagePsfMatchTask"""
1025 
1026  zogyConfig = pexConfig.ConfigField(
1027  dtype=ZogyConfig,
1028  doc='ZogyTask config to use when running on complete exposure (non spatially-varying)',
1029  )
1030 
1031  zogyMapReduceConfig = pexConfig.ConfigField(
1032  dtype=ZogyMapReduceConfig,
1033  doc='ZogyMapReduce config to use when running Zogy on each sub-image (spatially-varying)',
1034  )
1035 
1036  def setDefaults(self):
1037  self.zogyMapReduceConfig.gridStepX = self.zogyMapReduceConfig.gridStepY = 40
1038  self.zogyMapReduceConfig.cellSizeX = self.zogyMapReduceConfig.cellSizeY = 41
1039  self.zogyMapReduceConfig.borderSizeX = self.zogyMapReduceConfig.borderSizeY = 8
1040  self.zogyMapReduceConfig.reducer.reduceOperation = 'average'
1041  self.zogyConfig.inImageSpace = False
1042 
1043 
1045  """Task to perform Zogy PSF matching and image subtraction.
1046 
1047  This class inherits from ImagePsfMatchTask to contain the _warper
1048  subtask and related methods.
1049  """
1050 
1051  ConfigClass = ZogyImagePsfMatchConfig
1052 
1053  def __init__(self, *args, **kwargs):
1054  ImagePsfMatchTask.__init__(self, *args, **kwargs)
1055 
1056  def _computeImageMean(self, exposure):
1057  """Compute the sigma-clipped mean of the pixels image of `exposure`.
1058  """
1059  statsControl = afwMath.StatisticsControl()
1060  statsControl.setNumSigmaClip(3.)
1061  statsControl.setNumIter(3)
1062  ignoreMaskPlanes = ("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE")
1063  statsControl.setAndMask(afwImage.Mask.getPlaneBitMask(ignoreMaskPlanes))
1064  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getImage(),
1065  exposure.getMaskedImage().getMask(),
1066  afwMath.MEANCLIP | afwMath.MEDIAN, statsControl)
1067  mn = statObj.getValue(afwMath.MEANCLIP)
1068  med = statObj.getValue(afwMath.MEDIAN)
1069  return mn, med
1070 
1071  def subtractExposures(self, templateExposure, scienceExposure,
1072  doWarping=True, spatiallyVarying=True, inImageSpace=False,
1073  doPreConvolve=False):
1074  """Register, PSF-match, and subtract two Exposures using the ZOGY algorithm.
1075 
1076  Do the following, in order:
1077  - Warp templateExposure to match scienceExposure, if their WCSs do not already match
1078  - Compute subtracted exposure ZOGY image subtraction algorithm on the two exposures
1079 
1080  Parameters
1081  ----------
1082  templateExposure : `lsst.afw.image.Exposure`
1083  exposure to PSF-match to scienceExposure. The exposure's mean value is subtracted
1084  in-place.
1085  scienceExposure : `lsst.afw.image.Exposure`
1086  reference Exposure. The exposure's mean value is subtracted in-place.
1087  doWarping : `bool`
1088  what to do if templateExposure's and scienceExposure's WCSs do not match:
1089  - if True then warp templateExposure to match scienceExposure
1090  - if False then raise an Exception
1091  spatiallyVarying : bool
1092  If True, perform the operation over a grid of patches across the two exposures
1093  inImageSpace : `bool`
1094  If True, perform the Zogy convolutions in image space rather than in frequency space.
1095  doPreConvolve : `bool`
1096  ***Currently not implemented.*** If True assume we are to compute the match filter-convolved
1097  exposure which can be thresholded for detection. In the case of Zogy this would mean
1098  we compute the Scorr image.
1099 
1100  Returns
1101  -------
1102  A `lsst.pipe.base.Struct` containing these fields:
1103  - subtractedExposure: subtracted Exposure
1104  - warpedExposure: templateExposure after warping to match scienceExposure (if doWarping true)
1105  """
1106 
1107  mn1 = self._computeImageMean(templateExposure)
1108  mn2 = self._computeImageMean(scienceExposure)
1109  self.log.info("Exposure means=%f, %f; median=%f, %f:" % (mn1[0], mn2[0], mn1[1], mn2[1]))
1110  if not np.isnan(mn1[0]) and np.abs(mn1[0]) > 1:
1111  mi = templateExposure.getMaskedImage()
1112  mi -= mn1[0]
1113  if not np.isnan(mn2[0]) and np.abs(mn2[0]) > 1:
1114  mi = scienceExposure.getMaskedImage()
1115  mi -= mn2[0]
1116 
1117  self.log.info('Running Zogy algorithm: spatiallyVarying=%r' % spatiallyVarying)
1118 
1119  if not self._validateWcs(templateExposure, scienceExposure):
1120  if doWarping:
1121  self.log.info("Astrometrically registering template to science image")
1122  # Also warp the PSF
1123  xyTransform = afwGeom.makeWcsPairTransform(templateExposure.getWcs(),
1124  scienceExposure.getWcs())
1125  psfWarped = measAlg.WarpedPsf(templateExposure.getPsf(), xyTransform)
1126  templateExposure = self._warper.warpExposure(scienceExposure.getWcs(),
1127  templateExposure,
1128  destBBox=scienceExposure.getBBox())
1129 
1130  templateExposure.setPsf(psfWarped)
1131  else:
1132  self.log.error("ERROR: Input images not registered")
1133  raise RuntimeError("Input images not registered")
1134 
1135  def gm(exp):
1136  return exp.getMaskedImage().getMask()
1137 
1138  def ga(exp):
1139  return exp.getMaskedImage().getImage().getArray()
1140 
1141  if self.config.zogyConfig.inImageSpace:
1142  inImageSpace = True # Override
1143  self.log.info('Running Zogy algorithm: inImageSpace=%r' % inImageSpace)
1144  if spatiallyVarying:
1145  config = self.config.zogyMapReduceConfig
1146  task = ImageMapReduceTask(config=config)
1147  results = task.run(scienceExposure, template=templateExposure, inImageSpace=inImageSpace,
1148  doScorr=doPreConvolve, forceEvenSized=False)
1149  results.D = results.exposure
1150  # The CoaddPsf, when used for detection does not utilize its spatially-varying
1151  # properties; it simply computes the PSF at its getAveragePosition().
1152  # TODO: we need to get it to return the matchedExposure (convolved template)
1153  # too, for dipole fitting; but the imageMapReduce task might need to be engineered
1154  # for this purpose.
1155  else:
1156  config = self.config.zogyConfig
1157  task = ZogyTask(scienceExposure=scienceExposure, templateExposure=templateExposure,
1158  config=config)
1159  if not doPreConvolve:
1160  results = task.computeDiffim(inImageSpace=inImageSpace)
1161  results.matchedExposure = results.R
1162  else:
1163  results = task.computeScorr(inImageSpace=inImageSpace)
1164  results.D = results.S
1165 
1166  # Make sure masks of input images are propagated to diffim
1167  mask = results.D.getMaskedImage().getMask()
1168  mask |= scienceExposure.getMaskedImage().getMask()
1169  mask |= templateExposure.getMaskedImage().getMask()
1170  results.D.getMaskedImage().getMask()[:, :] = mask
1171  badBitsNan = mask.addMaskPlane('UNMASKEDNAN')
1172  resultsArr = results.D.getMaskedImage().getMask().getArray()
1173  resultsArr[np.isnan(resultsArr)] |= badBitsNan
1174  resultsArr[np.isnan(scienceExposure.getMaskedImage().getImage().getArray())] |= badBitsNan
1175  resultsArr[np.isnan(templateExposure.getMaskedImage().getImage().getArray())] |= badBitsNan
1176 
1177  results.subtractedExposure = results.D
1178  results.warpedExposure = templateExposure
1179  return results
1180 
1181  def subtractMaskedImages(self, templateExposure, scienceExposure,
1182  doWarping=True, spatiallyVarying=True, inImageSpace=False,
1183  doPreConvolve=False):
1184  raise NotImplementedError
1185 
1186 
1187 subtractAlgorithmRegistry.register('zogy', ZogyImagePsfMatchTask)
def _computeImageMean(self, exposure)
Definition: zogy.py:1056
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:646
Configuration for image-to-image Psf matching.
def computeDiffim(self, inImageSpace=None, padSize=None, returnMatchedTemplate=False, kwargs)
Definition: zogy.py:563
def computeDiffimFourierSpace(self, debug=False, returnMatchedTemplate=False, kwargs)
Definition: zogy.py:378
def __init__(self, args, kwargs)
Definition: zogy.py:882
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
Statistics makeStatistics(lsst::afw::math::MaskedVector< EntryT > const &mv, std::vector< WeightPixel > const &vweights, int const flags, StatisticsControl const &sctrl=StatisticsControl())
def run(self, subExposure, expandedSubExposure, fullBBox, template, kwargs)
Definition: zogy.py:886
def subtractExposures(self, templateExposure, scienceExposure, doWarping=True, spatiallyVarying=True, inImageSpace=False, doPreConvolve=False)
Definition: zogy.py:1073
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:699
std::shared_ptr< TransformPoint2ToPoint2 > makeWcsPairTransform(SkyWcs const &src, SkyWcs const &dst)
def _computeVarianceMean(self, exposure)
Definition: zogy.py:297
def computeDiffimPsf(self, padSize=0, keepFourier=False, psf1=None, psf2=None)
Definition: zogy.py:613
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:549
def subtractMaskedImages(self, templateExposure, scienceExposure, doWarping=True, spatiallyVarying=True, inImageSpace=False, doPreConvolve=False)
Definition: zogy.py:1183
def computeScorr(self, xVarAst=0., yVarAst=0., inImageSpace=None, padSize=0, kwargs)
Definition: zogy.py:839
def computeDiffimImageSpace(self, padSize=None, debug=False, kwargs)
Definition: zogy.py:494
void convolve(OutImageT &convolvedImage, InImageT const &inImage, KernelT const &kernel, bool doNormalize, bool doCopyEdge=false)
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:777
def __init__(self, args, kwargs)
Definition: zogy.py:1053