Hide keyboard shortcuts

Hot-keys on this page

r m x p   toggle line displays

j k   next/prev highlighted chunk

0   (zero) top of page

1   (one) first highlighted chunk

1# 

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 

136class LineCollection: 

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: 150 ↛ 153line 150 didn't jump to line 153, because the condition on line 150 was never false

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 

188class LineProfile: 

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: 234 ↛ 241line 234 didn't jump to line 241, because the condition on line 234 was never false

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): 292 ↛ 293line 292 didn't jump to line 293, because the condition on line 292 was never true

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: 406 ↛ 407line 406 didn't jump to line 407, because the condition on line 406 was never true

407 break 

408 dChi2 = oldChi2 - chi2 

409 cholesky = scipy.linalg.cho_factor(A) 

410 dx = scipy.linalg.cho_solve(cholesky, b) 

411 

412 factor, fmin, _, _ = scipy.optimize.brent(line_search, args=(dx,), full_output=True, tol=0.05) 

413 x -= factor * dx 

414 if (x[0] > 1.5 * self._rhoMax) or (iter > maxIter): 414 ↛ 415line 414 didn't jump to line 415, because the condition on line 414 was never true

415 fitFailure = True 

416 break 

417 oldChi2 = chi2 

418 iter += 1 

419 

420 outline = Line(x[0], x[1], abs(x[2])**-1) 

421 

422 return outline, chi2, fitFailure 

423 

424 

425class MaskStreaksConfig(pexConfig.Config): 

426 """Configuration parameters for `MaskStreaksTask` 

427 """ 

428 minimumKernelHeight = pexConfig.Field( 

429 doc="Minimum height of the streak-finding kernel relative to the tallest kernel", 

430 dtype=float, 

431 default=0.0, 

432 ) 

433 absMinimumKernelHeight = pexConfig.Field( 

434 doc="Minimum absolute height of the streak-finding kernel", 

435 dtype=float, 

436 default=5, 

437 ) 

438 clusterMinimumSize = pexConfig.Field( 

439 doc="Minimum size in pixels of detected clusters", 

440 dtype=int, 

441 default=50, 

442 ) 

443 clusterMinimumDeviation = pexConfig.Field( 

444 doc="Allowed deviation (in pixels) from a straight line for a detected " 

445 "line", 

446 dtype=int, 

447 default=2, 

448 ) 

449 delta = pexConfig.Field( 

450 doc="Stepsize in angle-radius parameter space", 

451 dtype=float, 

452 default=0.2, 

453 ) 

454 nSigma = pexConfig.Field( 

455 doc="Number of sigmas from center of kernel to include in voting " 

456 "procedure", 

457 dtype=float, 

458 default=2, 

459 ) 

460 rhoBinSize = pexConfig.Field( 

461 doc="Binsize in pixels for position parameter rho when finding " 

462 "clusters of detected lines", 

463 dtype=float, 

464 default=30, 

465 ) 

466 thetaBinSize = pexConfig.Field( 

467 doc="Binsize in degrees for angle parameter theta when finding " 

468 "clusters of detected lines", 

469 dtype=float, 

470 default=2, 

471 ) 

472 invSigma = pexConfig.Field( 

473 doc="Inverse of the Moffat sigma parameter (in units of pixels)" 

474 "describing the profile of the streak", 

475 dtype=float, 

476 default=10.**-1, 

477 ) 

478 footprintThreshold = pexConfig.Field( 

479 doc="Threshold at which to determine edge of line, in units of the line" 

480 "profile maximum", 

481 dtype=float, 

482 default=0.01 

483 ) 

484 dChi2Tolerance = pexConfig.Field( 

485 doc="Absolute difference in Chi2 between iterations of line profile" 

486 "fitting that is acceptable for convergence", 

487 dtype=float, 

488 default=0.1 

489 ) 

490 detectedMaskPlane = pexConfig.Field( 

491 doc="Name of mask with pixels above detection threshold, used for first" 

492 "estimate of streak locations", 

493 dtype=str, 

494 default="DETECTED" 

495 ) 

496 streaksMaskPlane = pexConfig.Field( 

497 doc="Name of mask plane holding detected streaks", 

498 dtype=str, 

499 default="STREAK" 

500 ) 

501 

502 

503class MaskStreaksTask(pipeBase.Task): 

504 """Find streaks or other straight lines in image data. 

505 

506 Nearby objects passing through the field of view of the telescope leave a 

507 bright trail in images. This class uses the Kernel Hough Transform (KHT) 

508 (Fernandes and Oliveira, 2007), implemented in `lsst.houghtransform`. The 

509 procedure works by taking a binary image, either provided as put or produced 

510 from the input data image, using a Canny filter to make an image of the 

511 edges in the original image, then running the KHT on the edge image. The KHT 

512 identifies clusters of non-zero points, breaks those clusters of points into 

513 straight lines, keeps clusters with a size greater than the user-set 

514 threshold, then performs a voting procedure to find the best-fit coordinates 

515 of any straight lines. Given the results of the KHT algorithm, clusters of 

516 lines are identified and grouped (generally these correspond to the two 

517 edges of a strea) and a profile is fit to the streak in the original 

518 (non-binary) image. 

519 """ 

