LSST Applications g070148d5b3+33e5256705,g0d53e28543+25c8b88941,g0da5cf3356+2dd1178308,g1081da9e2a+62d12e78cb,g17e5ecfddb+7e422d6136,g1c76d35bf8+ede3a706f7,g295839609d+225697d880,g2e2c1a68ba+cc1f6f037e,g2ffcdf413f+853cd4dcde,g38293774b4+62d12e78cb,g3b44f30a73+d953f1ac34,g48ccf36440+885b902d19,g4b2f1765b6+7dedbde6d2,g5320a0a9f6+0c5d6105b6,g56b687f8c9+ede3a706f7,g5c4744a4d9+ef6ac23297,g5ffd174ac0+0c5d6105b6,g6075d09f38+66af417445,g667d525e37+2ced63db88,g670421136f+2ced63db88,g71f27ac40c+2ced63db88,g774830318a+463cbe8d1f,g7876bc68e5+1d137996f1,g7985c39107+62d12e78cb,g7fdac2220c+0fd8241c05,g96f01af41f+368e6903a7,g9ca82378b8+2ced63db88,g9d27549199+ef6ac23297,gabe93b2c52+e3573e3735,gb065e2a02a+3dfbe639da,gbc3249ced9+0c5d6105b6,gbec6a3398f+0c5d6105b6,gc9534b9d65+35b9f25267,gd01420fc67+0c5d6105b6,geee7ff78d7+a14128c129,gf63283c776+ede3a706f7,gfed783d017+0c5d6105b6,w.2022.47
LSST Data Management Base Package
Loading...
Searching...
No Matches
maskStreaks.py
Go to the documentation of this file.
1# This file is part of pipe_tasks.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22__all__ = ["MaskStreaksConfig", "MaskStreaksTask", "setDetectionMask"]
23
24import lsst.pex.config as pexConfig
25import lsst.pipe.base as pipeBase
26import lsst.kht
27from lsst.utils.timer import timeMethod
28
29import numpy as np
30import scipy
31import textwrap
32import copy
33from skimage.feature import canny
34from sklearn.cluster import KMeans
35import warnings
36from dataclasses import dataclass
37
38
39def setDetectionMask(maskedImage, forceSlowBin=False, binning=None, detectedPlane="DETECTED",
40 badMaskPlanes=("NO_DATA", "INTRP", "BAD", "SAT", "EDGE"), detectionThreshold=5):
41 """Make detection mask and set the mask plane.
42
43 Creat a binary image from a masked image by setting all data with signal-to-
44 noise below some threshold to zero, and all data above the threshold to one.
45 If the binning parameter has been set, this procedure will be preceded by a
46 weighted binning of the data in order to smooth the result, after which the
47 result is scaled back to the original dimensions. Set the detection mask
48 plane with this binary image.
49
50 Parameters
51 ----------
52 maskedImage : `lsst.afw.image.maskedImage`
53 Image to be (optionally) binned and converted.
54 forceSlowBin : `bool`, optional
55 Force usage of slower binning method to check that the two methods
56 give the same result.
57 binning : `int`, optional
58 Number of pixels by which to bin image.
59 detectedPlane : `str`, optional
60 Name of mask with pixels that were detected above threshold in image.
61 badMaskPlanes : `set`, optional
62 Names of masks with pixels that are rejected.
63 detectionThreshold : `float`, optional
64 Boundary in signal-to-noise between non-detections and detections for
65 making a binary image from the original input image.
66 """
67 data = maskedImage.image.array
68 weights = 1 / maskedImage.variance.array
69 mask = maskedImage.getMask()
70
71 detectionMask = ((mask.array & mask.getPlaneBitMask(detectedPlane)))
72 badPixelMask = mask.getPlaneBitMask(badMaskPlanes)
73 badMask = (mask.array & badPixelMask) > 0
74 fitMask = detectionMask.astype(bool) & ~badMask
75
76 fitData = np.copy(data)
77 fitData[~fitMask] = 0
78 fitWeights = np.copy(weights)
79 fitWeights[~fitMask] = 0
80
81 if binning:
82 # Do weighted binning:
83 ymax, xmax = fitData.shape
84 if (ymax % binning == 0) and (xmax % binning == 0) and (not forceSlowBin):
85 # Faster binning method
86 binNumeratorReshape = (fitData * fitWeights).reshape(ymax // binning, binning,
87 xmax // binning, binning)
88 binDenominatorReshape = fitWeights.reshape(binNumeratorReshape.shape)
89 binnedNumerator = binNumeratorReshape.sum(axis=3).sum(axis=1)
90 binnedDenominator = binDenominatorReshape.sum(axis=3).sum(axis=1)
91 else:
92 # Slower binning method when (image shape mod binsize) != 0
93 warnings.warn('Using slow binning method--consider choosing a binsize that evenly divides '
94 f'into the image size, so that {ymax} mod binning == 0 '
95 f'and {xmax} mod binning == 0')
96 xarray = np.arange(xmax)
97 yarray = np.arange(ymax)
98 xmesh, ymesh = np.meshgrid(xarray, yarray)
99 xbins = np.arange(0, xmax + binning, binning)
100 ybins = np.arange(0, ymax + binning, binning)
101 numerator = fitWeights * fitData
102 binnedNumerator, *_ = scipy.stats.binned_statistic_2d(ymesh.ravel(), xmesh.ravel(),
103 numerator.ravel(), statistic='sum',
104 bins=(ybins, xbins))
105 binnedDenominator, *_ = scipy.stats.binned_statistic_2d(ymesh.ravel(), xmesh.ravel(),
106 fitWeights.ravel(), statistic='sum',
107 bins=(ybins, xbins))
108 binnedData = np.zeros(binnedNumerator.shape)
109 ind = binnedDenominator != 0
110 np.divide(binnedNumerator, binnedDenominator, out=binnedData, where=ind)
111 binnedWeight = binnedDenominator
112 binMask = (binnedData * binnedWeight**0.5) > detectionThreshold
113 tmpOutputMask = binMask.repeat(binning, axis=0)[:ymax]
114 outputMask = tmpOutputMask.repeat(binning, axis=1)[:, :xmax]
115 else:
116 outputMask = (fitData * fitWeights**0.5) > detectionThreshold
117
118 # Clear existing Detected Plane:
119 maskedImage.mask.array &= ~maskedImage.mask.getPlaneBitMask(detectedPlane)
120
121 # Set Detected Plane with the binary detection mask:
122 maskedImage.mask.array[outputMask] |= maskedImage.mask.getPlaneBitMask(detectedPlane)
123
124
125@dataclass
126class Line:
127 """A simple data class to describe a line profile. The parameter `rho`
128 describes the distance from the center of the image, `theta` describes
129 the angle, and `sigma` describes the width of the line.
130 """
131
132 rho: float
133 theta: float
134 sigma: float = 0
135
136
138 """Collection of `Line` objects.
139
140 Parameters
141 ----------
142 rhos : `np.ndarray`
143 Array of `Line` rho parameters.
144 thetas : `np.ndarray`
145 Array of `Line` theta parameters.
146 sigmas : `np.ndarray`, optional
147 Array of `Line` sigma parameters.
148 """
149
150 def __init__(self, rhos, thetas, sigmas=None):
151 if sigmas is None:
152 sigmas = np.zeros(len(rhos))
153
154 self._lines = [Line(rho, theta, sigma) for (rho, theta, sigma) in
155 zip(rhos, thetas, sigmas)]
156
157 def __len__(self):
158 return len(self._lines)
159
160 def __getitem__(self, index):
161 return self._lines[index]
162
163 def __iter__(self):
164 return iter(self._lines)
165
166 def __repr__(self):
167 joinedString = ", ".join(str(line) for line in self._lines)
168 return textwrap.shorten(joinedString, width=160, placeholder="...")
169
170 @property
171 def rhos(self):
172 return np.array([line.rho for line in self._lines])
173
174 @property
175 def thetas(self):
176 return np.array([line.theta for line in self._lines])
177
178 def append(self, newLine):
179 """Add line to current collection of lines.
180
181 Parameters
182 ----------
183 newLine : `Line`
184 `Line` to add to current collection of lines
185 """
186 self._lines.append(copy.copy(newLine))
187
188
190 """Construct and/or fit a model for a linear streak.
191
192 This assumes a simple model for a streak, in which the streak
193 follows a straight line in pixels space, with a Moffat-shaped profile. The
194 model is fit to data using a Newton-Raphson style minimization algorithm.
195 The initial guess for the line parameters is assumed to be fairly accurate,
196 so only a narrow band of pixels around the initial line estimate is used in
197 fitting the model, which provides a significant speed-up over using all the
198 data. The class can also be used just to construct a model for the data with
199 a line following the given coordinates.
200
201 Parameters
202 ----------
203 data : `np.ndarray`
204 2d array of data.
205 weights : `np.ndarray`
206 2d array of weights.
207 line : `Line`, optional
208 Guess for position of line. Data far from line guess is masked out.
209 Defaults to None, in which case only data with `weights` = 0 is masked
210 out.
211 """
212
213 def __init__(self, data, weights, line=None):
214 self.data = data
215 self.weights = weights
216 self._ymax, self._xmax = data.shape
217 self._dtype = data.dtype
218 xrange = np.arange(self._xmax) - self._xmax / 2.
219 yrange = np.arange(self._ymax) - self._ymax / 2.
220 self._rhoMax = ((0.5 * self._ymax)**2 + (0.5 * self._xmax)**2)**0.5
221 self._xmesh, self._ymesh = np.meshgrid(xrange, yrange)
222 self.mask = (weights != 0)
223
224 self._initLine = line
225 self.setLineMask(line)
226
227 def setLineMask(self, line):
228 """Set mask around the image region near the line.
229
230 Parameters
231 ----------
232 line : `Line`
233 Parameters of line in the image.
234 """
235 if line:
236 # Only fit pixels within 5 sigma of the estimated line
237 radtheta = np.deg2rad(line.theta)
238 distance = (np.cos(radtheta) * self._xmesh + np.sin(radtheta) * self._ymesh - line.rho)
239 m = (abs(distance) < 5 * line.sigma)
240 self.lineMask = self.mask & m
241 else:
242 self.lineMask = np.copy(self.mask)
243
244 self.lineMaskSize = self.lineMask.sum()
245 self._maskData = self.data[self.lineMask]
246 self._maskWeights = self.weights[self.lineMask]
247 self._mxmesh = self._xmesh[self.lineMask]
248 self._mymesh = self._ymesh[self.lineMask]
249
250 def _makeMaskedProfile(self, line, fitFlux=True):
251 """Construct the line model in the masked region and calculate its
252 derivatives.
253
254 Parameters
255 ----------
256 line : `Line`
257 Parameters of line profile for which to make profile in the masked
258 region.
259 fitFlux : `bool`
260 Fit the amplitude of the line profile to the data.
261
262 Returns
263 -------
264 model : `np.ndarray`
265 Model in the masked region.
266 dModel : `np.ndarray`
267 Derivative of the model in the masked region.
268 """
269 invSigma = line.sigma**-1
270 # Calculate distance between pixels and line
271 radtheta = np.deg2rad(line.theta)
272 costheta = np.cos(radtheta)
273 sintheta = np.sin(radtheta)
274 distance = (costheta * self._mxmesh + sintheta * self._mymesh - line.rho)
275 distanceSquared = distance**2
276
277 # Calculate partial derivatives of distance
278 drad = np.pi / 180
279 dDistanceSqdRho = 2 * distance * (-np.ones_like(self._mxmesh))
280 dDistanceSqdTheta = (2 * distance * (-sintheta * self._mxmesh + costheta * self._mymesh) * drad)
281
282 # Use pixel-line distances to make Moffat profile
283 profile = (1 + distanceSquared * invSigma**2)**-2.5
284 dProfile = -2.5 * (1 + distanceSquared * invSigma**2)**-3.5
285
286 if fitFlux:
287 # Calculate line flux from profile and data
288 flux = ((self._maskWeights * self._maskData * profile).sum()
289 / (self._maskWeights * profile**2).sum())
290 else:
291 # Approximately normalize the line
292 flux = invSigma**-1
293 if np.isnan(flux):
294 flux = 0
295
296 model = flux * profile
297
298 # Calculate model derivatives
299 fluxdProfile = flux * dProfile
300 fluxdProfileInvSigma = fluxdProfile * invSigma**2
301 dModeldRho = fluxdProfileInvSigma * dDistanceSqdRho
302 dModeldTheta = fluxdProfileInvSigma * dDistanceSqdTheta
303 dModeldInvSigma = fluxdProfile * distanceSquared * 2 * invSigma
304
305 dModel = np.array([dModeldRho, dModeldTheta, dModeldInvSigma])
306 return model, dModel
307
308 def makeProfile(self, line, fitFlux=True):
309 """Construct the line profile model.
310
311 Parameters
312 ----------
313 line : `Line`
314 Parameters of the line profile to model.
315 fitFlux : `bool`, optional
316 Fit the amplitude of the line profile to the data.
317
318 Returns
319 -------
320 finalModel : `np.ndarray`
321 Model for line profile.
322 """
323 model, _ = self._makeMaskedProfile(line, fitFlux=fitFlux)
324 finalModel = np.zeros((self._ymax, self._xmax), dtype=self._dtype)
325 finalModel[self.lineMask] = model
326 return finalModel
327
328 def _lineChi2(self, line, grad=True):
329 """Construct the chi2 between the data and the model.
330
331 Parameters
332 ----------
333 line : `Line`
334 `Line` parameters for which to build model and calculate chi2.
335 grad : `bool`, optional
336 Whether or not to return the gradient and hessian.
337
338 Returns
339 -------
340 reducedChi : `float`
341 Reduced chi2 of the model.
342 reducedDChi : `np.ndarray`
343 Derivative of the chi2 with respect to rho, theta, invSigma.
344 reducedHessianChi : `np.ndarray`
345 Hessian of the chi2 with respect to rho, theta, invSigma.
346 """
347 # Calculate chi2
348 model, dModel = self._makeMaskedProfile(line)
349 chi2 = (self._maskWeights * (self._maskData - model)**2).sum()
350 if not grad:
351 return chi2.sum() / self.lineMaskSize
352
353 # Calculate derivative and Hessian of chi2
354 derivChi2 = ((-2 * self._maskWeights * (self._maskData - model))[None, :] * dModel).sum(axis=1)
355 hessianChi2 = (2 * self._maskWeights * dModel[:, None, :] * dModel[None, :, :]).sum(axis=2)
356
357 reducedChi = chi2 / self.lineMaskSize
358 reducedDChi = derivChi2 / self.lineMaskSize
359 reducedHessianChi = hessianChi2 / self.lineMaskSize
360 return reducedChi, reducedDChi, reducedHessianChi
361
362 def fit(self, dChi2Tol=0.1, maxIter=100):
363 """Perform Newton-Raphson minimization to find line parameters.
364
365 This method takes advantage of having known derivative and Hessian of
366 the multivariate function to quickly and efficiently find the minimum.
367 This is more efficient than the scipy implementation of the Newton-
368 Raphson method, which doesn't take advantage of the Hessian matrix. The
369 method here also performs a line search in the direction of the steepest
370 derivative at each iteration, which reduces the number of iterations
371 needed.
372
373 Parameters
374 ----------
375 dChi2Tol : `float`, optional
376 Change in Chi2 tolerated for fit convergence.
377 maxIter : `int`, optional
378 Maximum number of fit iterations allowed. The fit should converge in
379 ~10 iterations, depending on the value of dChi2Tol, but this
380 maximum provides a backup.
381
382 Returns
383 -------
384 outline : `np.ndarray`
385 Coordinates and inverse width of fit line.
386 chi2 : `float`
387 Reduced Chi2 of model fit to data.
388 fitFailure : `bool`
389 Boolean where `False` corresponds to a successful fit.
390 """
391 # Do minimization on inverse of sigma to simplify derivatives:
392 x = np.array([self._initLine.rho, self._initLine.theta, self._initLine.sigma**-1])
393
394 dChi2 = 1
395 iter = 0
396 oldChi2 = 0
397 fitFailure = False
398
399 def line_search(c, dx):
400 testx = x - c * dx
401 testLine = Line(testx[0], testx[1], testx[2]**-1)
402 return self._lineChi2(testLine, grad=False)
403
404 while abs(dChi2) > dChi2Tol:
405 line = Line(x[0], x[1], x[2]**-1)
406 chi2, b, A = self._lineChi2(line)
407 if chi2 == 0:
408 break
409 if not np.isfinite(A).all():
410 # TODO: DM-30797 Add warning here.
411 fitFailure = True
412 break
413 dChi2 = oldChi2 - chi2
414 cholesky = scipy.linalg.cho_factor(A)
415 dx = scipy.linalg.cho_solve(cholesky, b)
416
417 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05)
418 x -= factor * dx
419 if (abs(x[0]) > 1.5 * self._rhoMax) or (iter > maxIter):
420 fitFailure = True
421 break
422 oldChi2 = chi2
423 iter += 1
424
425 outline = Line(x[0], x[1], abs(x[2])**-1)
426
427 return outline, chi2, fitFailure
428
429
430class MaskStreaksConfig(pexConfig.Config):
431 """Configuration parameters for `MaskStreaksTask`.
432 """
433
434 minimumKernelHeight = pexConfig.Field(
435 doc="Minimum height of the streak-finding kernel relative to the tallest kernel",
436 dtype=float,
437 default=0.0,
438 )
439 absMinimumKernelHeight = pexConfig.Field(
440 doc="Minimum absolute height of the streak-finding kernel",
441 dtype=float,
442 default=5,
443 )
444 clusterMinimumSize = pexConfig.Field(
445 doc="Minimum size in pixels of detected clusters",
446 dtype=int,
447 default=50,
448 )
449 clusterMinimumDeviation = pexConfig.Field(
450 doc="Allowed deviation (in pixels) from a straight line for a detected "
451 "line",
452 dtype=int,
453 default=2,
454 )
455 delta = pexConfig.Field(
456 doc="Stepsize in angle-radius parameter space",
457 dtype=float,
458 default=0.2,
459 )
460 nSigma = pexConfig.Field(
461 doc="Number of sigmas from center of kernel to include in voting "
462 "procedure",
463 dtype=float,
464 default=2,
465 )
466 rhoBinSize = pexConfig.Field(
467 doc="Binsize in pixels for position parameter rho when finding "
468 "clusters of detected lines",
469 dtype=float,
470 default=30,
471 )
472 thetaBinSize = pexConfig.Field(
473 doc="Binsize in degrees for angle parameter theta when finding "
474 "clusters of detected lines",
475 dtype=float,
476 default=2,
477 )
478 invSigma = pexConfig.Field(
479 doc="Inverse of the Moffat sigma parameter (in units of pixels)"
480 "describing the profile of the streak",
481 dtype=float,
482 default=10.**-1,
483 )
484 footprintThreshold = pexConfig.Field(
485 doc="Threshold at which to determine edge of line, in units of "
486 "nanoJanskys",
487 dtype=float,
488 default=0.01
489 )
490 dChi2Tolerance = pexConfig.Field(
491 doc="Absolute difference in Chi2 between iterations of line profile"
492 "fitting that is acceptable for convergence",
493 dtype=float,
494 default=0.1
495 )
496 detectedMaskPlane = pexConfig.Field(
497 doc="Name of mask with pixels above detection threshold, used for first"
498 "estimate of streak locations",
499 dtype=str,
500 default="DETECTED"
501 )
502 streaksMaskPlane = pexConfig.Field(
503 doc="Name of mask plane holding detected streaks",
504 dtype=str,
505 default="STREAK"
506 )
507
508
509class MaskStreaksTask(pipeBase.Task):
510 """Find streaks or other straight lines in image data.
511
512 Nearby objects passing through the field of view of the telescope leave a
513 bright trail in images. This class uses the Kernel Hough Transform (KHT)
514 (Fernandes and Oliveira, 2007), implemented in `lsst.houghtransform`. The
515 procedure works by taking a binary image, either provided as put or produced
516 from the input data image, using a Canny filter to make an image of the
517 edges in the original image, then running the KHT on the edge image. The KHT
518 identifies clusters of non-zero points, breaks those clusters of points into
519 straight lines, keeps clusters with a size greater than the user-set
520 threshold, then performs a voting procedure to find the best-fit coordinates
521 of any straight lines. Given the results of the KHT algorithm, clusters of
522 lines are identified and grouped (generally these correspond to the two
523 edges of a strea) and a profile is fit to the streak in the original
524 (non-binary) image.
525 """
526
527 ConfigClass = MaskStreaksConfig
528 _DefaultName = "maskStreaks"
529
530 @timeMethod
531 def find(self, maskedImage):
532 """Find streaks in a masked image.
533
534 Parameters
535 ----------
536 maskedImage : `lsst.afw.image.maskedImage`
537 The image in which to search for streaks.
538
539 Returns
540 -------
541 result : `lsst.pipe.base.Struct`
542 Results as a struct with attributes:
543
544 ``originalLines``
545 Lines identified by kernel hough transform.
546 ``lineClusters``
547 Lines grouped into clusters in rho-theta space.
548 ``lines``
549 Final result for lines after line-profile fit.
550 ``mask``
551 2-d boolean mask where detected lines are True.
552 """
553 mask = maskedImage.getMask()
554 detectionMask = (mask.array & mask.getPlaneBitMask(self.config.detectedMaskPlane))
555
556 self.edges = self._cannyFilter(detectionMask)
557 self.lines = self._runKHT(self.edges)
558
559 if len(self.lines) == 0:
560 lineMask = np.zeros(detectionMask.shape, dtype=bool)
561 fitLines = LineCollection([], [])
562 clusters = LineCollection([], [])
563 else:
564 clusters = self._findClusters(self.lines)
565 fitLines, lineMask = self._fitProfile(clusters, maskedImage)
566
567 # The output mask is the intersection of the fit streaks and the image detections
568 outputMask = lineMask & detectionMask.astype(bool)
569
570 return pipeBase.Struct(
571 lines=fitLines,
572 lineClusters=clusters,
573 originalLines=self.lines,
574 mask=outputMask,
575 )
576
577 @timeMethod
578 def run(self, maskedImage):
579 """Find and mask streaks in a masked image.
580
581 Finds streaks in the image and modifies maskedImage in place by adding a
582 mask plane with any identified streaks.
583
584 Parameters
585 ----------
586 maskedImage : `lsst.afw.image.maskedImage`
587 The image in which to search for streaks. The mask detection plane
588 corresponding to `config.detectedMaskPlane` must be set with the
589 detected pixels.
590
591 Returns
592 -------
593 result : `lsst.pipe.base.Struct`
594 Results as a struct with attributes:
595
596 ``originalLines``
597 Lines identified by kernel hough transform.
598 ``lineClusters``
599 Lines grouped into clusters in rho-theta space.
600 ``lines``
601 Final result for lines after line-profile fit.
602 """
603 streaks = self.find(maskedImage)
604
605 maskedImage.mask.addMaskPlane(self.config.streaksMaskPlane)
606 maskedImage.mask.array[streaks.mask] |= maskedImage.mask.getPlaneBitMask(self.config.streaksMaskPlane)
607
608 return pipeBase.Struct(
609 lines=streaks.lines,
610 lineClusters=streaks.lineClusters,
611 originalLines=streaks.originalLines,
612 )
613
614 def _cannyFilter(self, image):
615 """Apply a canny filter to the data in order to detect edges.
616
617 Parameters
618 ----------
619 image : `np.ndarray`
620 2-d image data on which to run filter.
621
622 Returns
623 -------
624 cannyData : `np.ndarray`
625 2-d image of edges found in input image.
626 """
627 filterData = image.astype(int)
628 return canny(filterData, low_threshold=0, high_threshold=1, sigma=0.1)
629
630 def _runKHT(self, image):
631 """Run Kernel Hough Transform on image.
632
633 Parameters
634 ----------
635 image : `np.ndarray`
636 2-d image data on which to detect lines.
637
638 Returns
639 -------
640 result : `LineCollection`
641 Collection of detected lines, with their detected rho and theta
642 coordinates.
643 """
644 lines = lsst.kht.find_lines(image, self.config.clusterMinimumSize,
645 self.config.clusterMinimumDeviation, self.config.delta,
646 self.config.minimumKernelHeight, self.config.nSigma,
647 self.config.absMinimumKernelHeight)
648
649 return LineCollection(lines.rho, lines.theta)
650
651 def _findClusters(self, lines):
652 """Group lines that are close in parameter space and likely describe
653 the same streak.
654
655 Parameters
656 ----------
657 lines : `LineCollection`
658 Collection of lines to group into clusters.
659
660 Returns
661 -------
662 result : `LineCollection`
663 Average `Line` for each cluster of `Line`s in the input
664 `LineCollection`.
665 """
666 # Scale variables by threshold bin-size variable so that rho and theta
667 # are on the same scale. Since the clustering algorithm below stops when
668 # the standard deviation <= 1, after rescaling each cluster will have a
669 # standard deviation at or below the bin-size.
670 x = lines.rhos / self.config.rhoBinSize
671 y = lines.thetas / self.config.thetaBinSize
672 X = np.array([x, y]).T
673 nClusters = 1
674
675 # Put line parameters in clusters by starting with all in one, then
676 # subdividing until the parameters of each cluster have std dev=1.
677 # If nClusters == len(lines), each line will have its own 'cluster', so
678 # the standard deviations of each cluster must be zero and the loop
679 # is guaranteed to stop.
680 while True:
681 kmeans = KMeans(n_clusters=nClusters).fit(X)
682 clusterStandardDeviations = np.zeros((nClusters, 2))
683 for c in range(nClusters):
684 inCluster = X[kmeans.labels_ == c]
685 clusterStandardDeviations[c] = np.std(inCluster, axis=0)
686 # Are the rhos and thetas in each cluster all below the threshold?
687 if (clusterStandardDeviations <= 1).all():
688 break
689 nClusters += 1
690
691 # The cluster centers are final line estimates
692 finalClusters = kmeans.cluster_centers_.T
693
694 # Rescale variables:
695 finalRhos = finalClusters[0] * self.config.rhoBinSize
696 finalThetas = finalClusters[1] * self.config.thetaBinSize
697 result = LineCollection(finalRhos, finalThetas)
698
699 return result
700
701 def _fitProfile(self, lines, maskedImage):
702 """Fit the profile of the streak.
703
704 Given the initial parameters of detected lines, fit a model for the
705 streak to the original (non-binary image). The assumed model is a
706 straight line with a Moffat profile.
707
708 Parameters
709 ----------
710 lines : `LineCollection`
711 Collection of guesses for `Line`s detected in the image.
712 maskedImage : `lsst.afw.image.maskedImage`
713 Original image to be used to fit profile of streak.
714
715 Returns
716 -------
717 lineFits : `LineCollection`
718 Collection of `Line` profiles fit to the data.
719 finalMask : `np.ndarray`
720 2d mask array with detected streaks=1.
721 """
722 data = maskedImage.image.array
723 weights = maskedImage.variance.array**-1
724 # Mask out any pixels with non-finite weights
725 weights[~np.isfinite(weights) | ~np.isfinite(data)] = 0
726
727 lineFits = LineCollection([], [])
728 finalLineMasks = [np.zeros(data.shape, dtype=bool)]
729 for line in lines:
730 line.sigma = self.config.invSigma**-1
731 lineModel = LineProfile(data, weights, line=line)
732 # Skip any lines that do not cover any data (sometimes happens because of chip gaps)
733 if lineModel.lineMaskSize == 0:
734 continue
735
736 fit, chi2, fitFailure = lineModel.fit(dChi2Tol=self.config.dChi2Tolerance)
737
738 # Initial estimate should be quite close: fit is deemed unsuccessful if rho or theta
739 # change more than the allowed bin in rho or theta:
740 if ((abs(fit.rho - line.rho) > 2 * self.config.rhoBinSize)
741 or (abs(fit.theta - line.theta) > 2 * self.config.thetaBinSize)):
742 fitFailure = True
743
744 if fitFailure:
745 continue
746
747 # Make mask
748 lineModel.setLineMask(fit)
749 finalModel = lineModel.makeProfile(fit)
750 # Take absolute value, as streaks are allowed to be negative
751 finalModelMax = abs(finalModel).max()
752 finalLineMask = abs(finalModel) > self.config.footprintThreshold
753 # Drop this line if the model profile is below the footprint threshold
754 if not finalLineMask.any():
755 continue
756 fit.chi2 = chi2
757 fit.finalModelMax = finalModelMax
758 lineFits.append(fit)
759 finalLineMasks.append(finalLineMask)
760
761 finalMask = np.array(finalLineMasks).any(axis=0)
762
763 return lineFits, finalMask
int max
table::Key< int > to
table::Key< int > a
def __init__(self, rhos, thetas, sigmas=None)
Definition: maskStreaks.py:150
def __init__(self, data, weights, line=None)
Definition: maskStreaks.py:213
def makeProfile(self, line, fitFlux=True)
Definition: maskStreaks.py:308
def _makeMaskedProfile(self, line, fitFlux=True)
Definition: maskStreaks.py:250
def _lineChi2(self, line, grad=True)
Definition: maskStreaks.py:328
def fit(self, dChi2Tol=0.1, maxIter=100)
Definition: maskStreaks.py:362
def _fitProfile(self, lines, maskedImage)
Definition: maskStreaks.py:701
def setDetectionMask(maskedImage, forceSlowBin=False, binning=None, detectedPlane="DETECTED", badMaskPlanes=("NO_DATA", "INTRP", "BAD", "SAT", "EDGE"), detectionThreshold=5)
Definition: maskStreaks.py:40