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