LSST Applications g0f08755f38+82efc23009,g12f32b3c4e+e7bdf1200e,g1653933729+a8ce1bb630,g1a0ca8cf93+50eff2b06f,g28da252d5a+52db39f6a5,g2bbee38e9b+37c5a29d61,g2bc492864f+37c5a29d61,g2cdde0e794+c05ff076ad,g3156d2b45e+41e33cbcdc,g347aa1857d+37c5a29d61,g35bb328faa+a8ce1bb630,g3a166c0a6a+37c5a29d61,g3e281a1b8c+fb992f5633,g414038480c+7f03dfc1b0,g41af890bb2+11b950c980,g5fbc88fb19+17cd334064,g6b1c1869cb+12dd639c9a,g781aacb6e4+a8ce1bb630,g80478fca09+72e9651da0,g82479be7b0+04c31367b4,g858d7b2824+82efc23009,g9125e01d80+a8ce1bb630,g9726552aa6+8047e3811d,ga5288a1d22+e532dc0a0b,gae0086650b+a8ce1bb630,gb58c049af0+d64f4d3760,gc28159a63d+37c5a29d61,gcf0d15dbbd+2acd6d4d48,gd7358e8bfb+778a810b6e,gda3e153d99+82efc23009,gda6a2b7d83+2acd6d4d48,gdaeeff99f8+1711a396fd,ge2409df99d+6b12de1076,ge79ae78c31+37c5a29d61,gf0baf85859+d0a5978c5a,gf3967379c6+4954f8c433,gfb92a5be7c+82efc23009,gfec2e1e490+2aaed99252,w.2024.46
LSST Data Management Base Package
Loading...
Searching...
No Matches
maskStreaks.py
Go to the documentation of this file.
1# This file is part of meas_algorithms.
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"]
23
24import numpy as np
25import scipy
26import textwrap
27import copy
28from skimage.feature import canny
29from sklearn.cluster import KMeans
30from dataclasses import dataclass
31import astropy.units as u
32
33from lsst.afw.geom import SpanSet
34import lsst.pex.config as pexConfig
35import lsst.pipe.base as pipeBase
36import lsst.kht
37from lsst.utils.timer import timeMethod
38
39
40@dataclass
41class Line:
42 """A simple data class to describe a line profile. The parameter `rho`
43 describes the distance from the center of the image, `theta` describes
44 the angle, `sigma` describes the width of the line, `reducedChi2` gives
45 the reduced chi2 of the fit, and `modelMaximum` gives the peak value of the
46 fit line profile.
47 """
48
49 rho: float
50 theta: float
51 sigma: float = 0
52 reducedChi2: float = 0
53 modelMaximum: float = 0
54
55
57 """Collection of `Line` objects.
58
59 Parameters
60 ----------
61 rhos : `np.ndarray`
62 Array of `Line` rho parameters.
63 thetas : `np.ndarray`
64 Array of `Line` theta parameters.
65 sigmas : `np.ndarray`, optional
66 Array of `Line` sigma parameters.
67 """
68
69 def __init__(self, rhos, thetas, sigmas=None):
70 if sigmas is None:
71 sigmas = np.zeros(len(rhos))
72
73 self._lines = [Line(rho, theta, sigma) for (rho, theta, sigma) in
74 zip(rhos, thetas, sigmas)]
75
76 def __len__(self):
77 return len(self._lines)
78
79 def __getitem__(self, index):
80 return self._lines[index]
81
82 def __iter__(self):
83 return iter(self._lines)
84
85 def __repr__(self):
86 joinedString = ", ".join(str(line) for line in self._lines)
87 return textwrap.shorten(joinedString, width=160, placeholder="...")
88
89 @property
90 def rhos(self):
91 return np.array([line.rho for line in self._lines])
92
93 @property
94 def thetas(self):
95 return np.array([line.theta for line in self._lines])
96
97 @property
98 def sigmas(self):
99 return np.array([line.sigma for line in self._lines])
100
101 def append(self, newLine):
102 """Add line to current collection of lines.
103
104 Parameters
105 ----------
106 newLine : `Line`
107 `Line` to add to current collection of lines
108 """
109 self._lines.append(copy.copy(newLine))
110
111
113 """Construct and/or fit a model for a linear streak.
114
115 This assumes a simple model for a streak, in which the streak
116 follows a straight line in pixels space, with a Moffat-shaped profile. The
117 model is fit to data using a Newton-Raphson style minimization algorithm.
118 The initial guess for the line parameters is assumed to be fairly accurate,
119 so only a narrow band of pixels around the initial line estimate is used in
120 fitting the model, which provides a significant speed-up over using all the
121 data. The class can also be used just to construct a model for the data with
122 a line following the given coordinates.
123
124 Parameters
125 ----------
126 data : `np.ndarray`
127 2d array of data.
128 weights : `np.ndarray`
129 2d array of weights.
130 line : `Line`, optional
131 Guess for position of line. Data far from line guess is masked out.
132 Defaults to None, in which case only data with `weights` = 0 is masked
133 out.
134 """
135
136 def __init__(self, data, weights, line=None):
137 self.data = data
138 self.weights = weights
139 self._ymax, self._xmax = data.shape
140 self._dtype = data.dtype
141 xrange = np.arange(self._xmax) - self._xmax / 2.
142 yrange = np.arange(self._ymax) - self._ymax / 2.
143 self._rhoMax = ((0.5 * self._ymax)**2 + (0.5 * self._xmax)**2)**0.5
144 self._xmesh, self._ymesh = np.meshgrid(xrange, yrange)
145 self.mask = (weights != 0)
146
147 self._initLine = line
148 self.setLineMask(line)
149
150 def getLineXY(self, line):
151 """Return the pixel coordinates of the ends of the line.
152
153 Parameters
154 ----------
155 line : `Line`
156 Line for which to find the endpoints.
157
158 Returns
159 -------
160 boxIntersections : `np.ndarray`
161 (x, y) coordinates of the start and endpoints of the line.
162 """
163 theta = line.theta * u.deg
164 # Determine where the line intersects with each edge of the bounding
165 # box.
166 # Bottom:
167 yA = -self._ymax / 2.
168 xA = (line.rho - yA * np.sin(theta)) / np.cos(theta)
169 # Left:
170 xB = -self._xmax / 2.
171 yB = (line.rho - xB * np.cos(theta)) / np.sin(theta)
172 # Top:
173 yC = self._ymax / 2.
174 xC = (line.rho - yC * np.sin(theta)) / np.cos(theta)
175 # Right:
176 xD = self._xmax / 2.
177 yD = (line.rho - xD * np.cos(theta)) / np.sin(theta)
178
179 lineIntersections = np.array([[xA, yA],
180 [xB, yB],
181 [xC, yC],
182 [xD, yD]])
183 lineIntersections[:, 0] += self._xmax / 2.
184 lineIntersections[:, 1] += self._ymax / 2.
185
186 # The line will necessarily intersect with exactly two edges of the
187 # bounding box itself.
188 inBox = ((lineIntersections[:, 0] >= 0) & (lineIntersections[:, 0] <= self._xmax)
189 & (lineIntersections[:, 1] >= 0) & (lineIntersections[:, 1] <= self._ymax))
190 boxIntersections = lineIntersections[inBox]
191
192 return boxIntersections
193
194 def setLineMask(self, line):
195 """Set mask around the image region near the line.
196
197 Parameters
198 ----------
199 line : `Line`
200 Parameters of line in the image.
201 """
202 if line:
203 # Only fit pixels within 5 sigma of the estimated line
204 radtheta = np.deg2rad(line.theta)
205 distance = (np.cos(radtheta) * self._xmesh + np.sin(radtheta) * self._ymesh - line.rho)
206 m = (abs(distance) < 5 * line.sigma)
207 self.lineMask = self.mask & m
208 else:
209 self.lineMask = np.copy(self.mask)
210
211 self.lineMaskSize = self.lineMask.sum()
212 self._maskData = self.data[self.lineMask]
213 self._maskWeights = self.weights[self.lineMask]
214 self._mxmesh = self._xmesh[self.lineMask]
215 self._mymesh = self._ymesh[self.lineMask]
216
217 def _makeMaskedProfile(self, line, fitFlux=True):
218 """Construct the line model in the masked region and calculate its
219 derivatives.
220
221 Parameters
222 ----------
223 line : `Line`
224 Parameters of line profile for which to make profile in the masked
225 region.
226 fitFlux : `bool`
227 Fit the amplitude of the line profile to the data.
228
229 Returns
230 -------
231 model : `np.ndarray`
232 Model in the masked region.
233 dModel : `np.ndarray`
234 Derivative of the model in the masked region.
235 """
236 invSigma = line.sigma**-1
237 # Calculate distance between pixels and line
238 radtheta = np.deg2rad(line.theta)
239 costheta = np.cos(radtheta)
240 sintheta = np.sin(radtheta)
241 distance = (costheta * self._mxmesh + sintheta * self._mymesh - line.rho)
242 distanceSquared = distance**2
243
244 # Calculate partial derivatives of distance
245 drad = np.pi / 180
246 dDistanceSqdRho = 2 * distance * (-np.ones_like(self._mxmesh))
247 dDistanceSqdTheta = (2 * distance * (-sintheta * self._mxmesh + costheta * self._mymesh) * drad)
248
249 # Use pixel-line distances to make Moffat profile
250 profile = (1 + distanceSquared * invSigma**2)**-2.5
251 dProfile = -2.5 * (1 + distanceSquared * invSigma**2)**-3.5
252
253 if fitFlux:
254 # Calculate line flux from profile and data
255 flux = ((self._maskWeights * self._maskData * profile).sum()
256 / (self._maskWeights * profile**2).sum())
257 else:
258 # Approximately normalize the line
259 flux = invSigma**-1
260 if np.isnan(flux):
261 flux = 0
262
263 model = flux * profile
264
265 # Calculate model derivatives
266 fluxdProfile = flux * dProfile
267 fluxdProfileInvSigma = fluxdProfile * invSigma**2
268 dModeldRho = fluxdProfileInvSigma * dDistanceSqdRho
269 dModeldTheta = fluxdProfileInvSigma * dDistanceSqdTheta
270 dModeldInvSigma = fluxdProfile * distanceSquared * 2 * invSigma
271
272 dModel = np.array([dModeldRho, dModeldTheta, dModeldInvSigma])
273 return model, dModel
274
275 def makeProfile(self, line, fitFlux=True):
276 """Construct the line profile model.
277
278 Parameters
279 ----------
280 line : `Line`
281 Parameters of the line profile to model.
282 fitFlux : `bool`, optional
283 Fit the amplitude of the line profile to the data.
284
285 Returns
286 -------
287 finalModel : `np.ndarray`
288 Model for line profile.
289 """
290 model, _ = self._makeMaskedProfile(line, fitFlux=fitFlux)
291 finalModel = np.zeros((self._ymax, self._xmax), dtype=self._dtype)
292 finalModel[self.lineMask] = model
293 return finalModel
294
295 def _lineChi2(self, line, grad=True):
296 """Construct the chi2 between the data and the model.
297
298 Parameters
299 ----------
300 line : `Line`
301 `Line` parameters for which to build model and calculate chi2.
302 grad : `bool`, optional
303 Whether or not to return the gradient and hessian.
304
305 Returns
306 -------
307 reducedChi : `float`
308 Reduced chi2 of the model.
309 reducedDChi : `np.ndarray`
310 Derivative of the chi2 with respect to rho, theta, invSigma.
311 reducedHessianChi : `np.ndarray`
312 Hessian of the chi2 with respect to rho, theta, invSigma.
313 """
314 # Calculate chi2
315 model, dModel = self._makeMaskedProfile(line)
316 chi2 = (self._maskWeights * (self._maskData - model)**2).sum()
317 if not grad:
318 return chi2.sum() / self.lineMaskSize
319
320 # Calculate derivative and Hessian of chi2
321 derivChi2 = ((-2 * self._maskWeights * (self._maskData - model))[None, :] * dModel).sum(axis=1)
322 hessianChi2 = (2 * self._maskWeights * dModel[:, None, :] * dModel[None, :, :]).sum(axis=2)
323
324 reducedChi = chi2 / self.lineMaskSize
325 reducedDChi = derivChi2 / self.lineMaskSize
326 reducedHessianChi = hessianChi2 / self.lineMaskSize
327 return reducedChi, reducedDChi, reducedHessianChi
328
329 def fit(self, dChi2Tol=0.1, maxIter=100, log=None):
330 """Perform Newton-Raphson minimization to find line parameters.
331
332 This method takes advantage of having known derivative and Hessian of
333 the multivariate function to quickly and efficiently find the minimum.
334 This is more efficient than the scipy implementation of the Newton-
335 Raphson method, which doesn't take advantage of the Hessian matrix. The
336 method here also performs a line search in the direction of the steepest
337 derivative at each iteration, which reduces the number of iterations
338 needed.
339
340 Parameters
341 ----------
342 dChi2Tol : `float`, optional
343 Change in Chi2 tolerated for fit convergence.
344 maxIter : `int`, optional
345 Maximum number of fit iterations allowed. The fit should converge in
346 ~10 iterations, depending on the value of dChi2Tol, but this
347 maximum provides a backup.
348 log : `lsst.utils.logging.LsstLogAdapter`, optional
349 Logger to use for reporting more details for fitting failures.
350
351 Returns
352 -------
353 outline : `Line`
354 Coordinates, inverse width, and chi2 of fit line.
355 fitFailure : `bool`
356 Boolean where `False` corresponds to a successful fit.
357 """
358 # Do minimization on inverse of sigma to simplify derivatives:
359 x = np.array([self._initLine.rho, self._initLine.theta, self._initLine.sigma**-1])
360
361 dChi2 = 1
362 iter = 0
363 oldChi2 = 0
364 fitFailure = False
365
366 def line_search(c, dx):
367 testx = x - c * dx
368 testLine = Line(testx[0], testx[1], testx[2]**-1)
369 return self._lineChi2(testLine, grad=False)
370
371 while abs(dChi2) > dChi2Tol:
372 line = Line(x[0], x[1], x[2]**-1)
373 chi2, b, A = self._lineChi2(line)
374 if chi2 == 0:
375 break
376 if not np.isfinite(A).all():
377 fitFailure = True
378 if log is not None:
379 log.warning("Hessian matrix has non-finite elements.")
380 break
381 dChi2 = oldChi2 - chi2
382 try:
383 cholesky = scipy.linalg.cho_factor(A)
384 except np.linalg.LinAlgError:
385 fitFailure = True
386 if log is not None:
387 log.warning("Hessian matrix is not invertible.")
388 break
389 dx = scipy.linalg.cho_solve(cholesky, b)
390
391 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05)
392 x -= factor * dx
393 if (abs(x[0]) > 1.5 * self._rhoMax) or (iter > maxIter):
394 fitFailure = True
395 break
396 oldChi2 = chi2
397 iter += 1
398
399 outline = Line(x[0], x[1], abs(x[2])**-1, chi2)
400
401 return outline, fitFailure
402
403
404class MaskStreaksConfig(pexConfig.Config):
405 """Configuration parameters for `MaskStreaksTask`.
406 """
407
408 minimumKernelHeight = pexConfig.Field(
409 doc="Minimum height of the streak-finding kernel relative to the tallest kernel",
410 dtype=float,
411 default=0.0,
412 )
413 absMinimumKernelHeight = pexConfig.Field(
414 doc="Minimum absolute height of the streak-finding kernel",
415 dtype=float,
416 default=5,
417 )
418 clusterMinimumSize = pexConfig.Field(
419 doc="Minimum size in pixels of detected clusters",
420 dtype=int,
421 default=50,
422 )
423 clusterMinimumDeviation = pexConfig.Field(
424 doc="Allowed deviation (in pixels) from a straight line for a detected "
425 "line",
426 dtype=int,
427 default=2,
428 )
429 delta = pexConfig.Field(
430 doc="Stepsize in angle-radius parameter space",
431 dtype=float,
432 default=0.2,
433 )
434 nSigma = pexConfig.Field(
435 doc="Number of sigmas from center of kernel to include in voting "
436 "procedure",
437 dtype=float,
438 default=2,
439 )
440 rhoBinSize = pexConfig.Field(
441 doc="Binsize in pixels for position parameter rho when finding "
442 "clusters of detected lines",
443 dtype=float,
444 default=30,
445 )
446 thetaBinSize = pexConfig.Field(
447 doc="Binsize in degrees for angle parameter theta when finding "
448 "clusters of detected lines",
449 dtype=float,
450 default=2,
451 )
452 invSigma = pexConfig.Field(
453 doc="Inverse of the Moffat sigma parameter (in units of pixels)"
454 "describing the profile of the streak",
455 dtype=float,
456 default=10.**-1,
457 )
458 footprintThreshold = pexConfig.Field(
459 doc="Threshold at which to determine edge of line, in units of "
460 "nanoJanskys",
461 dtype=float,
462 default=0.01
463 )
464 dChi2Tolerance = pexConfig.Field(
465 doc="Absolute difference in Chi2 between iterations of line profile"
466 "fitting that is acceptable for convergence",
467 dtype=float,
468 default=0.1
469 )
470 detectedMaskPlane = pexConfig.Field(
471 doc="Name of mask with pixels above detection threshold, used for first"
472 "estimate of streak locations",
473 dtype=str,
474 default="DETECTED"
475 )
476 onlyMaskDetected = pexConfig.Field(
477 doc=("If true, only propagate the part of the streak mask that "
478 "overlaps with the detection mask."),
479 dtype=bool,
480 default=True,
481 )
482 streaksMaskPlane = pexConfig.Field(
483 doc="Name of mask plane holding detected streaks",
484 dtype=str,
485 default="STREAK"
486 )
487 badMaskPlanes = pexConfig.ListField(
488 doc=("Names of mask plane regions to ignore entirely when doing streak"
489 " detection"),
490 dtype=str,
491 default=("NO_DATA", "INTRP", "BAD", "SAT", "EDGE"),
492 )
493
494
495class MaskStreaksTask(pipeBase.Task):
496 """Find streaks or other straight lines in image data.
497
498 Nearby objects passing through the field of view of the telescope leave a
499 bright trail in images. This class uses the Kernel Hough Transform (KHT)
500 (Fernandes and Oliveira, 2007), implemented in `lsst.houghtransform`. The
501 procedure works by taking a binary image, either provided as put or produced
502 from the input data image, using a Canny filter to make an image of the
503 edges in the original image, then running the KHT on the edge image. The KHT
504 identifies clusters of non-zero points, breaks those clusters of points into
505 straight lines, keeps clusters with a size greater than the user-set
506 threshold, then performs a voting procedure to find the best-fit coordinates
507 of any straight lines. Given the results of the KHT algorithm, clusters of
508 lines are identified and grouped (generally these correspond to the two
509 edges of a strea) and a profile is fit to the streak in the original
510 (non-binary) image.
511 """
512
513 ConfigClass = MaskStreaksConfig
514 _DefaultName = "maskStreaks"
515
516 @timeMethod
517 def find(self, maskedImage):
518 """Find streaks in a masked image.
519
520 Parameters
521 ----------
522 maskedImage : `lsst.afw.image.maskedImage`
523 The image in which to search for streaks.
524
525 Returns
526 -------
527 result : `lsst.pipe.base.Struct`
528 Results as a struct with attributes:
529
530 ``originalLines``
531 Lines identified by kernel hough transform.
532 ``lineClusters``
533 Lines grouped into clusters in rho-theta space.
534 ``lines``
535 Final result for lines after line-profile fit.
536 ``mask``
537 2-d boolean mask where detected lines are True.
538 """
539 mask = maskedImage.mask
540 detectionMask = (mask.array & mask.getPlaneBitMask(self.config.detectedMaskPlane))
541
542 initEdges = self._cannyFilter(detectionMask)
543 # Ignore regions with known bad masks, adding a one-pixel buffer around
544 # each to ensure that the edges of bad regions are also ignored.
545 badPixelMask = mask.getPlaneBitMask(self.config.badMaskPlanes)
546 badMaskSpanSet = SpanSet.fromMask(mask, badPixelMask).split()
547 for sset in badMaskSpanSet:
548 sset_dilated = sset.dilated(1)
549 sset_dilated.clippedTo(mask.getBBox()).setMask(mask, mask.getPlaneBitMask("BAD"))
550 dilatedBadMask = (mask.array & badPixelMask) > 0
551 self.edges = initEdges & ~dilatedBadMask
552 self.lines = self._runKHT(self.edges)
553
554 if len(self.lines) == 0:
555 lineMask = np.zeros(detectionMask.shape, dtype=bool)
556 fitLines = LineCollection([], [])
557 clusters = LineCollection([], [])
558 else:
559 clusters = self._findClusters(self.lines)
560 fitLines, lineMask = self._fitProfile(clusters, maskedImage)
561
562 if self.config.onlyMaskDetected:
563 # The output mask is the intersection of the fit streaks and the image detections
564 lineMask &= detectionMask.astype(bool)
565
566 return pipeBase.Struct(
567 lines=fitLines,
568 lineClusters=clusters,
569 originalLines=self.lines,
570 mask=lineMask,
571 )
572
573 @timeMethod
574 def run(self, maskedImage):
575 """Find and mask streaks in a masked image.
576
577 Finds streaks in the image and modifies maskedImage in place by adding a
578 mask plane with any identified streaks.
579
580 Parameters
581 ----------
582 maskedImage : `lsst.afw.image.Exposure` or `lsst.afw.image.maskedImage`
583 The image in which to search for streaks. The mask detection plane
584 corresponding to `config.detectedMaskPlane` must be set with the
585 detected pixels. The mask will have a plane added with any detected
586 streaks, and with the mask plane name set by
587 self.config.streaksMaskPlane.
588
589 Returns
590 -------
591 result : `lsst.pipe.base.Struct`
592 Results as a struct with attributes:
593
594 ``originalLines``
595 Lines identified by kernel hough transform.
596 ``lineClusters``
597 Lines grouped into clusters in rho-theta space.
598 ``lines``
599 Final result for lines after line-profile fit.
600 """
601 streaks = self.find(maskedImage)
602
603 maskedImage.mask.addMaskPlane(self.config.streaksMaskPlane)
604 maskedImage.mask.array[streaks.mask] |= maskedImage.mask.getPlaneBitMask(self.config.streaksMaskPlane)
605
606 return pipeBase.Struct(
607 lines=streaks.lines,
608 lineClusters=streaks.lineClusters,
609 originalLines=streaks.originalLines,
610 )
611
612 def _cannyFilter(self, image):
613 """Apply a canny filter to the data in order to detect edges.
614
615 Parameters
616 ----------
617 image : `np.ndarray`
618 2-d image data on which to run filter.
619
620 Returns
621 -------
622 cannyData : `np.ndarray`
623 2-d image of edges found in input image.
624 """
625 # Ensure that the pixels are zero or one. Change the datatype to
626 # np.float64 to be compatible with the Canny filter routine.
627 filterData = (image > 0).astype(np.float64)
628 return canny(filterData, use_quantiles=True, 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 self.log.info("The Kernel Hough Transform detected %s line(s)", len(lines))
649
650 return LineCollection(lines.rho, lines.theta)
651
652 def _findClusters(self, lines):
653 """Group lines that are close in parameter space and likely describe
654 the same streak.
655
656 Parameters
657 ----------
658 lines : `LineCollection`
659 Collection of lines to group into clusters.
660
661 Returns
662 -------
663 result : `LineCollection`
664 Average `Line` for each cluster of `Line`s in the input
665 `LineCollection`.
666 """
667 # Scale variables by threshold bin-size variable so that rho and theta
668 # are on the same scale. Since the clustering algorithm below stops when
669 # the standard deviation <= 1, after rescaling each cluster will have a
670 # standard deviation at or below the bin-size.
671 x = lines.rhos / self.config.rhoBinSize
672 y = lines.thetas / self.config.thetaBinSize
673 X = np.array([x, y]).T
674 nClusters = 1
675
676 # Put line parameters in clusters by starting with all in one, then
677 # subdividing until the parameters of each cluster have std dev=1.
678 # If nClusters == len(lines), each line will have its own 'cluster', so
679 # the standard deviations of each cluster must be zero and the loop
680 # is guaranteed to stop.
681 while True:
682 kmeans = KMeans(n_clusters=nClusters, n_init='auto').fit(X)
683 clusterStandardDeviations = np.zeros((nClusters, 2))
684 for c in range(nClusters):
685 inCluster = X[kmeans.labels_ == c]
686 clusterStandardDeviations[c] = np.std(inCluster, axis=0)
687 # Are the rhos and thetas in each cluster all below the threshold?
688 if (clusterStandardDeviations <= 1).all():
689 break
690 nClusters += 1
691
692 # The cluster centers are final line estimates
693 finalClusters = kmeans.cluster_centers_.T
694
695 # Rescale variables:
696 finalRhos = finalClusters[0] * self.config.rhoBinSize
697 finalThetas = finalClusters[1] * self.config.thetaBinSize
698 result = LineCollection(finalRhos, finalThetas)
699 self.log.info("Lines were grouped into %s potential streak(s)", len(finalRhos))
700
701 return result
702
703 def _fitProfile(self, lines, maskedImage):
704 """Fit the profile of the streak.
705
706 Given the initial parameters of detected lines, fit a model for the
707 streak to the original (non-binary image). The assumed model is a
708 straight line with a Moffat profile.
709
710 Parameters
711 ----------
712 lines : `LineCollection`
713 Collection of guesses for `Line`s detected in the image.
714 maskedImage : `lsst.afw.image.maskedImage`
715 Original image to be used to fit profile of streak.
716
717 Returns
718 -------
719 lineFits : `LineCollection`
720 Collection of `Line` profiles fit to the data.
721 finalMask : `np.ndarray`
722 2d mask array with detected streaks=1.
723 """
724 data = maskedImage.image.array
725 weights = maskedImage.variance.array**-1
726 mask = maskedImage.mask
727 badPixelMask = mask.getPlaneBitMask(self.config.badMaskPlanes)
728 badMask = (mask.array & badPixelMask) > 0
729 # Mask out any pixels with non-finite weights
730 weights[~np.isfinite(weights) | ~np.isfinite(data)] = 0
731 weights[badMask] = 0
732
733 lineFits = LineCollection([], [])
734 finalLineMasks = [np.zeros(data.shape, dtype=bool)]
735 nFinalLines = 0
736 nFitFailures = 0
737 for line in lines:
738 line.sigma = self.config.invSigma**-1
739 lineModel = LineProfile(data, weights, line=line)
740 # Skip any lines that do not cover any data (sometimes happens because of chip gaps)
741 if lineModel.lineMaskSize == 0:
742 continue
743
744 fit, fitFailure = lineModel.fit(dChi2Tol=self.config.dChi2Tolerance, log=self.log)
745
746 # Initial estimate should be quite close: fit is deemed unsuccessful if rho or theta
747 # change more than the allowed bin in rho or theta:
748 if ((abs(fit.rho - line.rho) > 2 * self.config.rhoBinSize)
749 or (abs(fit.theta - line.theta) > 2 * self.config.thetaBinSize)):
750 fitFailure = True
751 self.log.debug("Streak fit moved too far from initial estimate. Line will be dropped.")
752
753 if fitFailure:
754 nFitFailures += 1
755 continue
756
757 # Make mask
758 lineModel.setLineMask(fit)
759 finalModel = lineModel.makeProfile(fit)
760 # Take absolute value, as streaks are allowed to be negative
761 finalModelMax = abs(finalModel).max()
762 finalLineMask = abs(finalModel) > self.config.footprintThreshold
763 # Drop this line if the model profile is below the footprint
764 # threshold
765 if not finalLineMask.any():
766 continue
767 fit.modelMaximum = finalModelMax
768 lineFits.append(fit)
769 finalLineMasks.append(finalLineMask)
770 nFinalLines += 1
771
772 if nFitFailures > 0:
773 self.log.info("Streak profile could not be fit for %d out of %d detected lines.", nFitFailures,
774 len(lines))
775 finalMask = np.array(finalLineMasks).any(axis=0)
776
777 return lineFits, finalMask
int max
__init__(self, rhos, thetas, sigmas=None)
__init__(self, data, weights, line=None)