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