lsst.ip.diffim  21.0.0-15-g8b2b7d1+05ea3e5ea0
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 from lsst.geom import Box2I, Point2I, Extent2I
26 import lsst.afw.geom as afwGeom
27 import lsst.afw.image as afwImage
28 from lsst.afw.image import ImageOrigin
29 import lsst.afw.table as afwTable
30 import lsst.afw.math as afwMath
31 import lsst.meas.algorithms as measAlg
32 import lsst.pipe.base as pipeBase
33 import lsst.pex.config as pexConfig
34 
35 from .imagePsfMatch import (ImagePsfMatchTask, ImagePsfMatchConfig,
36  subtractAlgorithmRegistry)
37 
38 __all__ = ["ZogyTask", "ZogyConfig",
39  "ZogyImagePsfMatchConfig", "ZogyImagePsfMatchTask"]
40 
41 
42 """Tasks for performing the "Proper image subtraction" algorithm of
43 Zackay, et al. (2016), hereafter simply referred to as 'ZOGY (2016)'.
44 
45 `ZogyTask` contains methods to perform the basic estimation of the
46 ZOGY diffim ``D``, its updated PSF, and the variance-normalized
47 likelihood image ``S_corr``. We have implemented ZOGY using the
48 proscribed methodology, computing all convolutions in Fourier space,
49 and also variants in which the convolutions are performed in real
50 (image) space. The former is faster and results in fewer artifacts
51 when the PSFs are noisy (i.e., measured, for example, via
52 `PsfEx`). The latter is presumed to be preferred as it can account for
53 masks correctly with fewer "ringing" artifacts from edge effects or
54 saturated stars, but noisy PSFs result in their own smaller
55 artifacts. Removal of these artifacts is a subject of continuing
56 research. Currently, we "pad" the PSFs when performing the
57 subtractions in real space, which reduces, but does not entirely
58 eliminate these artifacts.
59 
60 All methods in `ZogyTask` assume template and science images are
61 already accurately photometrically and astrometrically registered.
62 
63 `ZogyMapper` is a wrapper which runs `ZogyTask` in the
64 `ImageMapReduce` framework, computing of ZOGY diffim's on small,
65 overlapping sub-images, thereby enabling complete ZOGY diffim's which
66 account for spatially-varying noise and PSFs across the two input
67 exposures. An example of the use of this task is in the `testZogy.py`
68 unit test.
69 """
70 
71 
72 class ZogyConfig(pexConfig.Config):
73  """Configuration parameters for the ZogyTask
74  """
75 
76  templateFluxScaling = pexConfig.Field(
77  dtype=float,
78  default=1.,
79  doc="Template flux scaling factor (Fr in ZOGY paper)"
80  )
81 
82  scienceFluxScaling = pexConfig.Field(
83  dtype=float,
84  default=1.,
85  doc="Science flux scaling factor (Fn in ZOGY paper)"
86  )
87 
88  scaleByCalibration = pexConfig.Field(
89  dtype=bool,
90  default=True,
91  doc="Compute the flux normalization scaling based on the image calibration."
92  "This overrides 'templateFluxScaling' and 'scienceFluxScaling'."
93  )
94 
95  correctBackground = pexConfig.Field(
96  dtype=bool,
97  default=False,
98  doc="Subtract exposure background mean to have zero expectation value."
99  )
100 
101  ignoreMaskPlanes = pexConfig.ListField(
102  dtype=str,
103  default=("INTRP", "EDGE", "DETECTED", "SAT", "CR", "BAD", "NO_DATA", "DETECTED_NEGATIVE"),
104  doc="Mask planes to ignore for statistics"
105  )
106  maxPsfCentroidDist = pexConfig.Field(
107  dtype=float,
108  default=0.2,
109  doc="Maximum centroid difference allowed between the two exposure PSFs (pixels)."
110  )
111  doSpatialGrid = pexConfig.Field(
112  dtype=bool,
113  default=False,
114  doc="Split the exposure and perform matching with the spatially varying PSF."
115  )
116  gridInnerSize = pexConfig.Field(
117  dtype=float,
118  default=8,
119  doc="Approximate useful inner size of the grid cells in units of the "
120  "estimated matching kernel size (doSpatialGrid=True only)."
121  )
122 
123 
124 class ZogyTask(pipeBase.Task):
125  """Task to perform ZOGY proper image subtraction. See module-level documentation for
126  additional details.
127 
128  """
129  ConfigClass = ZogyConfig
130  _DefaultName = "imageDifferenceZogy"
131 
132  def _computeVarianceMean(self, exposure):
133  """Compute the sigma-clipped mean of the variance image of ``exposure``.
134  """
135  statObj = afwMath.makeStatistics(exposure.getMaskedImage().getVariance(),
136  exposure.getMaskedImage().getMask(),
137  afwMath.MEANCLIP, self.statsControlstatsControl)
138  var = statObj.getValue(afwMath.MEANCLIP)
139  return var
140 
141  @staticmethod
142  def padCenterOriginArray(A, newShape, useInverse=False, dtype=None):
143  """Zero pad an image where the origin is at the center and replace the
144  origin to the corner as required by the periodic input of FFT.
145 
146  Implement also the inverse operation, crop the padding and re-center data.
147 
148  Parameters
149  ----------
150  A : `numpy.ndarray`
151  An array to copy from.
152  newShape : `tuple` of `int`
153  The dimensions of the resulting array. For padding, the resulting array
154  must be larger than A in each dimension. For the inverse operation this
155  must be the original, before padding dimensions of the array.
156  useInverse : bool, optional
157  Selector of forward, add padding, operation (False)
158  or its inverse, crop padding, operation (True).
159  dtype: `numpy.dtype`, optional
160  Dtype of output array. Values must be implicitly castable to this type.
161  Use to get expected result type, e.g. single float (nympy.float32).
162  If not specified, dtype is inherited from ``A``.
163 
164  Returns
165  -------
166  R : `numpy.ndarray`
167  The padded or unpadded array with shape of `newShape` and dtype of ``dtype``.
168 
169  Notes
170  -----
171  For odd dimensions, the splitting is rounded to
172  put the center pixel into the new corner origin (0,0). This is to be consistent
173  e.g. for a dirac delta kernel that is originally located at the center pixel.
174 
175 
176  Raises
177  ------
178  ValueError : ``newShape`` dimensions must be greater than or equal to the
179  dimensions of ``A`` for the forward operation and less than or equal to
180  for the inverse operation.
181  """
182 
183  # The forward and inverse operations should round odd dimension halves at the opposite
184  # sides to get the pixels back to their original positions.
185  if not useInverse:
186  # Forward operation: First and second halves with respect to the axes of A.
187  firstHalves = [x//2 for x in A.shape]
188  secondHalves = [x-y for x, y in zip(A.shape, firstHalves)]
189  for d1, d2 in zip(newShape, A.shape):
190  if d1 < d2:
191  raise ValueError("Newshape dimensions must be greater or equal")
192  else:
193  # Inverse operation: Opposite rounding
194  secondHalves = [x//2 for x in newShape]
195  firstHalves = [x-y for x, y in zip(newShape, secondHalves)]
196  for d1, d2 in zip(newShape, A.shape):
197  if d1 > d2:
198  raise ValueError("Newshape dimensions must be smaller or equal")
199 
200  if dtype is None:
201  dtype = A.dtype
202 
203  R = np.zeros(newShape, dtype=dtype)
204  R[-firstHalves[0]:, -firstHalves[1]:] = A[:firstHalves[0], :firstHalves[1]]
205  R[:secondHalves[0], -firstHalves[1]:] = A[-secondHalves[0]:, :firstHalves[1]]
206  R[:secondHalves[0], :secondHalves[1]] = A[-secondHalves[0]:, -secondHalves[1]:]
207  R[-firstHalves[0]:, :secondHalves[1]] = A[:firstHalves[0], -secondHalves[1]:]
208  return R
209 
210  def initializeSubImage(self, fullExp, innerBox, outerBox, noiseMeanVar, useNoise=True):
211  """Initializes a sub image.
212 
213  Parameters
214  ----------
215  fullExp : `lsst.afw.image.Exposure`
216  The full exposure to cut sub image from.
217  innerBox : `lsst.geom.Box2I`
218  The useful area of the calculation up to the whole bounding box of
219  ``fullExp``. ``fullExp`` must contain this box.
220  outerBox : `lsst.geom.Box2I`
221  The overall cutting area. ``outerBox`` must be at least 1 pixel larger
222  than ``inneBox`` in all directions and may not be fully contained by
223  ``fullExp``.
224  noiseMeanVar : `float` > 0.
225  The noise variance level to initialize variance plane and to generate
226  white noise for the non-overlapping region.
227  useNoise : `bool`, optional
228  If True, generate white noise for non-overlapping region. Otherwise,
229  zero padding will be used in the non-overlapping region.
230 
231  Notes
232  -----
233  ``innerBox``, ``outerBox`` must be in the PARENT system of ``fullExp``.
234 
235  Returns
236  -------
237  result : `lsst.pipe.base.Struct`
238  - ``subImg``, ``subVarImg`` : `lsst.afw.image.ImageD`
239  The new sub image and its sub variance plane.
240 
241  Notes
242  -----
243  Supports the non-grid option when ``innerBox`` equals to the
244  bounding box of ``fullExp``.
245  """
246  fullBox = fullExp.getBBox()
247  subImg = afwImage.ImageD(outerBox, 0)
248  subVarImg = afwImage.ImageD(outerBox, noiseMeanVar)
249  borderBoxes = self.splitBordersplitBorder(innerBox, outerBox)
250  # Initialize the border region that are not fully within the exposure
251  if useNoise:
252  noiseSig = np.sqrt(noiseMeanVar)
253  for box in borderBoxes:
254  if not fullBox.contains(box):
255  R = subImg[box].array
256  R[...] = self.rngrng.normal(scale=noiseSig, size=R.shape)
257  # Copy data to the fully contained inner region, allowing type conversion
258  subImg[innerBox].array[...] = fullExp.image[innerBox].array
259  subVarImg[innerBox].array[...] = fullExp.variance[innerBox].array
260  # Copy data to border regions that have at least a partial overlap
261  for box in borderBoxes:
262  overlapBox = box.clippedTo(fullBox)
263  if not overlapBox.isEmpty():
264  subImg[overlapBox].array[...] = fullExp.image[overlapBox].array
265  subVarImg[overlapBox].array[...] = fullExp.variance[overlapBox].array
266  return pipeBase.Struct(image=subImg, variance=subVarImg)
267 
268  @staticmethod
269  def estimateMatchingKernelSize(psf1, psf2):
270  """Estimate the image space size of the matching kernels.
271 
272  Return ten times the larger Gaussian sigma estimate but at least
273  the largest of the original psf dimensions.
274 
275  Parameters
276  ----------
277  psf1, psf2 : `lsst.afw.detection.Psf`
278  The PSFs of the two input exposures.
279 
280  Returns
281  -------
282  size : `int`
283  Conservative estimate for matching kernel size in pixels.
284  This is the minimum padding around the inner region at each side.
285 
286  Notes
287  -----
288  """
289  sig1 = psf1.computeShape().getDeterminantRadius()
290  sig2 = psf2.computeShape().getDeterminantRadius()
291  sig = max(sig1, sig2)
292  psfBBox1 = psf1.computeBBox()
293  psfBBox2 = psf2.computeBBox()
294  return max(10 * sig, psfBBox1.getWidth(), psfBBox1.getHeight(),
295  psfBBox2.getWidth(), psfBBox2.getHeight())
296 
297  @staticmethod
298  def splitBorder(innerBox, outerBox):
299  """Split the border area around the inner box into 8 disjunct boxes.
300 
301  Parameters
302  ----------
303  innerBox : `lsst.geom.Box2I`
304  The inner box.
305  outerBox : `lsst.geom.Box2I`
306  The outer box. It must be at least 1 pixel larger in each direction than the inner box.
307 
308  Returns
309  -------
310  resultBoxes : `list` of 8 boxes covering the edge around innerBox
311 
312  Notes
313  -----
314  The border boxes do not overlap. The border is covered counter clockwise
315  starting from lower left corner.
316 
317  Raises
318  ------
319  ValueError : If ``outerBox`` is not larger than ``innerBox``.
320  """
321  innerBox = innerBox.dilatedBy(1)
322  if not outerBox.contains(innerBox):
323  raise ValueError("OuterBox must be larger by at least 1 pixel in all directions")
324 
325  # ccw sequence of corners
326  o1, o2, o3, o4 = outerBox.getCorners()
327  i1, i2, i3, i4 = innerBox.getCorners()
328  p1 = Point2I(outerBox.minX, innerBox.minY)
329  p2 = Point2I(innerBox.maxX, outerBox.minY)
330  p3 = Point2I(outerBox.maxX, innerBox.maxY)
331  p4 = Point2I(innerBox.minX, outerBox.maxY)
332 
333  # The 8 border boxes ccw starting from lower left
334  pointPairs = ((o1, i1), (i1 + Extent2I(1, 0), p2 + Extent2I(-1, 0)), (o2, i2),
335  (i2 + Extent2I(0, 1), p3 + Extent2I(0, -1)), (o3, i3),
336  (i3 + Extent2I(-1, 0), p4 + Extent2I(1, 0)), (o4, i4),
337  (i4 + Extent2I(0, -1), p1 + Extent2I(0, 1)))
338  return [Box2I(x, y, invert=True) for (x, y) in pointPairs]
339 
340  @staticmethod
341  def generateGrid(imageBox, minEdgeDims, innerBoxDims, minTotalDims=None, powerOfTwo=False):
342  """Generate a splitting grid for an image.
343 
344  The inner boxes cover the input image without overlap, the edges around the inner boxes do overlap
345  and go beyond the image at the image edges.
346 
347  Parameters
348  ----------
349  imageBox : `lsst.geom.Box2I`
350  Bounding box of the exposure to split.
351  minEdgeDims : `lsst.geom.Extent2I`
352  Minimum edge width in (x,y) directions each side.
353  innerBoxDims : `lsst.geom.Extent2I`
354  Minimum requested inner box dimensions (x,y).
355  The actual dimensions can be larger due to rounding.
356  minTotalDims: `lsst.geom.Extent2I`, optional
357  If provided, minimum total outer dimensions (x,y). The edge will be increased until satisfied.
358  powerOfTwo : `bool`, optional
359  If True, the outer box dimensions should be rounded up to a power of 2
360  by increasing the border size. This is up to 8192, above this size,
361  rounding up is disabled.
362 
363  Notes
364  -----
365  Inner box dimensions are chosen to be as uniform as they can, remainder pixels at the edge of the
366  input will be appended to the last column/row boxes.
367 
368  See diffimTests/tickets/DM-28928_spatial_grid notebooks for demonstration of this code.
369 
370  This method can be used for both PARENT and LOCAL bounding boxes.
371 
372  The outerBox dimensions are always even.
373 
374  Returns
375  -------
376  boxList : `list` of `lsst.pipe.base.Struct`
377  ``innerBox``, ``outerBox`` : `lsst.geom.Box2I`, inner boxes and overlapping border around them.
378 
379  """
380  powersOf2 = np.array([16, 32, 64, 128, 256, 512, 1024, 2048, 4096, 8192])
381  doubleEdgeDims = minEdgeDims * 2
382  width, height = imageBox.getDimensions()
383  nX = width // innerBoxDims.x # Round down
384  if nX > 0:
385  innerWidth = width // nX # Round down
386  else:
387  innerWidth = width
388  nX = 1
389  xCorners = np.zeros(nX + 1)
390  xCorners[:-1] = np.arange(nX)*innerWidth + imageBox.minX
391  xCorners[-1] = imageBox.endX
392 
393  nY = height // innerBoxDims.y # Round down
394  if nY > 0:
395  innerHeight = height // nY # Round down
396  else:
397  innerHeight = height
398  nY = 1
399  yCorners = np.zeros(nY + 1)
400  yCorners[:-1] = np.arange(nY)*innerHeight + imageBox.minY
401  yCorners[-1] = imageBox.endY
402 
403  boxes = []
404 
405  for i_y in range(nY):
406  for i_x in range(nX):
407  innerBox = Box2I(Point2I(xCorners[i_x], yCorners[i_y]),
408  Point2I(xCorners[i_x + 1] - 1, yCorners[i_y + 1] - 1))
409 
410  paddedWidth = innerBox.width + doubleEdgeDims.x
411  if minTotalDims is not None and paddedWidth < minTotalDims.width:
412  paddedWidth = minTotalDims.width
413  if powerOfTwo:
414  i2x = np.searchsorted(powersOf2, paddedWidth, side='left')
415  if i2x < len(powersOf2):
416  paddedWidth = powersOf2[i2x]
417  if paddedWidth % 2 == 1:
418  paddedWidth += 1 # Ensure total width is even
419 
420  totalXedge = paddedWidth - innerBox.width
421 
422  paddedHeight = innerBox.height + doubleEdgeDims.y
423  if minTotalDims is not None and paddedHeight < minTotalDims.height:
424  paddedHeight = minTotalDims.height
425  if powerOfTwo:
426  i2y = np.searchsorted(powersOf2, paddedHeight, side='left')
427  if i2y < len(powersOf2):
428  paddedHeight = powersOf2[i2y]
429  if paddedHeight % 2 == 1:
430  paddedHeight += 1 # Ensure total height is even
431  totalYedge = paddedHeight - innerBox.height
432  outerBox = Box2I(Point2I(innerBox.minX - totalXedge//2, innerBox.minY - totalYedge//2),
433  Extent2I(paddedWidth, paddedHeight))
434  boxes.append(pipeBase.Struct(innerBox=innerBox, outerBox=outerBox))
435  return boxes
436 
437  def makeSpatialPsf(self, gridPsfs):
438  """Construct a CoaddPsf based on PSFs from individual sub image solutions.
439 
440  Parameters
441  ----------
442  gridPsfs : iterable of `lsst.pipe.base.Struct`
443  Iterable of bounding boxes (``bbox``) and Psf solutions (``psf``).
444 
445  Returns
446  -------
447  psf : `lsst.meas.algorithms.CoaddPsf`
448  A psf constructed from the PSFs of the individual subExposures.
449  """
450  schema = afwTable.ExposureTable.makeMinimalSchema()
451  schema.addField("weight", type="D", doc="Coadd weight")
452  mycatalog = afwTable.ExposureCatalog(schema)
453 
454  # We're just using the exposure's WCS (assuming that the subExposures'
455  # WCSs are the same, which they better be!).
456  wcsref = self.fullExp1fullExp1.getWcs()
457  for i, res in enumerate(gridPsfs):
458  record = mycatalog.getTable().makeRecord()
459  record.setPsf(res.psf)
460  record.setWcs(wcsref)
461  record.setBBox(res.bbox)
462  record['weight'] = 1.0
463  record['id'] = i
464  mycatalog.append(record)
465 
466  # create the CoaddPsf
467  psf = measAlg.CoaddPsf(mycatalog, wcsref, 'weight')
468  return psf
469 
470  def padAndFftImage(self, imgArr):
471  """Prepare and forward FFT an image array.
472 
473  Parameters
474  ----------
475  imgArr : `numpy.ndarray` of `float`
476  Original array. In-place modified as `numpy.nan` and `numpy.inf` are replaced by
477  array mean.
478 
479  Returns
480  -------
481  result : `lsst.pipe.base.Struct`
482  - ``imFft`` : `numpy.ndarray` of `numpy.complex`.
483  FFT of image.
484  - ``filtInf``, ``filtNaN`` : `numpy.ndarray` of `bool`
485 
486  Notes
487  -----
488  Save location of non-finite values for restoration, and replace them
489  with image mean values. Re-center and zero pad array by `padCenterOriginArray`.
490  """
491  filtInf = np.isinf(imgArr)
492  filtNaN = np.isnan(imgArr)
493  imgArr[filtInf] = np.nan
494  imgArr[filtInf | filtNaN] = np.nanmean(imgArr)
495  self.log.debug("Replacing {} Inf and {} NaN values.".format(
496  np.sum(filtInf), np.sum(filtNaN)))
497  imgArr = self.padCenterOriginArraypadCenterOriginArray(imgArr, self.freqSpaceShapefreqSpaceShape)
498  imgArr = np.fft.fft2(imgArr)
499  return pipeBase.Struct(imFft=imgArr, filtInf=filtInf, filtNaN=filtNaN)
500 
501  def removeNonFinitePixels(self, imgArr):
502  """Replace non-finite pixel values in-place.
503 
504  Save the locations of non-finite values for restoration, and replace them
505  with image mean values.
506 
507  Parameters
508  ----------
509  imgArr : `numpy.ndarray` of `float`
510  The image array. Non-finite values are replaced in-place in this array.
511 
512  Returns
513  -------
514  result : `lsst.pipe.base.Struct`
515  - ``filtInf``, ``filtNaN`` : `numpy.ndarray` of `bool`
516  The filter of the pixel values that were inf or nan.
517  """
518  filtInf = np.isinf(imgArr)
519  filtNaN = np.isnan(imgArr)
520  # Masked edge and bad pixels could also be removed here in the same way
521  # in the future
522  imgArr[filtInf] = np.nan
523  imgArr[filtInf | filtNaN] = np.nanmean(imgArr)
524  self.log.debugf("Replacing {} Inf and {} NaN values.",
525  np.sum(filtInf), np.sum(filtNaN))
526  return pipeBase.Struct(filtInf=filtInf, filtNaN=filtNaN)
527 
528  def inverseFftAndCropImage(self, imgArr, origSize, filtInf=None, filtNaN=None, dtype=None):
529  """Inverse FFT and crop padding from image array.
530 
531  Parameters
532  ----------
533  imgArr : `numpy.ndarray` of `numpy.complex`
534  Fourier space array representing a real image.
535 
536  origSize : `tuple` of `int`
537  Original unpadded shape tuple of the image to be cropped to.
538 
539  filtInf, filtNan : `numpy.ndarray` of bool or int, optional
540  If specified, they are used as index arrays for ``result`` to set values to
541  `numpy.inf` and `numpy.nan` respectively at these positions.
542 
543  dtype : `numpy.dtype`, optional
544  Dtype of result array to cast return values to implicitly. This is to
545  spare one array copy operation at reducing double precision to single.
546  If `None` result inherits dtype of `imgArr`.
547 
548  Returns
549  -------
550  result : `numpy.ndarray` of `dtype`
551  """
552  imgNew = np.fft.ifft2(imgArr)
553  imgNew = imgNew.real
554  imgNew = self.padCenterOriginArraypadCenterOriginArray(imgNew, origSize, useInverse=True, dtype=dtype)
555  if filtInf is not None:
556  imgNew[filtInf] = np.inf
557  if filtNaN is not None:
558  imgNew[filtNaN] = np.nan
559  return imgNew
560 
561  @staticmethod
562  def computePsfAtCenter(exposure):
563  """Computes the PSF image at the bbox center point.
564 
565  This may be at a fractional pixel position.
566 
567  Parameters
568  ----------
569  exposure : `lsst.afw.image.Exposure`
570  Exposure with psf.
571 
572  Returns
573  -------
574  psfImg : `lsst.afw.image.Image`
575  Calculated psf image.
576  """
577  pBox = exposure.getBBox()
578  cen = pBox.getCenter()
579  psf = exposure.getPsf()
580  psfImg = psf.computeKernelImage(cen) # Centered and normed
581  return psfImg
582 
583  @staticmethod
584  def subtractImageMean(image, mask, statsControl):
585  """In-place subtraction of sigma-clipped mean of the image.
586 
587  Parameters
588  ----------
589  image : `lsst.afw.image.Image`
590  Image to manipulate. Its sigma clipped mean is in-place subtracted.
591 
592  mask : `lsst.afw.image.Mask`
593  Mask to use for ignoring pixels.
594 
595  statsControl : `lsst.afw.math.StatisticsControl`
596  Config of sigma clipped mean statistics calculation.
597 
598  Returns
599  -------
600  None
601 
602  Raises
603  ------
604  ValueError : If image mean is nan.
605  """
606  statObj = afwMath.makeStatistics(image, mask,
607  afwMath.MEANCLIP, statsControl)
608  mean = statObj.getValue(afwMath.MEANCLIP)
609  if not np.isnan(mean):
610  image -= mean
611  else:
612  raise ValueError("Image mean is NaN.")
613 
614  def prepareFullExposure(self, exposure1, exposure2, correctBackground=False):
615  """Performs calculations that apply to the full exposures once only.
616 
617  Parameters
618  ----------
619 
620  exposure1, exposure2 : `lsst.afw.image.Exposure`
621  The input exposures. Copies are made for internal calculations.
622 
623  correctBackground : `bool`, optional
624  If True, subtracts sigma-clipped mean of exposures. The algorithm
625  assumes zero expectation value at background pixels.
626 
627  Returns
628  -------
629  None
630 
631  Notes
632  -----
633  Set a number of instance fields with pre-calculated values.
634 
635  Raises
636  ------
637  ValueError : If photometric calibrations are not available while
638  ``config.scaleByCalibration`` equals True.
639  """
641  self.statsControlstatsControl.setNumSigmaClip(3.)
642  self.statsControlstatsControl.setNumIter(3)
643  self.statsControlstatsControl.setAndMask(afwImage.Mask.getPlaneBitMask(
644  self.config.ignoreMaskPlanes))
645 
646  exposure1 = exposure1.clone()
647  exposure2 = exposure2.clone()
648  # Fallback values if sub exposure variance calculation is problematic
649  sig1 = np.sqrt(self._computeVarianceMean_computeVarianceMean(exposure1))
650  self.fullExpVar1fullExpVar1 = sig1*sig1
651  sig2 = np.sqrt(self._computeVarianceMean_computeVarianceMean(exposure2))
652  self.fullExpVar2fullExpVar2 = sig2*sig2
653 
654  # If 'scaleByCalibration' is True then these norms are overwritten
655  if self.config.scaleByCalibration:
656  calibObj1 = exposure1.getPhotoCalib()
657  calibObj2 = exposure2.getPhotoCalib()
658  if calibObj1 is None or calibObj2 is None:
659  raise ValueError("Photometric calibrations are not available for both exposures.")
660  mImg1 = calibObj1.calibrateImage(exposure1.maskedImage)
661  mImg2 = calibObj2.calibrateImage(exposure2.maskedImage)
662  self.F1F1 = 1.
663  self.F2F2 = 1.
664  else:
665  self.F1F1 = self.config.templateFluxScaling # default is 1
666  self.F2F2 = self.config.scienceFluxScaling # default is 1
667  mImg1 = exposure1.maskedImage
668  mImg2 = exposure2.maskedImage
669 
670  # mImgs can be in-place modified
671  if correctBackground:
672  self.subtractImageMeansubtractImageMean(mImg1.image, mImg1.mask, self.statsControlstatsControl)
673  self.subtractImageMeansubtractImageMean(mImg2.image, mImg2.mask, self.statsControlstatsControl)
674 
675  # Determine border size
676  self.borderSizeborderSize = self.estimateMatchingKernelSizeestimateMatchingKernelSize(exposure1.getPsf(), exposure2.getPsf())
677  self.log.debugf("Minimum padding border size: {} pixels", self.borderSizeborderSize)
678  # Remove non-finite values from the images in-place
679  self.filtsImg1filtsImg1 = self.removeNonFinitePixelsremoveNonFinitePixels(mImg1.image.array)
680  self.filtsImg2filtsImg2 = self.removeNonFinitePixelsremoveNonFinitePixels(mImg2.image.array)
681  self.filtsVar1filtsVar1 = self.removeNonFinitePixelsremoveNonFinitePixels(mImg1.variance.array)
682  self.filtsVar2filtsVar2 = self.removeNonFinitePixelsremoveNonFinitePixels(mImg2.variance.array)
683 
684  exposure1.maskedImage = mImg1
685  exposure2.maskedImage = mImg2
686 
687  self.fullExp1fullExp1 = exposure1
688  self.fullExp2fullExp2 = exposure2
689 
690  def prepareSubExposure(self, localCutout, psf1=None, psf2=None, sig1=None, sig2=None):
691  """Perform per-sub exposure preparations.
692 
693  Parameters
694  ----------
695  sig1, sig2 : `float`, optional
696  For debug purposes only, copnsider that the image
697  may already be rescaled by the photometric calibration.
698  localCutout : `lsst.pipe.base.Struct`
699  - innerBox, outerBox: `lsst.geom.Box2I` LOCAL inner and outer boxes
700  psf1, psf2 : `lsst.afw.detection.Psf`, optional
701  If specified, use given psf as the sub exposure psf. For debug purposes.
702  sig1, sig2 : `float`, optional
703  If specified, use value as the sub-exposures' background noise sigma value.
704 
705  Returns
706  -------
707  None
708 
709  """
710  self.log.debugf("Processing LOCAL cell w/ inner box:{}, outer box:{}",
711  localCutout.innerBox, localCutout.outerBox)
712  # The PARENT origin cutout boxes for the two exposures
713  self.cutBoxes1cutBoxes1 = pipeBase.Struct(
714  innerBox=localCutout.innerBox.shiftedBy(Extent2I(self.fullExp1fullExp1.getXY0())),
715  outerBox=localCutout.outerBox.shiftedBy(Extent2I(self.fullExp1fullExp1.getXY0())))
716  self.cutBoxes2cutBoxes2 = pipeBase.Struct(
717  innerBox=localCutout.innerBox.shiftedBy(Extent2I(self.fullExp2fullExp2.getXY0())),
718  outerBox=localCutout.outerBox.shiftedBy(Extent2I(self.fullExp2fullExp2.getXY0())))
719  # The sub-exposure views of the useful inner area of this grid cell
720  innerSubExp1 = self.fullExp1fullExp1[self.cutBoxes1cutBoxes1.innerBox]
721  innerSubExp2 = self.fullExp2fullExp2[self.cutBoxes2cutBoxes2.innerBox]
722  if psf1 is None:
723  self.subExpPsf1subExpPsf1 = self.computePsfAtCentercomputePsfAtCenter(innerSubExp1)
724  else:
725  self.subExpPsf1subExpPsf1 = psf1
726  if psf2 is None:
727  self.subExpPsf2subExpPsf2 = self.computePsfAtCentercomputePsfAtCenter(innerSubExp2)
728  else:
729  self.subExpPsf2subExpPsf2 = psf2
730  self.checkCentroidscheckCentroids(self.subExpPsf1subExpPsf1.array, self.subExpPsf2subExpPsf2.array)
731  psfBBox1 = self.subExpPsf1subExpPsf1.getBBox()
732  psfBBox2 = self.subExpPsf2subExpPsf2.getBBox()
733  self.psfShape1psfShape1 = (psfBBox1.getHeight(), psfBBox1.getWidth())
734  self.psfShape2psfShape2 = (psfBBox2.getHeight(), psfBBox2.getWidth())
735  # sig1 and sig2 should not be set externally, just for debug purpose
736  if sig1 is None:
737  sig1 = np.sqrt(self._computeVarianceMean_computeVarianceMean(innerSubExp1))
738  if sig1 > 0.: # Not zero and not nan
739  self.subExpVar1subExpVar1 = sig1*sig1
740  else:
741  self.subExpVar1subExpVar1 = self.fullExpVar1fullExpVar1
742  if sig2 is None:
743  sig2 = np.sqrt(self._computeVarianceMean_computeVarianceMean(innerSubExp2))
744  if sig2 > 0.: # Not zero and not nan
745  self.subExpVar2subExpVar2 = sig2*sig2
746  else:
747  self.subExpVar2subExpVar2 = self.fullExpVar2fullExpVar2
748  # Initialize random number generator to a deterministic state
749  self.rngrng = np.random.default_rng(seed=np.array([self.subExpVar1subExpVar1]).view(int))
750  self.freqSpaceShapefreqSpaceShape = (localCutout.outerBox.getHeight(), localCutout.outerBox.getWidth())
751 
752  self.subImg1subImg1 = self.initializeSubImageinitializeSubImage(
753  self.fullExp1fullExp1, self.cutBoxes1cutBoxes1.innerBox, self.cutBoxes1cutBoxes1.outerBox,
754  self.subExpVar1subExpVar1, useNoise=True)
755  self.subImg2subImg2 = self.initializeSubImageinitializeSubImage(
756  self.fullExp2fullExp2, self.cutBoxes2cutBoxes2.innerBox, self.cutBoxes2cutBoxes2.outerBox,
757  self.subExpVar2subExpVar2, useNoise=True)
758 
759  D = self.padCenterOriginArraypadCenterOriginArray(self.subImg1subImg1.image.array, self.freqSpaceShapefreqSpaceShape)
760  self.subImgFft1subImgFft1 = np.fft.fft2(D)
761  D = self.padCenterOriginArraypadCenterOriginArray(self.subImg1subImg1.variance.array, self.freqSpaceShapefreqSpaceShape)
762  self.subVarImgFft1subVarImgFft1 = np.fft.fft2(D)
763 
764  D = self.padCenterOriginArraypadCenterOriginArray(self.subImg2subImg2.image.array, self.freqSpaceShapefreqSpaceShape)
765  self.subImgFft2subImgFft2 = np.fft.fft2(D)
766  D = self.padCenterOriginArraypadCenterOriginArray(self.subImg2subImg2.variance.array, self.freqSpaceShapefreqSpaceShape)
767  self.subVarImgFft2subVarImgFft2 = np.fft.fft2(D)
768 
769  D = self.padCenterOriginArraypadCenterOriginArray(self.subExpPsf1subExpPsf1.array, self.freqSpaceShapefreqSpaceShape)
770  self.psfFft1psfFft1 = np.fft.fft2(D)
771  D = self.padCenterOriginArraypadCenterOriginArray(self.subExpPsf2subExpPsf2.array, self.freqSpaceShapefreqSpaceShape)
772  self.psfFft2psfFft2 = np.fft.fft2(D)
773 
774  @staticmethod
776  """Square the argument in pixel space.
777 
778  Parameters
779  ----------
780  D : 2D `numpy.ndarray` of `numpy.complex`
781  Fourier transform of a real valued array.
782 
783  Returns
784  -------
785  R : `numpy.ndarray` of `numpy.complex`
786 
787  Notes
788  -----
789  ``D`` is to be inverse Fourier transformed, squared and then
790  forward Fourier transformed again, i.e. an autoconvolution in Fourier space.
791  This operation is not distributive over multiplication.
792  ``pixelSpaceSquare(A*B) != pixelSpaceSquare(A)*pixelSpaceSquare(B)``
793  """
794  R = np.real(np.fft.ifft2(D))
795  R *= R
796  R = np.fft.fft2(R)
797  return R
798 
799  @staticmethod
800  def getCentroid(A):
801  """Calculate the centroid coordinates of a 2D array.
802 
803  Parameters
804  ----------
805  A : 2D `numpy.ndarray` of `float`
806  The input array. Must not be all exact zero.
807 
808  Notes
809  -----
810  Calculates the centroid as if the array represented a 2D geometrical shape with
811  weights per cell, allowing for "negative" weights. If sum equals to exact (float) zero,
812  calculates centroid of absolute value array.
813 
814  The geometrical center is defined as (0,0), independently of the array shape.
815  For an odd dimension, this is the center of the center pixel,
816  for an even dimension, this is between the two center pixels.
817 
818  Returns
819  -------
820  ycen, xcen : `tuple` of `float`
821 
822  """
823  s = np.sum(A)
824  if s == 0.:
825  A = np.fabs(A)
826  s = np.sum(A)
827  w = np.arange(A.shape[0], dtype=float) - (A.shape[0] - 1.)/2
828  ycen = np.sum(w[:, np.newaxis]*A)/s
829  w = np.arange(A.shape[1], dtype=float) - (A.shape[1] - 1.)/2
830  xcen = np.sum(w[np.newaxis, :]*A)/s
831 
832  return ycen, xcen
833 
834  def checkCentroids(self, psfArr1, psfArr2):
835  """Check whether two PSF array centroids' distance is within tolerance.
836 
837  Parameters
838  ----------
839  psfArr1, psfArr2 : `numpy.ndarray` of `float`
840  Input PSF arrays to check.
841 
842  Returns
843  -------
844  None
845 
846  Raises
847  ------
848  ValueError:
849  Centroid distance exceeds `config.maxPsfCentroidDist` pixels.
850  """
851  yc1, xc1 = self.getCentroidgetCentroid(psfArr1)
852  yc2, xc2 = self.getCentroidgetCentroid(psfArr2)
853  dy = yc2 - yc1
854  dx = xc2 - xc1
855  if dy*dy + dx*dx > self.config.maxPsfCentroidDist*self.config.maxPsfCentroidDist:
856  raise ValueError(
857  f"PSF centroids are offset by more than {self.config.maxPsfCentroidDist:.2f} pixels.")
858 
859  def calculateFourierDiffim(self, psf1, im1, varPlane1, F1, varMean1,
860  psf2, im2, varPlane2, F2, varMean2, calculateScore=True):
861  """Convolve and subtract two images in Fourier space.
862 
863  Calculate the ZOGY proper difference image, score image and their PSFs.
864  All input and output arrays are in Fourier space.
865 
866  Parameters
867  ----------
868  psf1, psf2, im1, im2, varPlane1, varPlane2 : `numpy.ndarray` of `numpy.complex`,
869  shape ``self.freqSpaceShape``
870  Psf, image and variance plane arrays respectively.
871  All arrays must be already in Fourier space.
872 
873  varMean1, varMean2: `numpy.float` > 0.
874  Average per-pixel noise variance in im1, im2 respectively. Used as weighing
875  of input images. Must be greater than zero.
876 
877  F1, F2 : `numpy.float` > 0.
878  Photometric scaling of the images. See eqs. (5)--(9)
879 
880  calculateScore : `bool`, optional
881  If True (default), calculate and return the detection significance (score) image.
882  Otherwise, these return fields are `None`.
883 
884  Returns
885  -------
886  result : `pipe.base.Struct`
887  All arrays are in Fourier space and have shape ``self.freqSpaceShape``.
888  - ``Fd`` : `float`
889  Photometric level of ``D``.
890  - ``D`` : `numpy.ndarray` of `numpy.complex`
891  The difference image.
892  - ``varplaneD`` : `numpy.ndarray` of `numpy.complex`
893  Variance plane of ``D``.
894  - ``Pd`` : `numpy.ndarray` of `numpy.complex`
895  PSF of ``D``.
896  - ``S`` : `numpy.ndarray` of `numpy.complex` or `None`
897  Significance (score) image.
898  - ``varplaneS`` : `numpy.ndarray` of `numpy.complex` or `None`
899  Variance plane of ``S``.
900  - ``Ps`` : `numpy.ndarray` of `numpy.complex`
901  PSF of ``S``.
902 
903  Notes
904  -----
905  All array inputs and outputs are Fourier-space images with shape of
906  `self.freqSpaceShape` in this method.
907 
908  ``varMean1``, ``varMean2`` quantities are part of the noise model and not to be confused
909  with the variance of image frequency components or with ``varPlane1``, ``varPlane2`` that
910  are the Fourier transform of the variance planes.
911  """
912  var1F2Sq = varMean1*F2*F2
913  var2F1Sq = varMean2*F1*F1
914  # We need reals for comparison, also real operations are usually faster
915  psfAbsSq1 = np.real(np.conj(psf1)*psf1)
916  psfAbsSq2 = np.real(np.conj(psf2)*psf2)
917  FdDenom = np.sqrt(var1F2Sq + var2F1Sq) # one number
918 
919  # Secure positive limit to avoid floating point operations resulting in exact zero
920  tiny = np.finfo(psf1.dtype).tiny * 100
921  sDenom = var1F2Sq*psfAbsSq2 + var2F1Sq*psfAbsSq1 # array, eq. (12)
922  # Frequencies where both psfs are too close to zero.
923  # We expect this only in cases when psf1, psf2 are identical,
924  # and either having very well sampled Gaussian tails
925  # or having "edges" such that some sinc-like zero crossings are found at symmetry points
926  #
927  # if sDenom < tiny then it can be == 0. -> `denom` = 0. and 0/0 occur at `c1` , `c2`
928  # if we keep SDenom = tiny, denom ~ O(sqrt(tiny)), Pd ~ O(sqrt(tiny)), S ~ O(sqrt(tiny)*tiny) == 0
929  # Where S = 0 then Pd = 0 and D should still yield the same variance ~ O(1)
930  # For safety, we set S = 0 explicitly, too, though it should be unnecessary.
931  fltZero = sDenom < tiny
932  nZero = np.sum(fltZero)
933  self.log.debug(f"There are {nZero} frequencies where both FFTd PSFs are close to zero.")
934  if nZero > 0:
935  # We expect only a small fraction of such frequencies
936  fltZero = np.nonzero(fltZero) # Tuple of index arrays
937  sDenom[fltZero] = tiny # Avoid division problem but overwrite result anyway
938  denom = np.sqrt(sDenom) # array, eq. (13)
939 
940  c1 = F2*psf2/denom
941  c2 = F1*psf1/denom
942  if nZero > 0:
943  c1[fltZero] = F2/FdDenom
944  c2[fltZero] = F1/FdDenom
945  D = c1*im1 - c2*im2 # Difference image eq. (13)
946  varPlaneD = self.pixelSpaceSquarepixelSpaceSquare(c1)*varPlane1 + self.pixelSpaceSquarepixelSpaceSquare(c2)*varPlane2 # eq. (26)
947 
948  Pd = FdDenom*psf1*psf2/denom # Psf of D eq. (14)
949  if nZero > 0:
950  Pd[fltZero] = 0
951 
952  Fd = F1*F2/FdDenom # Flux scaling of D eq. (15)
953  if calculateScore:
954  c1 = F1*F2*F2*np.conj(psf1)*psfAbsSq2/sDenom
955  c2 = F2*F1*F1*np.conj(psf2)*psfAbsSq1/sDenom
956  if nZero > 0:
957  c1[fltZero] = 0
958  c2[fltZero] = 0
959  S = c1*im1 - c2*im2 # eq. (12)
960  varPlaneS = self.pixelSpaceSquarepixelSpaceSquare(c1)*varPlane1 + self.pixelSpaceSquarepixelSpaceSquare(c2)*varPlane2
961  Ps = np.conj(Pd)*Pd # eq. (17) Source detection expects a PSF
962  else:
963  S = None
964  Ps = None
965  varPlaneS = None
966  return pipeBase.Struct(D=D, Pd=Pd, varPlaneD=varPlaneD, Fd=Fd,
967  S=S, Ps=Ps, varPlaneS=varPlaneS)
968 
969  @staticmethod
970  def calculateMaskPlane(mask1, mask2, effPsf1=None, effPsf2=None):
971  """Calculate the mask plane of the difference image.
972 
973  Parameters
974  ----------
975  mask1, maks2 : `lsst.afw.image.Mask`
976  Mask planes of the two exposures.
977 
978 
979  Returns
980  -------
981  diffmask : `lsst.afw.image.Mask`
982  Mask plane for the subtraction result.
983 
984  Notes
985  -----
986  TODO DM-25174 : Specification of effPsf1, effPsf2 are not yet supported.
987  """
988 
989  # mask1 x effPsf2 | mask2 x effPsf1
990  if effPsf1 is not None or effPsf2 is not None:
991  # TODO: DM-25174 effPsf1, effPsf2: the effective psf for cross-blurring.
992  # We need a "size" approximation of the c1 and c2 coefficients to make effPsfs
993  # Also convolution not yet supports mask-only operation
994  raise NotImplementedError("Mask plane only 'convolution' operation is not yet supported")
995  R = mask1.clone()
996  R |= mask2
997  return R
998 
999  @staticmethod
1001  """Create a non spatially varying PSF from a `numpy.ndarray`.
1002 
1003  Parameters
1004  ----------
1005  A : `numpy.ndarray`
1006  2D array to use as the new psf image. The pixels are copied.
1007 
1008  Returns
1009  -------
1010  psfNew : `lsst.meas.algorithms.KernelPsf`
1011  The constructed PSF.
1012  """
1013  psfImg = afwImage.ImageD(A.astype(np.float64, copy=True), deep=False)
1014  psfNew = measAlg.KernelPsf(afwMath.FixedKernel(psfImg))
1015  return psfNew
1016 
1017  def pasteSubDiffImg(self, ftDiff, diffExp, scoreExp=None):
1018  """Paste sub image results back into result Exposure objects.
1019 
1020  Parameters
1021  ----------
1022  ftDiff : `lsst.pipe.base.Struct`
1023  Result struct by `calculateFourierDiffim`.
1024  diffExp : `lsst.afw.image.Exposure`
1025  The result exposure to paste into the sub image result.
1026  Must be dimensions and dtype of ``self.fullExp1``.
1027  scoreExp : `lsst.afw.image.Exposure` or `None`
1028  The result score exposure to paste into the sub image result.
1029  Must be dimensions and dtype of ``self.fullExp1``.
1030  If `None`, the score image results are disregarded.
1031 
1032  Returns
1033  -------
1034  None
1035 
1036  Notes
1037  -----
1038  The PSF of the score image is just to make the score image resemble a
1039  regular exposure and to study the algorithm performance.
1040 
1041  Add an entry to the ``self.gridPsfs`` list.
1042 
1043  gridPsfs : `list` of `lsst.pipe.base.Struct`
1044  - ``bbox`` : `lsst.geom.Box2I`
1045  The inner region of the grid cell.
1046  - ``Pd`` : `lsst.meas.algorithms.KernelPsf`
1047  The diffim PSF in this cell.
1048  - ``Ps`` : `lsst.meas.algorithms.KernelPsf` or `None`
1049  The score image PSF in this cell or `None` if the score
1050  image was not calculated.
1051  """
1052  D = self.inverseFftAndCropImageinverseFftAndCropImage(
1053  ftDiff.D, self.freqSpaceShapefreqSpaceShape, dtype=self.fullExp1fullExp1.image.dtype)
1054  varPlaneD = self.inverseFftAndCropImageinverseFftAndCropImage(
1055  ftDiff.varPlaneD, self.freqSpaceShapefreqSpaceShape, dtype=self.fullExp1fullExp1.variance.dtype)
1056  Pd = self.inverseFftAndCropImageinverseFftAndCropImage(
1057  ftDiff.Pd, self.psfShape1psfShape1, dtype=self.subExpPsf1subExpPsf1.dtype)
1058  sumPd = np.sum(Pd)
1059  # If this is smaller than 1. it is an indicator that it does not fit its original dimensions
1060  self.log.infof("Pd sum before normalization: {:.3f}", sumPd)
1061  Pd /= sumPd
1062  # Convert Pd into a Psf object
1063  Pd = self.makeKernelPsfFromArraymakeKernelPsfFromArray(Pd)
1064 
1065  xy0 = self.cutBoxes1cutBoxes1.outerBox.getMin()
1066  # D is already converted back to dtype of fullExp1
1067  # Encapsulate D simply into an afwImage.Image for correct inner-outer box handling
1068  imgD = afwImage.Image(D, deep=False, xy0=xy0, dtype=self.fullExp1fullExp1.image.dtype)
1069  diffExp.image[self.cutBoxes1cutBoxes1.innerBox] = imgD[self.cutBoxes1cutBoxes1.innerBox]
1070  imgVarPlaneD = afwImage.Image(varPlaneD, deep=False, xy0=xy0,
1071  dtype=self.fullExp1fullExp1.variance.dtype)
1072  diffExp.variance[self.cutBoxes1cutBoxes1.innerBox] = imgVarPlaneD[self.cutBoxes1cutBoxes1.innerBox]
1073  diffExp.mask[self.cutBoxes1cutBoxes1.innerBox] = self.calculateMaskPlanecalculateMaskPlane(
1074  self.fullExp1fullExp1.mask[self.cutBoxes1cutBoxes1.innerBox], self.fullExp2fullExp2.mask[self.cutBoxes2cutBoxes2.innerBox])
1075 
1076  # Calibrate the image; subimages on the grid must be on the same photometric scale
1077  # Now the calibration object will be 1. everywhere
1078  diffExp.maskedImage[self.cutBoxes1cutBoxes1.innerBox] /= ftDiff.Fd
1079 
1080  if ftDiff.S is not None and scoreExp is not None:
1081  S = self.inverseFftAndCropImageinverseFftAndCropImage(
1082  ftDiff.S, self.freqSpaceShapefreqSpaceShape, dtype=self.fullExp1fullExp1.image.dtype)
1083  varPlaneS = self.inverseFftAndCropImageinverseFftAndCropImage(
1084  ftDiff.varPlaneS, self.freqSpaceShapefreqSpaceShape, dtype=self.fullExp1fullExp1.variance.dtype)
1085 
1086  imgS = afwImage.Image(S, deep=False, xy0=xy0, dtype=self.fullExp1fullExp1.image.dtype)
1087  imgVarPlaneS = afwImage.Image(varPlaneS, deep=False, xy0=xy0,
1088  dtype=self.fullExp1fullExp1.variance.dtype)
1089  scoreExp.image[self.cutBoxes1cutBoxes1.innerBox] = imgS[self.cutBoxes1cutBoxes1.innerBox]
1090  scoreExp.variance[self.cutBoxes1cutBoxes1.innerBox] = imgVarPlaneS[self.cutBoxes1cutBoxes1.innerBox]
1091 
1092  # PSF of S
1093  Ps = self.inverseFftAndCropImageinverseFftAndCropImage(ftDiff.Ps, self.psfShape1psfShape1, dtype=self.subExpPsf1subExpPsf1.dtype)
1094  sumPs = np.sum(Ps)
1095  self.log.infof("Ps sum before normalization: {:.3f}", sumPs)
1096  Ps /= sumPs
1097 
1098  # TODO DM-23855 : Additional score image corrections may be done here
1099 
1100  scoreExp.mask[self.cutBoxes1cutBoxes1.innerBox] = diffExp.mask[self.cutBoxes1cutBoxes1.innerBox]
1101  # Convert Ps into a Psf object
1102  Ps = self.makeKernelPsfFromArraymakeKernelPsfFromArray(Ps)
1103  else:
1104  Ps = None
1105  self.gridPsfsgridPsfs.append(pipeBase.Struct(bbox=self.cutBoxes1cutBoxes1.innerBox, Pd=Pd, Ps=Ps))
1106 
1107  def finishResultExposures(self, diffExp, scoreExp=None):
1108  """Perform final steps on the full difference exposure result.
1109 
1110  Set photometric calibration, psf properties of the exposures.
1111 
1112  Parameters
1113  ----------
1114  diffExp : `lsst.afw.image.Exposure`
1115  The result difference image exposure to finalize.
1116  scoreExp : `lsst.afw.image.Exposure` or `None`
1117  The result score exposure to finalize.
1118 
1119  Returns
1120  -------
1121  None.
1122  """
1123  # Set Calibration and PSF of the result exposures
1124  calibOne = afwImage.PhotoCalib(1.)
1125  diffExp.setPhotoCalib(calibOne)
1126  # Create the spatially varying PSF and set it for the diffExp
1127  # Set the PSF of this subExposure
1128  if len(self.gridPsfsgridPsfs) > 1:
1129  diffExp.setPsf(
1130  self.makeSpatialPsfmakeSpatialPsf(
1131  pipeBase.Struct(bbox=x.bbox, psf=x.Pd) for x in self.gridPsfsgridPsfs
1132  ))
1133  if scoreExp is not None:
1134  scoreExp.setPsf(
1135  self.makeSpatialPsfmakeSpatialPsf(
1136  pipeBase.Struct(bbox=x.bbox, psf=x.Ps) for x in self.gridPsfsgridPsfs
1137  ))
1138  else:
1139  # We did not have a grid, use the result psf without
1140  # making a CoaddPsf
1141  diffExp.setPsf(self.gridPsfsgridPsfs[0].Pd)
1142  if scoreExp is not None:
1143  scoreExp.setPsf(self.gridPsfsgridPsfs[0].Ps)
1144 
1145  # diffExp.setPsf(self.makeKernelPsfFromArray(Pd))
1146  if scoreExp is not None:
1147  scoreExp.setPhotoCalib(calibOne)
1148  # Zero score exposure where its variance is zero or the inputs are non-finite
1149  flt = (self.filtsImg1filtsImg1.filtInf | self.filtsImg2filtsImg2.filtInf
1150  | self.filtsImg1filtsImg1.filtNaN | self.filtsImg2filtsImg2.filtNaN
1151  | self.filtsVar1filtsVar1.filtInf | self.filtsVar2filtsVar2.filtInf
1152  | self.filtsVar1filtsVar1.filtNaN | self.filtsVar2filtsVar2.filtNaN)
1153  # Ensure that no division by 0 occurs in S/sigma(S).
1154  # S is set to be always finite, 0 where pixels non-finite
1155  tiny = np.finfo(scoreExp.variance.dtype).tiny * 100
1156  flt = np.logical_or(flt, scoreExp.variance.array < tiny)
1157  # Don't set variance to tiny.
1158  # It becomes 0 in case of conversion to single precision.
1159  # Set variance to 1, indicating that zero is in units of "sigmas" already.
1160  scoreExp.variance.array[flt] = 1
1161  scoreExp.image.array[flt] = 0
1162 
1163  def run(self, exposure1, exposure2, calculateScore=True):
1164  """Task entry point to perform the zogy subtraction
1165  of ``exposure1-exposure2``.
1166 
1167  Parameters
1168  ----------
1169  exposure1, exposure2 : `lsst.afw.image.Exposure`
1170  Two exposures warped and matched into matching pixel dimensions.
1171  calculateScore : `bool`, optional
1172  If True (default), calculate the score image and return in ``scoreExp``.
1173 
1174 
1175  Returns
1176  -------
1177  resultName : `lsst.pipe.base.Struct`
1178  - ``diffExp`` : `lsst.afw.image.Exposure`
1179  The Zogy difference exposure (``exposure1-exposure2``).
1180  - ``scoreExp`` : `lsst.afw.image.Exposure` or `None`
1181  The Zogy significance or score (S) exposure if ``calculateScore==True``.
1182  - ``ftDiff`` : `lsst.pipe.base.Struct`
1183  Lower level return struct by `calculateFourierDiffim` with added
1184  fields from the task instance. For debug purposes.
1185 
1186  Notes
1187  -----
1188 
1189  ``diffExp`` and ``scoreExp`` always inherit their metadata from
1190  ``exposure1`` (e.g. dtype, bbox, wcs).
1191 
1192  The score image (``S``) is defined in the ZOGY paper as the detection
1193  statistic value at each pixel. In the ZOGY image model, the input images
1194  have uniform variance noises and thus ``S`` has uniform per pixel
1195  variance (though it is not scaled to 1). In Section 3.3 of the paper,
1196  there are "corrections" defined to the score image to correct the
1197  significance values for some deviations from the image model. The first
1198  of these corrections is the calculation of the _variance plane_ of ``S``
1199  allowing for different per pixel variance values by following the
1200  overall convolution operation on the pixels of the input images. ``S``
1201  scaled (divided) by its corrected per pixel noise is referred as
1202  ``Scorr`` in the paper.
1203 
1204  In the current implementation, ``scoreExp`` contains ``S`` in its image
1205  plane and the calculated (non-uniform) variance plane of ``S`` in its
1206  variance plane. ``scoreExp`` can be used directly for source detection
1207  as a likelihood image by respecting its variance plane or can be divided
1208  by the square root of the variance plane to scale detection significance
1209  values into units of sigma. ``S`` should be interpreted as a detection
1210  likelihood directly on a per-pixel basis. The calculated PSF
1211  of ``S`` is merely an indication how much the input PSFs localize point
1212  sources.
1213 
1214  TODO DM-23855 : Implement further correction tags to the variance of
1215  ``scoreExp``. As of DM-25174 it is not determined how important these
1216  further correction tags are.
1217  """
1218  # We use the dimensions of the 1st image only in the code
1219  if exposure1.getDimensions() != exposure2.getDimensions():
1220  raise ValueError("Exposure dimensions do not match ({} != {} )".format(
1221  exposure1.getDimensions(), exposure2.getDimensions()))
1222 
1223  self.prepareFullExposureprepareFullExposure(exposure1, exposure2, correctBackground=self.config.correctBackground)
1224  # Do not use the exposure1, exposure2 input arguments from here
1225  exposure1 = None
1226  exposure2 = None
1227  if self.config.doSpatialGrid:
1228  gridBoxes = self.generateGridgenerateGrid(
1229  self.fullExp1fullExp1.getBBox(ImageOrigin.LOCAL), Extent2I(self.borderSizeborderSize, self.borderSizeborderSize),
1230  Extent2I(Extent2I(self.borderSizeborderSize, self.borderSizeborderSize) * self.config.gridInnerSize),
1231  powerOfTwo=True)
1232  else:
1233  gridBoxes = self.generateGridgenerateGrid(
1234  self.fullExp1fullExp1.getBBox(ImageOrigin.LOCAL), Extent2I(self.borderSizeborderSize, self.borderSizeborderSize),
1235  self.fullExp1fullExp1.getBBox().getDimensions(), powerOfTwo=True)
1236 
1237  diffExp = self.fullExp1fullExp1.clone()
1238  if calculateScore:
1239  scoreExp = self.fullExp1fullExp1.clone()
1240  else:
1241  scoreExp = None
1242  self.gridPsfsgridPsfs = []
1243  # Loop through grid boxes
1244  for boxPair in gridBoxes:
1245  self.prepareSubExposureprepareSubExposure(boxPair) # Extract sub images and fft
1246  ftDiff = self.calculateFourierDiffimcalculateFourierDiffim(
1247  self.psfFft1psfFft1, self.subImgFft1subImgFft1, self.subVarImgFft1subVarImgFft1, self.F1F1, self.subExpVar1subExpVar1,
1248  self.psfFft2psfFft2, self.subImgFft2subImgFft2, self.subVarImgFft2subVarImgFft2, self.F2F2, self.subExpVar2subExpVar2,
1249  calculateScore=calculateScore)
1250  self.pasteSubDiffImgpasteSubDiffImg(ftDiff, diffExp, scoreExp) # Paste back result
1251  self.finishResultExposuresfinishResultExposures(diffExp, scoreExp)
1252  # Add debug info from the task instance
1253  ftDiff.freqSpaceShape = self.freqSpaceShapefreqSpaceShape # The outer shape of the last grid cell
1254  ftDiff.psfShape1 = self.psfShape1psfShape1 # The psf image shape in exposure1
1255  ftDiff.psfShape2 = self.psfShape2psfShape2 # The psf image shape in exposure2
1256  ftDiff.borderSize = self.borderSizeborderSize # The requested padding around the inner region
1257  return pipeBase.Struct(diffExp=diffExp,
1258  scoreExp=scoreExp,
1259  ftDiff=ftDiff)
1260 
1261 
1263  """Config for the ZogyImagePsfMatchTask"""
1264 
1265  zogyConfig = pexConfig.ConfigField(
1266  dtype=ZogyConfig,
1267  doc='ZogyTask config to use',
1268  )
1269 
1270 
1272  """Task to perform Zogy PSF matching and image subtraction.
1273 
1274  This class inherits from ImagePsfMatchTask to contain the _warper
1275  subtask and related methods.
1276  """
1277 
1278  ConfigClass = ZogyImagePsfMatchConfig
1279 
1280  def __init__(self, *args, **kwargs):
1281  ImagePsfMatchTask.__init__(self, *args, **kwargs)
1282 
1283  def run(self, scienceExposure, templateExposure, doWarping=True):
1284  """Register, PSF-match, and subtract two Exposures, ``scienceExposure - templateExposure``
1285  using the ZOGY algorithm.
1286 
1287  Parameters
1288  ----------
1289  templateExposure : `lsst.afw.image.Exposure`
1290  exposure to be warped to scienceExposure.
1291  scienceExposure : `lsst.afw.image.Exposure`
1292  reference Exposure.
1293  doWarping : `bool`
1294  what to do if templateExposure's and scienceExposure's WCSs do not match:
1295  - if True then warp templateExposure to match scienceExposure
1296  - if False then raise an Exception
1297 
1298  Notes
1299  -----
1300  Do the following, in order:
1301  - Warp templateExposure to match scienceExposure, if their WCSs do not already match
1302  - Compute subtracted exposure ZOGY image subtraction algorithm on the two exposures
1303 
1304  This is the new entry point of the task as of DM-25115.
1305 
1306  Returns
1307  -------
1308  results : `lsst.pipe.base.Struct` containing these fields:
1309  - subtractedExposure: `lsst.afw.image.Exposure`
1310  The subtraction result.
1311  - warpedExposure: `lsst.afw.image.Exposure` or `None`
1312  templateExposure after warping to match scienceExposure
1313  """
1314 
1315  if not self._validateWcs_validateWcs(scienceExposure, templateExposure):
1316  if doWarping:
1317  self.log.info("Warping templateExposure to scienceExposure")
1318  xyTransform = afwGeom.makeWcsPairTransform(templateExposure.getWcs(),
1319  scienceExposure.getWcs())
1320  psfWarped = measAlg.WarpedPsf(templateExposure.getPsf(), xyTransform)
1321  templateExposure = self._warper_warper.warpExposure(
1322  scienceExposure.getWcs(), templateExposure, destBBox=scienceExposure.getBBox())
1323  templateExposure.setPsf(psfWarped)
1324  else:
1325  raise RuntimeError("Input images are not registered. Consider setting doWarping=True.")
1326 
1327  config = self.config.zogyConfig
1328  task = ZogyTask(config=config)
1329  results = task.run(scienceExposure, templateExposure)
1330  results.warpedExposure = templateExposure
1331  return results
1332 
1333  def subtractExposures(self, templateExposure, scienceExposure, *args):
1334  raise NotImplementedError
1335 
1336  def subtractMaskedImages(self, templateExposure, scienceExposure, *args):
1337  raise NotImplementedError
1338 
1339 
1340 subtractAlgorithmRegistry.register('zogy', ZogyImagePsfMatchTask)
def _validateWcs(self, templateExposure, scienceExposure)
def subtractMaskedImages(self, templateExposure, scienceExposure, *args)
Definition: zogy.py:1336
def __init__(self, *args, **kwargs)
Definition: zogy.py:1280
def run(self, scienceExposure, templateExposure, doWarping=True)
Definition: zogy.py:1283
def subtractExposures(self, templateExposure, scienceExposure, *args)
Definition: zogy.py:1333
def generateGrid(imageBox, minEdgeDims, innerBoxDims, minTotalDims=None, powerOfTwo=False)
Definition: zogy.py:341
def padCenterOriginArray(A, newShape, useInverse=False, dtype=None)
Definition: zogy.py:142
def subtractImageMean(image, mask, statsControl)
Definition: zogy.py:584
def calculateFourierDiffim(self, psf1, im1, varPlane1, F1, varMean1, psf2, im2, varPlane2, F2, varMean2, calculateScore=True)
Definition: zogy.py:860
def calculateMaskPlane(mask1, mask2, effPsf1=None, effPsf2=None)
Definition: zogy.py:970
def run(self, exposure1, exposure2, calculateScore=True)
Definition: zogy.py:1163
def prepareFullExposure(self, exposure1, exposure2, correctBackground=False)
Definition: zogy.py:614
def estimateMatchingKernelSize(psf1, psf2)
Definition: zogy.py:269
def removeNonFinitePixels(self, imgArr)
Definition: zogy.py:501
def splitBorder(innerBox, outerBox)
Definition: zogy.py:298
def initializeSubImage(self, fullExp, innerBox, outerBox, noiseMeanVar, useNoise=True)
Definition: zogy.py:210
def makeSpatialPsf(self, gridPsfs)
Definition: zogy.py:437
def padAndFftImage(self, imgArr)
Definition: zogy.py:470
def checkCentroids(self, psfArr1, psfArr2)
Definition: zogy.py:834
def pasteSubDiffImg(self, ftDiff, diffExp, scoreExp=None)
Definition: zogy.py:1017
def finishResultExposures(self, diffExp, scoreExp=None)
Definition: zogy.py:1107
def prepareSubExposure(self, localCutout, psf1=None, psf2=None, sig1=None, sig2=None)
Definition: zogy.py:690
def computePsfAtCenter(exposure)
Definition: zogy.py:562
def inverseFftAndCropImage(self, imgArr, origSize, filtInf=None, filtNaN=None, dtype=None)
Definition: zogy.py:528
def _computeVarianceMean(self, exposure)
Definition: zogy.py:132
std::shared_ptr< TransformPoint2ToPoint2 > makeWcsPairTransform(SkyWcs const &src, SkyWcs const &dst)
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())