25from lsst.geom import Box2I, Point2I, Extent2I
35from .imagePsfMatch
import (ImagePsfMatchTask, ImagePsfMatchConfig,
36 subtractAlgorithmRegistry)
38__all__ = [
"ZogyTask",
"ZogyConfig",
39 "ZogyImagePsfMatchConfig",
"ZogyImagePsfMatchTask"]
42"""Tasks for performing the "Proper image subtraction" algorithm of
43Zackay, et al. (2016), hereafter simply referred to as 'ZOGY (2016)'.
45`ZogyTask` contains methods to perform the basic estimation of the
46ZOGY diffim ``D``, its updated PSF, and the variance-normalized
47likelihood image ``S_corr``. We have implemented ZOGY using the
48proscribed methodology, computing all convolutions in Fourier space,
49and also variants in which the convolutions are performed in real
50(image) space. The former is faster and results in fewer artifacts
51when the PSFs are noisy (i.e., measured, for example, via
52`PsfEx`). The latter is presumed to be preferred as it can account for
53masks correctly with fewer "ringing" artifacts from edge effects or
54saturated stars, but noisy PSFs result in their own smaller
55artifacts. Removal of these artifacts is a subject of continuing
56research. Currently, we "pad" the PSFs when performing the
57subtractions in real space, which reduces, but does not entirely
58eliminate these artifacts.
60All methods in `ZogyTask` assume template and science images are
61already accurately photometrically and astrometrically registered.
63`ZogyMapper` is a wrapper which runs `ZogyTask` in the
64`ImageMapReduce` framework, computing of ZOGY diffim's on small,
65overlapping sub-images, thereby enabling complete ZOGY diffim's which
66account for spatially-varying noise and PSFs across the two input
67exposures. An example of the use of this task is in the `testZogy.py`
73 """Configuration parameters for the ZogyTask
76 templateFluxScaling = pexConfig.Field(
79 doc="Template flux scaling factor (Fr in ZOGY paper)"
82 scienceFluxScaling = pexConfig.Field(
85 doc=
"Science flux scaling factor (Fn in ZOGY paper)"
88 scaleByCalibration = pexConfig.Field(
91 doc=
"Compute the flux normalization scaling based on the image calibration."
92 "This overrides 'templateFluxScaling' and 'scienceFluxScaling'."
95 correctBackground = pexConfig.Field(
98 doc=
"Subtract exposure background mean to have zero expectation value."
101 ignoreMaskPlanes = pexConfig.ListField(
103 default=(
"INTRP",
"EDGE",
"DETECTED",
"SAT",
"CR",
"BAD",
"NO_DATA",
"DETECTED_NEGATIVE"),
104 doc=
"Mask planes to ignore for statistics"
106 maxPsfCentroidDist = pexConfig.Field(
109 doc=
"Maximum centroid difference allowed between the two exposure PSFs (pixels)."
111 doSpatialGrid = pexConfig.Field(
114 doc=
"Split the exposure and perform matching with the spatially varying PSF."
116 gridInnerSize = pexConfig.Field(
119 doc=
"Approximate useful inner size of the grid cells in units of the "
120 "estimated matching kernel size (doSpatialGrid=True only)."
125 """Task to perform ZOGY proper image subtraction. See module-level documentation for
129 ConfigClass = ZogyConfig
130 _DefaultName = "imageDifferenceZogy"
132 def _computeVarianceMean(self, exposure):
133 """Compute the sigma-clipped mean of the variance image of ``exposure``.
136 exposure.getMaskedImage().getMask(),
138 var = statObj.getValue(afwMath.MEANCLIP)
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.
146 Implement also the inverse operation, crop the padding
and re-center data.
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``.
167 The padded
or unpadded array
with shape of `newShape`
and dtype of ``dtype``.
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.
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.
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):
191 raise ValueError(
"Newshape dimensions must be greater or equal")
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):
198 raise ValueError(
"Newshape dimensions must be smaller or equal")
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]:]
211 """Initializes a sub image.
216 The full exposure to cut sub image from.
218 The useful area of the calculation up to the whole bounding box of
219 ``fullExp``. ``fullExp`` must contain this box.
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
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.
233 result : `lsst.pipe.base.Struct`
234 - ``subImg``, ``subVarImg`` : `lsst.afw.image.ImageD`
235 The new sub image
and its sub variance plane.
239 ``innerBox``, ``outerBox`` must be
in the PARENT system of ``fullExp``.
241 Supports the non-grid option when ``innerBox`` equals to the
242 bounding box of ``fullExp``.
244 fullBox = fullExp.getBBox()
245 subImg = afwImage.ImageD(outerBox, 0)
246 subVarImg = afwImage.ImageD(outerBox, noiseMeanVar)
247 borderBoxes = self.splitBordersplitBorder(innerBox, outerBox)
250 rng = np.random.Generator(
251 np.random.PCG64(seed=np.array([noiseMeanVar]).view(int)))
252 noiseSig = np.sqrt(noiseMeanVar)
253 for box
in borderBoxes:
254 if not fullBox.contains(box):
255 R = subImg[box].array
256 R[...] = rng.normal(scale=noiseSig, size=R.shape)
258 subImg[innerBox].array[...] = fullExp.image[innerBox].array
259 subVarImg[innerBox].array[...] = fullExp.variance[innerBox].array
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)
270 """Estimate the image space size of the matching kernels.
272 Return ten times the larger Gaussian sigma estimate but at least
273 the largest of the original psf dimensions.
278 The PSFs of the two input exposures.
283 Conservative estimate for matching kernel size
in pixels.
284 This
is the minimum padding around the inner region at each side.
289 sig1 = psf1.computeShape(psf1.getAveragePosition()).getDeterminantRadius()
290 sig2 = psf2.computeShape(psf2.getAveragePosition()).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())
299 """Split the border area around the inner box into 8 disjunct boxes.
306 The outer box. It must be at least 1 pixel larger in each direction than the inner box.
310 resultBoxes : `list` of 8 boxes covering the edge around innerBox
314 The border boxes do
not overlap. The border
is covered counter clockwise
315 starting
from lower left corner.
319 ValueError : If ``outerBox``
is not larger than ``innerBox``.
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")
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)
334 pointPairs = ((o1, i1), (i1 +
Extent2I(1, 0), p2 +
Extent2I(-1, 0)), (o2, i2),
338 return [
Box2I(x, y, invert=
True)
for (x, y)
in pointPairs]
341 def generateGrid(imageBox, minEdgeDims, innerBoxDims, minTotalDims=None, powerOfTwo=False):
342 """Generate a splitting grid for an image.
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.
350 Bounding box of the exposure to split.
352 Minimum edge width
in (x,y) directions each side.
354 Minimum requested inner box dimensions (x,y).
355 The actual dimensions can be larger due to rounding.
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.
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.
368 See diffimTests/tickets/DM-28928_spatial_grid notebooks
for demonstration of this code.
370 This method can be used
for both PARENT
and LOCAL bounding boxes.
372 The outerBox dimensions are always even.
376 boxList : `list` of `lsst.pipe.base.Struct`
377 ``innerBox``, ``outerBox`` : `
lsst.geom.Box2I`, inner boxes
and overlapping border around them.
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
385 innerWidth = width // nX
389 xCorners = np.zeros(nX + 1)
390 xCorners[:-1] = np.arange(nX)*innerWidth + imageBox.minX
391 xCorners[-1] = imageBox.endX
393 nY = height // innerBoxDims.y
395 innerHeight = height // nY
399 yCorners = np.zeros(nY + 1)
400 yCorners[:-1] = np.arange(nY)*innerHeight + imageBox.minY
401 yCorners[-1] = imageBox.endY
405 for i_y
in range(nY):
406 for i_x
in range(nX):
408 Point2I(xCorners[i_x + 1] - 1, yCorners[i_y + 1] - 1))
410 paddedWidth = innerBox.width + doubleEdgeDims.x
411 if minTotalDims
is not None and paddedWidth < minTotalDims.width:
412 paddedWidth = minTotalDims.width
414 i2x = np.searchsorted(powersOf2, paddedWidth, side=
'left')
415 if i2x < len(powersOf2):
416 paddedWidth = powersOf2[i2x]
417 if paddedWidth % 2 == 1:
420 totalXedge = paddedWidth - innerBox.width
422 paddedHeight = innerBox.height + doubleEdgeDims.y
423 if minTotalDims
is not None and paddedHeight < minTotalDims.height:
424 paddedHeight = minTotalDims.height
426 i2y = np.searchsorted(powersOf2, paddedHeight, side=
'left')
427 if i2y < len(powersOf2):
428 paddedHeight = powersOf2[i2y]
429 if paddedHeight % 2 == 1:
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))
438 """Construct a CoaddPsf based on PSFs from individual sub image solutions.
442 gridPsfs : iterable of `lsst.pipe.base.Struct`
443 Iterable of bounding boxes (``bbox``) and Psf solutions (``psf``).
448 A psf constructed
from the PSFs of the individual subExposures.
450 schema = afwTable.ExposureTable.makeMinimalSchema()
451 schema.addField("weight", type=
"D", doc=
"Coadd weight")
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
464 mycatalog.append(record)
467 psf = measAlg.CoaddPsf(mycatalog, wcsref,
'weight')
471 """Prepare and forward FFT an image array.
475 imgArr : `numpy.ndarray` of `float`
476 Original array. In-place modified as `numpy.nan`
and `numpy.inf` are replaced by
481 result : `lsst.pipe.base.Struct`
482 - ``imFft`` : `numpy.ndarray` of `numpy.complex`.
484 - ``filtInf``, ``filtNaN`` : `numpy.ndarray` of `bool`
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`.
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 %s Inf and %s NaN values.",
496 np.sum(filtInf), np.sum(filtNaN))
498 imgArr = np.fft.fft2(imgArr)
499 return pipeBase.Struct(imFft=imgArr, filtInf=filtInf, filtNaN=filtNaN)
502 """Replace non-finite pixel values in-place.
504 Save the locations of non-finite values for restoration,
and replace them
505 with image mean values.
509 imgArr : `numpy.ndarray` of `float`
510 The image array. Non-finite values are replaced
in-place
in this array.
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.
518 filtInf = np.isinf(imgArr)
519 filtNaN = np.isnan(imgArr)
522 imgArr[filtInf] = np.nan
523 imgArr[filtInf | filtNaN] = np.nanmean(imgArr)
524 self.log.
debug(
"Replacing %s Inf and %s NaN values.",
525 np.sum(filtInf), np.sum(filtNaN))
526 return pipeBase.Struct(filtInf=filtInf, filtNaN=filtNaN)
529 """Inverse FFT and crop padding from image array.
533 imgArr : `numpy.ndarray` of `numpy.complex`
534 Fourier space array representing a real image.
536 origSize : `tuple` of `int`
537 Original unpadded shape tuple of the image to be cropped to.
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.
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`.
550 result : `numpy.ndarray` of `dtype`
552 imgNew = np.fft.ifft2(imgArr)
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
563 """Computes the PSF image at the bbox center point.
565 This may be at a fractional pixel position.
575 Calculated psf image.
577 pBox = exposure.getBBox()
578 cen = pBox.getCenter()
579 psf = exposure.getPsf()
580 psfImg = psf.computeKernelImage(cen)
585 """In-place subtraction of sigma-clipped mean of the image.
590 Image to manipulate. Its sigma clipped mean is in-place subtracted.
593 Mask to use
for ignoring pixels.
596 Config of sigma clipped mean statistics calculation.
604 ValueError : If image mean
is nan.
607 afwMath.MEANCLIP, statsControl)
608 mean = statObj.getValue(afwMath.MEANCLIP)
609 if not np.isnan(mean):
612 raise ValueError(
"Image mean is NaN.")
615 """Performs calculations that apply to the full exposures once only.
621 The input exposures. Copies are made for internal calculations.
623 correctBackground : `bool`, optional
624 If
True, subtracts sigma-clipped mean of exposures. The algorithm
625 assumes zero expectation value at background pixels.
633 Set a number of instance fields
with pre-calculated values.
637 ValueError : If photometric calibrations are
not available
while
638 ``config.scaleByCalibration`` equals
True.
643 self.statsControlstatsControl.setAndMask(afwImage.Mask.getPlaneBitMask(
644 self.config.ignoreMaskPlanes))
646 exposure1 = exposure1.clone()
647 exposure2 = exposure2.clone()
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)
665 self.
F1F1 = self.config.templateFluxScaling
666 self.
F2F2 = self.config.scienceFluxScaling
667 mImg1 = exposure1.maskedImage
668 mImg2 = exposure2.maskedImage
671 if correctBackground:
677 self.log.
debug(
"Minimum padding border size: %s pixels", self.
borderSizeborderSize)
684 exposure1.maskedImage = mImg1
685 exposure2.maskedImage = mImg2
691 """Perform per-sub exposure preparations.
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`
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.
710 self.log.debug("Processing LOCAL cell w/ inner box:%s, outer box:%s",
711 localCutout.innerBox, localCutout.outerBox)
714 innerBox=localCutout.innerBox.shiftedBy(
Extent2I(self.
fullExp1fullExp1.getXY0())),
715 outerBox=localCutout.outerBox.shiftedBy(
Extent2I(self.
fullExp1fullExp1.getXY0())))
717 innerBox=localCutout.innerBox.shiftedBy(
Extent2I(self.
fullExp2fullExp2.getXY0())),
718 outerBox=localCutout.outerBox.shiftedBy(
Extent2I(self.
fullExp2fullExp2.getXY0())))
731 psfBBox1 = self.
subExpPsf1subExpPsf1.getBBox()
732 psfBBox2 = self.
subExpPsf2subExpPsf2.getBBox()
733 self.
psfShape1psfShape1 = (psfBBox1.getHeight(), psfBBox1.getWidth())
734 self.
psfShape2psfShape2 = (psfBBox2.getHeight(), psfBBox2.getWidth())
748 self.
freqSpaceShapefreqSpaceShape = (localCutout.outerBox.getHeight(), localCutout.outerBox.getWidth())
774 """Square the argument in pixel space.
778 D : 2D `numpy.ndarray` of `numpy.complex`
779 Fourier transform of a real valued array.
783 R : `numpy.ndarray` of `numpy.complex`
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.
792 R = np.real(np.fft.ifft2(D))
799 """Calculate the centroid coordinates of a 2D array.
803 A : 2D `numpy.ndarray` of `float`
804 The input array. Must not be all exact zero.
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.
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.
818 ycen, xcen : `tuple` of `float`
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
833 """Check whether two PSF array centroids' distance is within tolerance.
837 psfArr1, psfArr2 : `numpy.ndarray` of `float`
838 Input PSF arrays to check.
847 Centroid distance exceeds `config.maxPsfCentroidDist` pixels.
853 if dy*dy + dx*dx > self.config.maxPsfCentroidDist*self.config.maxPsfCentroidDist:
855 f
"PSF centroids are offset by more than {self.config.maxPsfCentroidDist:.2f} pixels.")
858 psf2, im2, varPlane2, F2, varMean2, calculateScore=True):
859 """Convolve and subtract two images in Fourier space.
861 Calculate the ZOGY proper difference image, score image and their PSFs.
862 All input
and output arrays are
in Fourier space.
866 psf1, psf2 : `numpy.ndarray`, (``self.
freqSpaceShapefreqSpaceShape``,)
867 Psf arrays. Must be already
in Fourier space.
868 im1, im2 : `numpy.ndarray`, (``self.
freqSpaceShapefreqSpaceShape``,)
869 Image arrays. Must be already
in Fourier space.
870 varPlane1, varPlane2 : `numpy.ndarray`, (``self.
freqSpaceShapefreqSpaceShape``,)
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`.
883 result : `pipe.base.Struct`
884 All arrays are
in Fourier space
and have shape ``self.
freqSpaceShapefreqSpaceShape``.
887 Photometric level of ``D`` (`float`).
889 The difference image (`numpy.ndarray` [`numpy.complex`]).
891 Variance plane of ``D`` (`numpy.ndarray` [`numpy.complex`]).
893 PSF of ``D`` (`numpy.ndarray` [`numpy.complex`]).
895 Significance (score) image (`numpy.ndarray` [`numpy.complex`]
or `
None`).
897 Variance plane of ``S`` ((`numpy.ndarray` [`numpy.complex`]
or `
None`).
899 PSF of ``S`` (`numpy.ndarray` [`numpy.complex`]).
903 All array inputs
and outputs are Fourier-space images
with shape of
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.
910 var1F2Sq = varMean1*F2*F2
911 var2F1Sq = varMean2*F1*F1
913 psfAbsSq1 = np.real(np.conj(psf1)*psf1)
914 psfAbsSq2 = np.real(np.conj(psf2)*psf2)
915 FdDenom = np.sqrt(var1F2Sq + var2F1Sq)
918 tiny = np.finfo(psf1.dtype).tiny * 100
919 sDenom = var1F2Sq*psfAbsSq2 + var2F1Sq*psfAbsSq1
929 fltZero = sDenom < tiny
930 nZero = np.sum(fltZero)
931 self.log.
debug(
"There are %s frequencies where both FFTd PSFs are close to zero.", nZero)
934 fltZero = np.nonzero(fltZero)
935 sDenom[fltZero] = tiny
936 denom = np.sqrt(sDenom)
941 c1[fltZero] = F2/FdDenom
942 c2[fltZero] = F1/FdDenom
946 Pd = FdDenom*psf1*psf2/denom
952 c1 = F1*F2*F2*np.conj(psf1)*psfAbsSq2/sDenom
953 c2 = F2*F1*F1*np.conj(psf2)*psfAbsSq1/sDenom
964 return pipeBase.Struct(D=D, Pd=Pd, varPlaneD=varPlaneD, Fd=Fd,
965 S=S, Ps=Ps, varPlaneS=varPlaneS)
969 """Calculate the mask plane of the difference image.
974 Mask planes of the two exposures.
980 Mask plane for the subtraction result.
984 TODO DM-25174 : Specification of effPsf1, effPsf2 are
not yet supported.
988 if effPsf1
is not None or effPsf2
is not None:
992 raise NotImplementedError(
"Mask plane only 'convolution' operation is not yet supported")
999 """Create a non spatially varying PSF from a `numpy.ndarray`.
1004 2D array to use as the new psf image. The pixels are copied.
1009 The constructed PSF.
1011 psfImg = afwImage.ImageD(A.astype(np.float64, copy=True), deep=
False)
1016 """Paste sub image results back into result Exposure objects.
1020 ftDiff : `lsst.pipe.base.Struct`
1021 Result struct by `calculateFourierDiffim`.
1023 The result exposure to paste into the sub image result.
1024 Must be dimensions and dtype of ``self.
fullExp1fullExp1``.
1026 The result score exposure to paste into the sub image result.
1027 Must be dimensions
and dtype of ``self.
fullExp1fullExp1``.
1028 If `
None`, the score image results are disregarded.
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.
1039 Add an entry to the ``self.
gridPsfsgridPsfs`` list.
1041 gridPsfs : `list` of `lsst.pipe.base.Struct`
1043 The inner region of the grid cell.
1045 The diffim PSF
in this cell.
1047 The score image PSF
in this cell
or `
None`
if the score
1048 image was
not calculated.
1058 self.log.
info(
"Pd sum before normalization: %.3f", sumPd)
1063 xy0 = self.
cutBoxes1cutBoxes1.outerBox.getMin()
1067 diffExp.image[self.
cutBoxes1cutBoxes1.innerBox] = imgD[self.
cutBoxes1cutBoxes1.innerBox]
1069 dtype=self.
fullExp1fullExp1.variance.dtype)
1070 diffExp.variance[self.
cutBoxes1cutBoxes1.innerBox] = imgVarPlaneD[self.
cutBoxes1cutBoxes1.innerBox]
1076 diffExp.maskedImage[self.
cutBoxes1cutBoxes1.innerBox] /= ftDiff.Fd
1078 if ftDiff.S
is not None and scoreExp
is not None:
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]
1093 self.log.
info(
"Ps sum before normalization: %.3f", sumPs)
1098 scoreExp.mask[self.
cutBoxes1cutBoxes1.innerBox] = diffExp.mask[self.
cutBoxes1cutBoxes1.innerBox]
1106 """Perform final steps on the full difference exposure result.
1108 Set photometric calibration, psf properties of the exposures.
1113 The result difference image exposure to finalize.
1115 The result score exposure to finalize.
1123 diffExp.setPhotoCalib(calibOne)
1129 pipeBase.Struct(bbox=x.bbox, psf=x.Pd)
for x
in self.
gridPsfsgridPsfs
1131 if scoreExp
is not None:
1134 pipeBase.Struct(bbox=x.bbox, psf=x.Ps)
for x
in self.
gridPsfsgridPsfs
1139 diffExp.setPsf(self.
gridPsfsgridPsfs[0].Pd)
1140 if scoreExp
is not None:
1141 scoreExp.setPsf(self.
gridPsfsgridPsfs[0].Ps)
1144 if scoreExp
is not None:
1145 scoreExp.setPhotoCalib(calibOne)
1153 tiny = np.finfo(scoreExp.variance.dtype).tiny * 100
1154 flt = np.logical_or(flt, scoreExp.variance.array < tiny)
1158 scoreExp.variance.array[flt] = 1
1159 scoreExp.image.array[flt] = 0
1161 def run(self, exposure1, exposure2, calculateScore=True):
1162 """Task entry point to perform the zogy subtraction
1163 of ``exposure1-exposure2``.
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``.
1175 resultName : `lsst.pipe.base.Struct`
1177 The Zogy difference exposure (``exposure1-exposure2``).
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.
1187 ``diffExp``
and ``scoreExp`` always inherit their metadata
from
1188 ``exposure1`` (e.g. dtype, bbox, wcs).
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.
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
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.
1217 if exposure1.getDimensions() != exposure2.getDimensions():
1218 raise ValueError(
"Exposure dimensions do not match ({} != {} )".
format(
1219 exposure1.getDimensions(), exposure2.getDimensions()))
1221 self.
prepareFullExposureprepareFullExposure(exposure1, exposure2, correctBackground=self.config.correctBackground)
1225 if self.config.doSpatialGrid:
1233 self.
fullExp1fullExp1.getBBox().getDimensions(), powerOfTwo=
True)
1242 for boxPair
in gridBoxes:
1247 calculateScore=calculateScore)
1252 ftDiff.psfShape1 = self.
psfShape1psfShape1
1253 ftDiff.psfShape2 = self.
psfShape2psfShape2
1254 ftDiff.borderSize = self.
borderSizeborderSize
1255 return pipeBase.Struct(diffExp=diffExp,
1261 """Config for the ZogyImagePsfMatchTask"""
1263 zogyConfig = pexConfig.ConfigField(
1265 doc=
'ZogyTask config to use',
1270 """Task to perform Zogy PSF matching and image subtraction.
1273 subtask
and related methods.
1276 ConfigClass = ZogyImagePsfMatchConfig
1279 ImagePsfMatchTask.__init__(self, *args, **kwargs)
1281 def run(self, scienceExposure, templateExposure, doWarping=True):
1282 """Register, PSF-match, and subtract two Exposures, ``scienceExposure - templateExposure``
1283 using the ZOGY algorithm.
1288 exposure to be warped to scienceExposure.
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
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
1302 This
is the new entry point of the task
as of DM-25115.
1306 results : `lsst.pipe.base.Struct` containing these fields:
1308 The subtraction result.
1310 templateExposure after warping to match scienceExposure
1313 if not self.
_validateWcs_validateWcs(scienceExposure, templateExposure):
1315 self.log.
info(
"Warping templateExposure to scienceExposure")
1317 scienceExposure.getWcs())
1318 psfWarped = measAlg.WarpedPsf(templateExposure.getPsf(), xyTransform)
1320 scienceExposure.getWcs(), templateExposure, destBBox=scienceExposure.getBBox())
1321 templateExposure.setPsf(psfWarped)
1323 raise RuntimeError(
"Input images are not registered. Consider setting doWarping=True.")
1325 config = self.config.zogyConfig
1327 results = task.run(scienceExposure, templateExposure)
1328 results.warpedExposure = templateExposure
1332 raise NotImplementedError
1335 raise NotImplementedError
1338subtractAlgorithmRegistry.register(
'zogy', ZogyImagePsfMatchTask)
A polymorphic base class for representing an image's Point Spread Function.
A class to contain the data, WCS, and other information needed to describe an image of the sky.
A class to represent a 2-dimensional array of pixels.
Represent a 2-dimensional array of bitmask pixels.
The photometric calibration of an exposure.
A kernel created from an Image.
Pass parameters to a Statistics object.
Custom catalog class for ExposureRecord/Table.
An integer coordinate rectangle.
def _validateWcs(self, templateExposure, scienceExposure)
def subtractMaskedImages(self, templateExposure, scienceExposure, *args)
def __init__(self, *args, **kwargs)
def run(self, scienceExposure, templateExposure, doWarping=True)
def subtractExposures(self, templateExposure, scienceExposure, *args)
def generateGrid(imageBox, minEdgeDims, innerBoxDims, minTotalDims=None, powerOfTwo=False)
def padCenterOriginArray(A, newShape, useInverse=False, dtype=None)
def subtractImageMean(image, mask, statsControl)
def calculateFourierDiffim(self, psf1, im1, varPlane1, F1, varMean1, psf2, im2, varPlane2, F2, varMean2, calculateScore=True)
def makeKernelPsfFromArray(A)
def calculateMaskPlane(mask1, mask2, effPsf1=None, effPsf2=None)
def run(self, exposure1, exposure2, calculateScore=True)
def prepareFullExposure(self, exposure1, exposure2, correctBackground=False)
def estimateMatchingKernelSize(psf1, psf2)
def removeNonFinitePixels(self, imgArr)
def splitBorder(innerBox, outerBox)
def initializeSubImage(self, fullExp, innerBox, outerBox, noiseMeanVar, useNoise=True)
def makeSpatialPsf(self, gridPsfs)
def padAndFftImage(self, imgArr)
def checkCentroids(self, psfArr1, psfArr2)
def pasteSubDiffImg(self, ftDiff, diffExp, scoreExp=None)
def finishResultExposures(self, diffExp, scoreExp=None)
def prepareSubExposure(self, localCutout, psf1=None, psf2=None, sig1=None, sig2=None)
def computePsfAtCenter(exposure)
def inverseFftAndCropImage(self, imgArr, origSize, filtInf=None, filtNaN=None, dtype=None)
def _computeVarianceMean(self, exposure)
CoaddPsf is the Psf derived to be used for non-PSF-matched Coadd images.
A Psf defined by a Kernel.
std::shared_ptr< FrameSet > append(FrameSet const &first, FrameSet const &second)
Construct a FrameSet that performs two transformations in series.
std::shared_ptr< TransformPoint2ToPoint2 > makeWcsPairTransform(SkyWcs const &src, SkyWcs const &dst)
A Transform obtained by putting two SkyWcs objects "back to back".
Backwards-compatibility support for depersisting the old Calib (FluxMag0/FluxMag0Err) objects.
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())
Handle a watered-down front-end to the constructor (no variance)
int warpExposure(DestExposureT &destExposure, SrcExposureT const &srcExposure, WarpingControl const &control, typename DestExposureT::MaskedImageT::SinglePixel padValue=lsst::afw::math::edgePixel< typename DestExposureT::MaskedImageT >(typename lsst::afw::image::detail::image_traits< typename DestExposureT::MaskedImageT >::image_category()))
Warp (remap) one exposure to another.
Extent< int, 2 > Extent2I
Fit spatial kernel using approximate fluxes for candidates, and solving a linear system of equations.
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)