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