520 

521 ConfigClass = MaskStreaksConfig 

522 _DefaultName = "maskStreaks" 

523 

524 @pipeBase.timeMethod 

525 def find(self, maskedImage): 

526 """Find streaks in a masked image 

527 

528 Parameters 

529 ---------- 

530 maskedImage : `lsst.afw.image.maskedImage` 

531 The image in which to search for streaks. 

532 

533 Returns 

534 ------- 

535 result : `lsst.pipe.base.Struct` 

536 Result struct with components: 

537 

538 - ``originalLines``: lines identified by kernel hough transform 

539 - ``lineClusters``: lines grouped into clusters in rho-theta space 

540 - ``lines``: final result for lines after line-profile fit 

541 - ``mask``: 2-d boolean mask where detected lines are True 

542 """ 

543 mask = maskedImage.getMask() 

544 detectionMask = (mask.array & mask.getPlaneBitMask(self.config.detectedMaskPlane)) 

545 

546 self.edges = self._cannyFilter(detectionMask) 

547 self.lines = self._runKHT(self.edges) 

548 

549 if len(self.lines) == 0: 

550 lineMask = np.zeros(detectionMask.shape, dtype=bool) 

551 fitLines = LineCollection([], []) 

552 clusters = LineCollection([], []) 

553 else: 

554 clusters = self._findClusters(self.lines) 

555 fitLines, lineMask = self._fitProfile(clusters, maskedImage) 

556 

557 # The output mask is the intersection of the fit streaks and the image detections 

558 outputMask = lineMask & detectionMask.astype(bool) 

559 

560 return pipeBase.Struct( 

561 lines=fitLines, 

562 lineClusters=clusters, 

563 originalLines=self.lines, 

564 mask=outputMask, 

565 ) 

566 

567 @pipeBase.timeMethod 

568 def run(self, maskedImage): 

569 """Find and mask streaks in a masked image. 

570 

571 Finds streaks in the image and modifies maskedImage in place by adding a 

572 mask plane with any identified streaks. 

573 

574 Parameters 

575 ---------- 

576 maskedImage : `lsst.afw.image.maskedImage` 

577 The image in which to search for streaks. The mask detection plane 

578 corresponding to `config.detectedMaskPlane` must be set with the 

579 detected pixels. 

580 

581 Returns 

582 ------- 

583 result : `lsst.pipe.base.Struct` 

584 Result struct with components: 

585 

586 - ``originalLines``: lines identified by kernel hough transform 

587 - ``lineClusters``: lines grouped into clusters in rho-theta space 

588 - ``lines``: final result for lines after line-profile fit 

589 """ 

590 streaks = self.find(maskedImage) 

591 

592 maskedImage.mask.addMaskPlane(self.config.streaksMaskPlane) 

593 maskedImage.mask.array[streaks.mask] |= maskedImage.mask.getPlaneBitMask(self.config.streaksMaskPlane) 

594 

595 return pipeBase.Struct( 

596 lines=streaks.lines, 

597 lineClusters=streaks.lineClusters, 

598 originalLines=streaks.originalLines, 

599 ) 

600 

601 def _cannyFilter(self, image): 

602 """Apply a canny filter to the data in order to detect edges 

603 

604 Parameters 

605 ---------- 

606 image : `np.ndarray` 

607 2-d image data on which to run filter 

608 

609 Returns 

610 ------- 

611 cannyData : `np.ndarray` 

612 2-d image of edges found in input image 

613 """ 

614 filterData = image.astype(int) 

615 return canny(filterData, low_threshold=0, high_threshold=1, sigma=0.1) 

616 

617 def _runKHT(self, image): 

618 """Run Kernel Hough Transform on image. 

619 

620 Parameters 

621 ---------- 

622 image : `np.ndarray` 

623 2-d image data on which to detect lines 

624 

625 Returns 

626 ------- 

627 result : `LineCollection` 

628 Collection of detected lines, with their detected rho and theta 

629 coordinates. 

630 """ 

631 lines = lsst.kht.find_lines(image, self.config.clusterMinimumSize, 

632 self.config.clusterMinimumDeviation, self.config.delta, 

633 self.config.minimumKernelHeight, self.config.nSigma, 

634 self.config.absMinimumKernelHeight) 

635 

636 return LineCollection(lines.rho, lines.theta) 

637 

638 def _findClusters(self, lines): 

