lsst.pipe.tasks g4a6547c0d5+de68eba77a
Loading...
Searching...
No Matches
maskStreaks.py
Go to the documentation of this file.
2# LSST Data Management System
3# Copyright 2014 LSST Corporation.
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 <http://www.lsstcorp.org/LegalNotices/>.
21#
22
23__all__ = ["MaskStreaksConfig", "MaskStreaksTask", "setDetectionMask"]
24
25import lsst.pex.config as pexConfig
26import lsst.pipe.base as pipeBase
27import lsst.kht
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 rho: float
132 theta: float
133 sigma: float = 0
134
135
137 """Collection of `Line` objects.
138
139 Parameters
140 ----------
141 rhos : np.ndarray
142 Array of `Line` rho parameters
143 thetas : np.ndarray
144 Array of `Line` theta parameters
145 sigmas : np.ndarray (optional)
146 Array of `Line` sigma parameters
147 """
148
149 def __init__(self, rhos, thetas, sigmas=None):
150 if sigmas is None:
151 sigmas = np.zeros(len(rhos))
152
153 self._lines = [Line(rho, theta, sigma) for (rho, theta, sigma) in
154 zip(rhos, thetas, sigmas)]
155
156 def __len__(self):
157 return len(self._lines)
158
159 def __getitem__(self, index):
160 return self._lines[index]
161
162 def __iter__(self):
163 return iter(self._lines)
164
165 def __repr__(self):
166 joinedString = ", ".join(str(line) for line in self._lines)
167 return textwrap.shorten(joinedString, width=160, placeholder="...")
168
169 @property
170 def rhos(self):
171 return np.array([line.rho for line in self._lines])
172
173 @property
174 def thetas(self):
175 return np.array([line.theta for line in self._lines])
176
177 def append(self, newLine):
178 """Add line to current collection of lines.
179
180 Parameters
181 ----------
182 newLine : `Line`
183 `Line` to add to current collection of lines
184 """
185 self._lines.append(copy.copy(newLine))
186
187
189 """Construct and/or fit a model for a linear streak.
190
191 This assumes a simple model for a streak, in which the streak
192 follows a straight line in pixels space, with a Moffat-shaped profile. The
193 model is fit to data using a Newton-Raphson style minimization algorithm.
194 The initial guess for the line parameters is assumed to be fairly accurate,
195 so only a narrow band of pixels around the initial line estimate is used in
196 fitting the model, which provides a significant speed-up over using all the
197 data. The class can also be used just to construct a model for the data with
198 a line following the given coordinates.
199
200 Parameters
201 ----------
202 data : np.ndarray
203 2d array of data
204 weights : np.ndarray
205 2d array of weights
206 line : `Line` (optional)
207 Guess for position of line. Data far from line guess is masked out.
208 Defaults to None, in which case only data with `weights` = 0 is masked
209 out.
210 """
211
212 def __init__(self, data, weights, line=None):
213 self.data = data
214 self.weights = weights
215 self._ymax, self._xmax = data.shape
216 self._dtype = data.dtype
217 xrange = np.arange(self._xmax) - self._xmax / 2.
218 yrange = np.arange(self._ymax) - self._ymax / 2.
219 self._rhoMax = ((0.5 * self._ymax)**2 + (0.5 * self._xmax)**2)**0.5
220 self._xmesh, self._ymesh = np.meshgrid(xrange, yrange)
221 self.mask = (weights != 0)
222
223 self._initLine = line
224 self.setLineMask(line)
225
226 def setLineMask(self, line):
227 """Set mask around the image region near the line
228
229 Parameters
230 ----------
231 line : `Line`
232 Parameters of line in the image
233 """
234 if line:
235 # Only fit pixels within 5 sigma of the estimated line
236 radtheta = np.deg2rad(line.theta)
237 distance = (np.cos(radtheta) * self._xmesh + np.sin(radtheta) * self._ymesh - line.rho)
238 m = (abs(distance) < 5 * line.sigma)
239 self.lineMask = self.mask & m
240 else:
241 self.lineMask = np.copy(self.mask)
242
243 self.lineMaskSize = self.lineMask.sum()
244 self._maskData = self.data[self.lineMask]
245 self._maskWeights = self.weights[self.lineMask]
246 self._mxmesh = self._xmesh[self.lineMask]
247 self._mymesh = self._ymesh[self.lineMask]
248
249 def _makeMaskedProfile(self, line, fitFlux=True):
250 """Construct the line model in the masked region and calculate its
251 derivatives
252
253 Parameters
254 ----------
255 line : `Line`
256 Parameters of line profile for which to make profile in the masked
257 region
258 fitFlux : bool
259 Fit the amplitude of the line profile to the data
260
261 Returns
262 -------
263 model : np.ndarray
264 Model in the masked region
265 dModel : np.ndarray
266 Derivative of the model in the masked region
267 """
268 invSigma = line.sigma**-1
269 # Calculate distance between pixels and line
270 radtheta = np.deg2rad(line.theta)
271 costheta = np.cos(radtheta)
272 sintheta = np.sin(radtheta)
273 distance = (costheta * self._mxmesh + sintheta * self._mymesh - line.rho)
274 distanceSquared = distance**2
275
276 # Calculate partial derivatives of distance
277 drad = np.pi / 180
278 dDistanceSqdRho = 2 * distance * (-np.ones_like(self._mxmesh))
279 dDistanceSqdTheta = (2 * distance * (-sintheta * self._mxmesh + costheta * self._mymesh) * drad)
280
281 # Use pixel-line distances to make Moffat profile
282 profile = (1 + distanceSquared * invSigma**2)**-2.5
283 dProfile = -2.5 * (1 + distanceSquared * invSigma**2)**-3.5
284
285 if fitFlux:
286 # Calculate line flux from profile and data
287 flux = ((self._maskWeights * self._maskData * profile).sum()
288 / (self._maskWeights * profile**2).sum())
289 else:
290 # Approximately normalize the line
291 flux = invSigma**-1
292 if np.isnan(flux):
293 flux = 0
294
295 model = flux * profile
296
297 # Calculate model derivatives
298 fluxdProfile = flux * dProfile
299 fluxdProfileInvSigma = fluxdProfile * invSigma**2
300 dModeldRho = fluxdProfileInvSigma * dDistanceSqdRho
301 dModeldTheta = fluxdProfileInvSigma * dDistanceSqdTheta
302 dModeldInvSigma = fluxdProfile * distanceSquared * 2 * invSigma
303
304 dModel = np.array([dModeldRho, dModeldTheta, dModeldInvSigma])
305 return model, dModel
306
307 def makeProfile(self, line, fitFlux=True):
308 """Construct the line profile model
309
310 Parameters
311 ----------
312 line : `Line`
313 Parameters of the line profile to model
314 fitFlux : bool (optional)
315 Fit the amplitude of the line profile to the data
316
317 Returns
318 -------
319 finalModel : np.ndarray
320 Model for line profile
321 """
322 model, _ = self._makeMaskedProfile(line, fitFlux=fitFlux)
323 finalModel = np.zeros((self._ymax, self._xmax), dtype=self._dtype)
324 finalModel[self.lineMask] = model
325 return finalModel
326
327 def _lineChi2(self, line, grad=True):
328 """Construct the chi2 between the data and the model
329
330 Parameters
331 ----------
332 line : `Line`
333 `Line` parameters for which to build model and calculate chi2
334 grad : bool (optional)
335 Whether or not to return the gradient and hessian
336
337 Returns
338 -------
339 reducedChi : float
340 Reduced chi2 of the model
341 reducedDChi : np.ndarray
342 Derivative of the chi2 with respect to rho, theta, invSigma
343 reducedHessianChi : np.ndarray
344 Hessian of the chi2 with respect to rho, theta, invSigma
345 """
346 # Calculate chi2
347 model, dModel = self._makeMaskedProfile(line)
348 chi2 = (self._maskWeights * (self._maskData - model)**2).sum()
349 if not grad:
350 return chi2.sum() / self.lineMaskSize
351
352 # Calculate derivative and Hessian of chi2
353 derivChi2 = ((-2 * self._maskWeights * (self._maskData - model))[None, :] * dModel).sum(axis=1)
354 hessianChi2 = (2 * self._maskWeights * dModel[:, None, :] * dModel[None, :, :]).sum(axis=2)
355
356 reducedChi = chi2 / self.lineMaskSize
357 reducedDChi = derivChi2 / self.lineMaskSize
358 reducedHessianChi = hessianChi2 / self.lineMaskSize
359 return reducedChi, reducedDChi, reducedHessianChi
360
361 def fit(self, dChi2Tol=0.1, maxIter=100):
362 """Perform Newton-Raphson minimization to find line parameters
363
364 This method takes advantage of having known derivative and Hessian of
365 the multivariate function to quickly and efficiently find the minimum.
366 This is more efficient than the scipy implementation of the Newton-
367 Raphson method, which doesn't take advantage of the Hessian matrix. The
368 method here also performs a line search in the direction of the steepest
369 derivative at each iteration, which reduces the number of iterations
370 needed.
371
372 Parameters
373 ----------
374 dChi2Tol : float (optional)
375 Change in Chi2 tolerated for fit convergence
376 maxIter : int (optional)
377 Maximum number of fit iterations allowed. The fit should converge in
378 ~10 iterations, depending on the value of dChi2Tol, but this
379 maximum provides a backup.
380
381 Returns
382 -------
383 outline : np.ndarray
384 Coordinates and inverse width of fit line
385 chi2 : float
386 Reduced Chi2 of model fit to data
387 fitFailure : bool
388 Boolean where `False` corresponds to a successful fit
389 """
390 # Do minimization on inverse of sigma to simplify derivatives:
391 x = np.array([self._initLine.rho, self._initLine.theta, self._initLine.sigma**-1])
392
393 dChi2 = 1
394 iter = 0
395 oldChi2 = 0
396 fitFailure = False
397
398 def line_search(c, dx):
399 testx = x - c * dx
400 testLine = Line(testx[0], testx[1], testx[2]**-1)
401 return self._lineChi2(testLine, grad=False)
402
403 while abs(dChi2) > dChi2Tol:
404 line = Line(x[0], x[1], x[2]**-1)
405 chi2, b, A = self._lineChi2(line)
406 if chi2 == 0:
407 break
408 if not np.isfinite(A).all():
409 # TODO: DM-30797 Add warning here.
410 fitFailure = True
411 break
412 dChi2 = oldChi2 - chi2
413 cholesky = scipy.linalg.cho_factor(A)
414 dx = scipy.linalg.cho_solve(cholesky, b)
415
416 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05)
417 x -= factor * dx
418 if (abs(x[0]) > 1.5 * self._rhoMax) or (iter > maxIter):
419 fitFailure = True
420 break
421 oldChi2 = chi2
422 iter += 1
423
424 outline = Line(x[0], x[1], abs(x[2])**-1)
425
426 return outline, chi2, fitFailure
427
428
429class MaskStreaksConfig(pexConfig.Config):
430 """Configuration parameters for `MaskStreaksTask`
431 """
432 minimumKernelHeight = pexConfig.Field(
433 doc="Minimum height of the streak-finding kernel relative to the tallest kernel",
434 dtype=float,
435 default=0.0,
436 )
437 absMinimumKernelHeight = pexConfig.Field(
438 doc="Minimum absolute height of the streak-finding kernel",
439 dtype=float,
440 default=5,
441 )
442 clusterMinimumSize = pexConfig.Field(
443 doc="Minimum size in pixels of detected clusters",
444 dtype=int,
445 default=50,
446 )
447 clusterMinimumDeviation = pexConfig.Field(
448 doc="Allowed deviation (in pixels) from a straight line for a detected "
449 "line",
450 dtype=int,
451 default=2,
452 )
453 delta = pexConfig.Field(
454 doc="Stepsize in angle-radius parameter space",
455 dtype=float,
456 default=0.2,
457 )
458 nSigma = pexConfig.Field(
459 doc="Number of sigmas from center of kernel to include in voting "
460 "procedure",
461 dtype=float,
462 default=2,
463 )
464 rhoBinSize = pexConfig.Field(
465 doc="Binsize in pixels for position parameter rho when finding "
466 "clusters of detected lines",
467 dtype=float,
468 default=30,
469 )
470 thetaBinSize = pexConfig.Field(
471 doc="Binsize in degrees for angle parameter theta when finding "
472 "clusters of detected lines",
473 dtype=float,
474 default=2,
475 )
476 invSigma = pexConfig.Field(
477 doc="Inverse of the Moffat sigma parameter (in units of pixels)"
478 "describing the profile of the streak",
479 dtype=float,
480 default=10.**-1,
481 )
482 footprintThreshold = pexConfig.Field(
483 doc="Threshold at which to determine edge of line, in units of the line"
484 "profile maximum",
485 dtype=float,
486 default=0.01
487 )
488 dChi2Tolerance = pexConfig.Field(
489 doc="Absolute difference in Chi2 between iterations of line profile"
490 "fitting that is acceptable for convergence",
491 dtype=float,
492 default=0.1
493 )
494 detectedMaskPlane = pexConfig.Field(
495 doc="Name of mask with pixels above detection threshold, used for first"
496 "estimate of streak locations",
497 dtype=str,
498 default="DETECTED"
499 )
500 streaksMaskPlane = pexConfig.Field(
501 doc="Name of mask plane holding detected streaks",
502 dtype=str,
503 default="STREAK"
504 )
505
506
507class MaskStreaksTask(pipeBase.Task):
508 """Find streaks or other straight lines in image data.
509
510 Nearby objects passing through the field of view of the telescope leave a
511 bright trail in images. This class uses the Kernel Hough Transform (KHT)
512 (Fernandes and Oliveira, 2007), implemented in `lsst.houghtransform`. The
513 procedure works by taking a binary image, either provided as put or produced
514 from the input data image, using a Canny filter to make an image of the
515 edges in the original image, then running the KHT on the edge image. The KHT
516 identifies clusters of non-zero points, breaks those clusters of points into
517 straight lines, keeps clusters with a size greater than the user-set
518 threshold, then performs a voting procedure to find the best-fit coordinates
519 of any straight lines. Given the results of the KHT algorithm, clusters of
520 lines are identified and grouped (generally these correspond to the two
521 edges of a strea) and a profile is fit to the streak in the original
522 (non-binary) image.
523 """
524
525 ConfigClass = MaskStreaksConfig
526 _DefaultName = "maskStreaks"
527
528 @pipeBase.timeMethod
529 def find(self, maskedImage):
530 """Find streaks in a masked image
531
532 Parameters
533 ----------
534 maskedImage : `lsst.afw.image.maskedImage`
535 The image in which to search for streaks.
536
537 Returns
538 -------
539 result : `lsst.pipe.base.Struct`
540 Result struct with components:
541
542 - ``originalLines``: lines identified by kernel hough transform
543 - ``lineClusters``: lines grouped into clusters in rho-theta space
544 - ``lines``: final result for lines after line-profile fit
545 - ``mask``: 2-d boolean mask where detected lines are True
546 """
547 mask = maskedImage.getMask()
548 detectionMask = (mask.array & mask.getPlaneBitMask(self.config.detectedMaskPlane))
549
550 self.edges = self._cannyFilter(detectionMask)
551 self.lines = self._runKHT(self.edges)
552
553 if len(self.lines) == 0:
554 lineMask = np.zeros(detectionMask.shape, dtype=bool)
555 fitLines = LineCollection([], [])
556 clusters = LineCollection([], [])
557 else:
558 clusters = self._findClusters(self.lines)
559 fitLines, lineMask = self._fitProfile(clusters, maskedImage)
560
561 # The output mask is the intersection of the fit streaks and the image detections
562 outputMask = lineMask & detectionMask.astype(bool)
563
564 return pipeBase.Struct(
565 lines=fitLines,
566 lineClusters=clusters,
567 originalLines=self.lines,
568 mask=outputMask,
569 )
570
571 @pipeBase.timeMethod
572 def run(self, maskedImage):
573 """Find and mask streaks in a masked image.
574
575 Finds streaks in the image and modifies maskedImage in place by adding a
576 mask plane with any identified streaks.
577
578 Parameters
579 ----------
580 maskedImage : `lsst.afw.image.maskedImage`
581 The image in which to search for streaks. The mask detection plane
582 corresponding to `config.detectedMaskPlane` must be set with the
583 detected pixels.
584
585 Returns
586 -------
587 result : `lsst.pipe.base.Struct`
588 Result struct with components:
589
590 - ``originalLines``: lines identified by kernel hough transform
591 - ``lineClusters``: lines grouped into clusters in rho-theta space
592 - ``lines``: final result for lines after line-profile fit
593 """
594 streaks = self.find(maskedImage)
595
596 maskedImage.mask.addMaskPlane(self.config.streaksMaskPlane)
597 maskedImage.mask.array[streaks.mask] |= maskedImage.mask.getPlaneBitMask(self.config.streaksMaskPlane)
598
599 return pipeBase.Struct(
600 lines=streaks.lines,
601 lineClusters=streaks.lineClusters,
602 originalLines=streaks.originalLines,
603 )
604
605 def _cannyFilter(self, image):
606 """Apply a canny filter to the data in order to detect edges
607
608 Parameters
609 ----------
610 image : `np.ndarray`
611 2-d image data on which to run filter
612
613 Returns
614 -------
615 cannyData : `np.ndarray`
616 2-d image of edges found in input image
617 """
618 filterData = image.astype(int)
619 return canny(filterData, low_threshold=0, high_threshold=1, sigma=0.1)
620
621 def _runKHT(self, image):
622 """Run Kernel Hough Transform on image.
623
624 Parameters
625 ----------
626 image : `np.ndarray`
627 2-d image data on which to detect lines
628
629 Returns
630 -------
631 result : `LineCollection`
632 Collection of detected lines, with their detected rho and theta
633 coordinates.
634 """
635 lines = lsst.kht.find_lines(image, self.config.clusterMinimumSize,
636 self.config.clusterMinimumDeviation, self.config.delta,
637 self.config.minimumKernelHeight, self.config.nSigma,
638 self.config.absMinimumKernelHeight)
639
640 return LineCollection(lines.rho, lines.theta)
641
642 def _findClusters(self, lines):
643 """Group lines that are close in parameter space and likely describe
644 the same streak.
645
646 Parameters
647 ----------
648 lines : `LineCollection`
649 Collection of lines to group into clusters
650
651 Returns
652 -------
653 result : `LineCollection`
654 Average `Line` for each cluster of `Line`s in the input
655 `LineCollection`
656 """
657 # Scale variables by threshold bin-size variable so that rho and theta
658 # are on the same scale. Since the clustering algorithm below stops when
659 # the standard deviation <= 1, after rescaling each cluster will have a
660 # standard deviation at or below the bin-size.
661 x = lines.rhos / self.config.rhoBinSize
662 y = lines.thetas / self.config.thetaBinSize
663 X = np.array([x, y]).T
664 nClusters = 1
665
666 # Put line parameters in clusters by starting with all in one, then
667 # subdividing until the parameters of each cluster have std dev=1.
668 # If nClusters == len(lines), each line will have its own 'cluster', so
669 # the standard deviations of each cluster must be zero and the loop
670 # is guaranteed to stop.
671 while True:
672 kmeans = KMeans(n_clusters=nClusters).fit(X)
673 clusterStandardDeviations = np.zeros((nClusters, 2))
674 for c in range(nClusters):
675 inCluster = X[kmeans.labels_ == c]
676 clusterStandardDeviations[c] = np.std(inCluster, axis=0)
677 # Are the rhos and thetas in each cluster all below the threshold?
678 if (clusterStandardDeviations <= 1).all():
679 break
680 nClusters += 1
681
682 # The cluster centers are final line estimates
683 finalClusters = kmeans.cluster_centers_.T
684
685 # Rescale variables:
686 finalRhos = finalClusters[0] * self.config.rhoBinSize
687 finalThetas = finalClusters[1] * self.config.thetaBinSize
688 result = LineCollection(finalRhos, finalThetas)
689
690 return result
691
692 def _fitProfile(self, lines, maskedImage):
693 """Fit the profile of the streak.
694
695 Given the initial parameters of detected lines, fit a model for the
696 streak to the original (non-binary image). The assumed model is a
697 straight line with a Moffat profile.
698
699 Parameters
700 ----------
701 lines : `LineCollection`
702 Collection of guesses for `Line`s detected in the image
703 maskedImage : `lsst.afw.image.maskedImage`
704 Original image to be used to fit profile of streak.
705
706 Returns
707 -------
708 lineFits : `LineCollection`
709 Collection of `Line` profiles fit to the data
710 finalMask : `np.ndarray`
711 2d mask array with detected streaks=1.
712 """
713 data = maskedImage.image.array
714 weights = maskedImage.variance.array**-1
715 # Mask out any pixels with non-finite weights
716 weights[~np.isfinite(weights) | ~np.isfinite(data)] = 0
717
718 lineFits = LineCollection([], [])
719 finalLineMasks = [np.zeros(data.shape, dtype=bool)]
720 for line in lines:
721 line.sigma = self.config.invSigma**-1
722 lineModel = LineProfile(data, weights, line=line)
723 # Skip any lines that do not cover any data (sometimes happens because of chip gaps)
724 if lineModel.lineMaskSize == 0:
725 continue
726
727 fit, chi2, fitFailure = lineModel.fit(dChi2Tol=self.config.dChi2Tolerance)
728
729 # Initial estimate should be quite close: fit is deemed unsuccessful if rho or theta
730 # change more than the allowed bin in rho or theta:
731 if ((abs(fit.rho - line.rho) > 2 * self.config.rhoBinSize)
732 or (abs(fit.theta - line.theta) > 2 * self.config.thetaBinSize)):
733 fitFailure = True
734
735 if fitFailure:
736 continue
737
738 # Make mask
739 lineModel.setLineMask(fit)
740 finalModel = lineModel.makeProfile(fit)
741 # Take absolute value, as streaks are allowed to be negative
742 finalModelMax = abs(finalModel).max()
743 finalLineMask = abs(finalModel) > self.config.footprintThreshold
744 # Drop this line if the model profile is below the footprint threshold
745 if not finalLineMask.any():
746 continue
747 fit.chi2 = chi2
748 fit.finalModelMax = finalModelMax
749 lineFits.append(fit)
750 finalLineMasks.append(finalLineMask)
751
752 finalMask = np.array(finalLineMasks).any(axis=0)
753
754 return lineFits, finalMask
def __init__(self, rhos, thetas, sigmas=None)
Definition: maskStreaks.py:149
def __init__(self, data, weights, line=None)
Definition: maskStreaks.py:212
def makeProfile(self, line, fitFlux=True)
Definition: maskStreaks.py:307
def _makeMaskedProfile(self, line, fitFlux=True)
Definition: maskStreaks.py:249
def _lineChi2(self, line, grad=True)
Definition: maskStreaks.py:327
def fit(self, dChi2Tol=0.1, maxIter=100)
Definition: maskStreaks.py:361
def _fitProfile(self, lines, maskedImage)
Definition: maskStreaks.py:692
def setDetectionMask(maskedImage, forceSlowBin=False, binning=None, detectedPlane="DETECTED", badMaskPlanes=("NO_DATA", "INTRP", "BAD", "SAT", "EDGE"), detectionThreshold=5)
Definition: maskStreaks.py:40