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