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