639 """Group lines that are close in parameter space and likely describe 

640 the same streak. 

641 

642 Parameters 

643 ---------- 

644 lines : `LineCollection` 

645 Collection of lines to group into clusters 

646 

647 Returns 

648 ------- 

649 result : `LineCollection` 

650 Average `Line` for each cluster of `Line`s in the input 

651 `LineCollection` 

652 """ 

653 # Scale variables by threshold bin-size variable so that rho and theta 

654 # are on the same scale. Since the clustering algorithm below stops when 

655 # the standard deviation <= 1, after rescaling each cluster will have a 

656 # standard deviation at or below the bin-size. 

657 x = lines.rhos / self.config.rhoBinSize 

658 y = lines.thetas / self.config.thetaBinSize 

659 X = np.array([x, y]).T 

660 nClusters = 1 

661 

662 # Put line parameters in clusters by starting with all in one, then 

663 # subdividing until the parameters of each cluster have std dev=1. 

664 # If nClusters == len(lines), each line will have its own 'cluster', so 

665 # the standard deviations of each cluster must be zero and the loop 

666 # is guaranteed to stop. 

667 while True: 

668 kmeans = KMeans(n_clusters=nClusters).fit(X) 

669 clusterStandardDeviations = np.zeros((nClusters, 2)) 

670 for c in range(nClusters): 

671 inCluster = X[kmeans.labels_ == c] 

672 clusterStandardDeviations[c] = np.std(inCluster, axis=0) 

673 # Are the rhos and thetas in each cluster all below the threshold? 

674 if (clusterStandardDeviations <= 1).all(): 674 ↛ 676line 674 didn't jump to line 676, because the condition on line 674 was never false

675 break 

676 nClusters += 1 

677 

678 # The cluster centers are final line estimates 

679 finalClusters = kmeans.cluster_centers_.T 

680 

681 # Rescale variables: 

682 finalRhos = finalClusters[0] * self.config.rhoBinSize 

683 finalThetas = finalClusters[1] * self.config.thetaBinSize 

684 result = LineCollection(finalRhos, finalThetas) 

685 

686 return result 

687 

688 def _fitProfile(self, lines, maskedImage): 

689 """Fit the profile of the streak. 

690 

691 Given the initial parameters of detected lines, fit a model for the 

692 streak to the original (non-binary image). The assumed model is a 

693 straight line with a Moffat profile. 

694 

695 Parameters 

696 ---------- 

697 lines : `LineCollection` 

698 Collection of guesses for `Line`s detected in the image 

699 maskedImage : `lsst.afw.image.maskedImage` 

700 Original image to be used to fit profile of streak. 

701 

702 Returns 

703 ------- 

704 lineFits : `LineCollection` 

705 Collection of `Line` profiles fit to the data 

706 finalMask : `np.ndarray` 

707 2d mask array with detected streaks=1. 

708 """ 

709 data = maskedImage.image.array 

710 weights = maskedImage.variance.array**-1 

711 

712 lineFits = LineCollection([], []) 

713 finalLineMasks = [np.zeros(data.shape, dtype=bool)] 

714 for line in lines: 

715 line.sigma = self.config.invSigma**-1 

716 lineModel = LineProfile(data, weights, line=line) 

717 # Skip any lines that do not cover any data (sometimes happens because of chip gaps) 

718 if lineModel.lineMaskSize == 0: 718 ↛ 719line 718 didn't jump to line 719, because the condition on line 718 was never true

719 continue 

720 

721 fit, chi2, fitFailure = lineModel.fit(dChi2Tol=self.config.dChi2Tolerance) 

722 

723 # Initial estimate should be quite close: fit is deemed unsuccessful if rho or theta 

724 # change more than the allowed bin in rho or theta: 

725 if ((abs(fit.rho - line.rho) > 2 * self.config.rhoBinSize) 725 ↛ 727line 725 didn't jump to line 727, because the condition on line 725 was never true

726 or (abs(fit.theta - line.theta) > 2 * self.config.thetaBinSize)): 

727 fitFailure = True 

728 

729 if fitFailure: 729 ↛ 730line 729 didn't jump to line 730, because the condition on line 729 was never true

730 continue 

731 

732 # Make mask 

733 lineModel.setLineMask(fit) 

734 finalModel = lineModel.makeProfile(fit) 

735 # Take absolute value, as streaks are allowed to be negative 

736 finalModelMax = abs(finalModel).max() 

737 finalLineMask = abs(finalModel) > self.config.footprintThreshold 

738 # Drop this line if the model profile is below the footprint threshold 

739 if not finalLineMask.any(): 739 ↛ 740line 739 didn't jump to line 740, because the condition on line 739 was never true

740 continue 

741 fit.chi2 = chi2 

742 fit.finalModelMax = finalModelMax 

743 lineFits.append(fit) 

744 finalLineMasks.append(finalLineMask) 

745 

746 finalMask = np.array(finalLineMasks).any(axis=0) 

747 

748 return lineFits, finalMask