Coverage for python/lsst/ip/isr/overscan.py: 11%

322 statements  

« prev     ^ index     » next       coverage.py v7.2.5, created at 2023-05-23 08:59 +0000

1# This file is part of ip_isr. 

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__ = ["OverscanCorrectionTaskConfig", "OverscanCorrectionTask"] 

23 

24import numpy as np 

25import lsst.afw.math as afwMath 

26import lsst.afw.image as afwImage 

27import lsst.geom as geom 

28import lsst.pipe.base as pipeBase 

29import lsst.pex.config as pexConfig 

30 

31from .isr import fitOverscanImage 

32from .isrFunctions import makeThresholdMask 

33 

34 

35class OverscanCorrectionTaskConfig(pexConfig.Config): 

36 """Overscan correction options. 

37 """ 

38 fitType = pexConfig.ChoiceField( 

39 dtype=str, 

40 doc="The method for fitting the overscan bias level.", 

41 default='MEDIAN', 

42 allowed={ 

43 "POLY": "Fit ordinary polynomial to the longest axis of the overscan region", 

44 "CHEB": "Fit Chebyshev polynomial to the longest axis of the overscan region", 

45 "LEG": "Fit Legendre polynomial to the longest axis of the overscan region", 

46 "NATURAL_SPLINE": "Fit natural spline to the longest axis of the overscan region", 

47 "CUBIC_SPLINE": "Fit cubic spline to the longest axis of the overscan region", 

48 "AKIMA_SPLINE": "Fit Akima spline to the longest axis of the overscan region", 

49 "MEAN": "Correct using the mean of the overscan region", 

50 "MEANCLIP": "Correct using a clipped mean of the overscan region", 

51 "MEDIAN": "Correct using the median of the overscan region", 

52 "MEDIAN_PER_ROW": "Correct using the median per row of the overscan region", 

53 }, 

54 ) 

55 order = pexConfig.Field( 

56 dtype=int, 

57 doc=("Order of polynomial to fit if overscan fit type is a polynomial, " 

58 "or number of spline knots if overscan fit type is a spline."), 

59 default=1, 

60 ) 

61 numSigmaClip = pexConfig.Field( 

62 dtype=float, 

63 doc="Rejection threshold (sigma) for collapsing overscan before fit", 

64 default=3.0, 

65 ) 

66 maskPlanes = pexConfig.ListField( 

67 dtype=str, 

68 doc="Mask planes to reject when measuring overscan", 

69 default=['BAD', 'SAT'], 

70 ) 

71 overscanIsInt = pexConfig.Field( 

72 dtype=bool, 

73 doc="Treat overscan as an integer image for purposes of fitType=MEDIAN" 

74 " and fitType=MEDIAN_PER_ROW.", 

75 default=True, 

76 ) 

77 

78 doParallelOverscan = pexConfig.Field( 

79 dtype=bool, 

80 doc="Correct using parallel overscan after serial overscan correction?", 

81 default=False, 

82 ) 

83 parallelOverscanMaskThreshold = pexConfig.Field( 

84 dtype=int, 

85 doc="Threshold above which pixels in the parallel overscan are masked as bleeds.", 

86 default=100000, 

87 ) 

88 parallelOverscanMaskGrowSize = pexConfig.Field( 

89 dtype=int, 

90 doc="Masks created from saturated bleeds should be grown by this many " 

91 "pixels during construction of the parallel overscan mask. " 

92 "This value determined from the ITL chip in the LATISS camera", 

93 default=7, 

94 ) 

95 

96 leadingColumnsToSkip = pexConfig.Field( 

97 dtype=int, 

98 doc="Number of leading columns to skip in serial overscan correction.", 

99 default=0, 

100 ) 

101 trailingColumnsToSkip = pexConfig.Field( 

102 dtype=int, 

103 doc="Number of trailing columns to skip in serial overscan correction.", 

104 default=0, 

105 ) 

106 leadingRowsToSkip = pexConfig.Field( 

107 dtype=int, 

108 doc="Number of leading rows to skip in parallel overscan correction.", 

109 default=0, 

110 ) 

111 trailingRowsToSkip = pexConfig.Field( 

112 dtype=int, 

113 doc="Number of trailing rows to skip in parallel overscan correction.", 

114 default=0, 

115 ) 

116 

117 maxDeviation = pexConfig.Field( 117 ↛ exitline 117 didn't jump to the function exit

118 dtype=float, 

119 doc="Maximum deviation from median (in ADU) to mask in overscan correction.", 

120 default=1000.0, check=lambda x: x > 0, 

121 ) 

122 

123 

124class OverscanCorrectionTask(pipeBase.Task): 

125 """Correction task for overscan. 

126 

127 This class contains a number of utilities that are easier to 

128 understand and use when they are not embedded in nested if/else 

129 loops. 

130 

131 Parameters 

132 ---------- 

133 statControl : `lsst.afw.math.StatisticsControl`, optional 

134 Statistics control object. 

135 """ 

136 ConfigClass = OverscanCorrectionTaskConfig 

137 _DefaultName = "overscan" 

138 

139 def __init__(self, statControl=None, **kwargs): 

140 super().__init__(**kwargs) 

141 self.allowDebug = True 

142 

143 if statControl: 

144 self.statControl = statControl 

145 else: 

146 self.statControl = afwMath.StatisticsControl() 

147 self.statControl.setNumSigmaClip(self.config.numSigmaClip) 

148 self.statControl.setAndMask(afwImage.Mask.getPlaneBitMask(self.config.maskPlanes)) 

149 

150 def run(self, exposure, amp, isTransposed=False): 

151 """Measure and remove an overscan from an amplifier image. 

152 

153 Parameters 

154 ---------- 

155 exposure : `lsst.afw.image.Exposure` 

156 Image data that will have the overscan corrections applied. 

157 amp : `lsst.afw.cameraGeom.Amplifier` 

158 Amplifier to use for debugging purposes. 

159 isTransposed : `bool`, optional 

160 Is the image transposed, such that serial and parallel 

161 overscan regions are reversed? Default is False. 

162 

163 Returns 

164 ------- 

165 overscanResults : `lsst.pipe.base.Struct` 

166 Result struct with components: 

167 

168 ``imageFit`` 

169 Value or fit subtracted from the amplifier image data 

170 (scalar or `lsst.afw.image.Image`). 

171 ``overscanFit`` 

172 Value or fit subtracted from the serial overscan image 

173 data (scalar or `lsst.afw.image.Image`). 

174 ``overscanImage`` 

175 Image of the serial overscan region with the serial 

176 overscan correction applied 

177 (`lsst.afw.image.Image`). This quantity is used to 

178 estimate the amplifier read noise empirically. 

179 ``parallelOverscanFit`` 

180 Value or fit subtracted from the parallel overscan 

181 image data (scalar, `lsst.afw.image.Image`, or None). 

182 ``parallelOverscanImage`` 

183 Image of the parallel overscan region with the 

184 parallel overscan correction applied 

185 (`lsst.afw.image.Image` or None). 

186 

187 Raises 

188 ------ 

189 RuntimeError 

190 Raised if an invalid overscan type is set. 

191 """ 

192 # Do Serial overscan first. 

193 serialOverscanBBox = amp.getRawSerialOverscanBBox() 

194 imageBBox = amp.getRawDataBBox() 

195 

196 if self.config.doParallelOverscan: 

197 # We need to extend the serial overscan BBox to the full 

198 # size of the detector. 

199 parallelOverscanBBox = amp.getRawParallelOverscanBBox() 

200 imageBBox = imageBBox.expandedTo(parallelOverscanBBox) 

201 

202 serialOverscanBBox = geom.Box2I(geom.Point2I(serialOverscanBBox.getMinX(), 

203 imageBBox.getMinY()), 

204 geom.Extent2I(serialOverscanBBox.getWidth(), 

205 imageBBox.getHeight())) 

206 serialResults = self.correctOverscan(exposure, amp, 

207 imageBBox, serialOverscanBBox, isTransposed=isTransposed) 

208 overscanMean = serialResults.overscanMean 

209 overscanMedian = serialResults.overscanMedian 

210 overscanSigma = serialResults.overscanSigma 

211 residualMean = serialResults.overscanMeanResidual 

212 residualMedian = serialResults.overscanMedianResidual 

213 residualSigma = serialResults.overscanSigmaResidual 

214 

215 # Do Parallel Overscan 

216 parallelResults = None 

217 if self.config.doParallelOverscan: 

218 # This does not need any extensions, as we'll only 

219 # subtract it from the data region. 

220 parallelOverscanBBox = amp.getRawParallelOverscanBBox() 

221 imageBBox = amp.getRawDataBBox() 

222 

223 maskIm = exposure.getMaskedImage() 

224 maskIm = maskIm.Factory(maskIm, parallelOverscanBBox) 

225 

226 # The serial overscan correction has removed some signal 

227 # from the parallel overscan region, but that is largely a 

228 # constant offset. The collapseArray method now attempts 

229 # to fill fully masked columns with the median of 

230 # neighboring values, with a fallback to the median of the 

231 # correction in all other columns. Filling with neighbor 

232 # values ensures that large variations in the parallel 

233 # overscan do not create new outlier points. The 

234 # MEDIAN_PER_ROW method does this filling as a separate 

235 # operation, using the same method. 

236 parallelResults = self.correctOverscan(exposure, amp, 

237 imageBBox, parallelOverscanBBox, 

238 isTransposed=not isTransposed) 

239 overscanMean = (overscanMean, parallelResults.overscanMean) 

240 overscanMedian = (overscanMedian, parallelResults.overscanMedian) 

241 overscanSigma = (overscanSigma, parallelResults.overscanSigma) 

242 residualMean = (residualMean, parallelResults.overscanMeanResidual) 

243 residualMedian = (residualMedian, parallelResults.overscanMedianResidual) 

244 residualSigma = (residualSigma, parallelResults.overscanSigmaResidual) 

245 

246 parallelOverscanFit = parallelResults.overscanOverscanModel if parallelResults else None 

247 parallelOverscanImage = parallelResults.overscanImage if parallelResults else None 

248 

249 return pipeBase.Struct(imageFit=serialResults.ampOverscanModel, 

250 overscanFit=serialResults.overscanOverscanModel, 

251 overscanImage=serialResults.overscanImage, 

252 

253 parallelOverscanFit=parallelOverscanFit, 

254 parallelOverscanImage=parallelOverscanImage, 

255 overscanMean=overscanMean, 

256 overscanMedian=overscanMedian, 

257 overscanSigma=overscanSigma, 

258 residualMean=residualMean, 

259 residualMedian=residualMedian, 

260 residualSigma=residualSigma) 

261 

262 def correctOverscan(self, exposure, amp, imageBBox, overscanBBox, isTransposed=True): 

263 """Trim the exposure, fit the overscan, subtract the fit, and 

264 calculate statistics. 

265 

266 Parameters 

267 ---------- 

268 exposure : `lsst.afw.image.Exposure` 

269 Exposure containing the data. 

270 amp : `lsst.afw.cameraGeom.Amplifier` 

271 The amplifier that is to be corrected. 

272 imageBBox: `lsst.geom.Box2I` 

273 Bounding box of the image data that will have the overscan 

274 subtracted. If parallel overscan will be performed, that 

275 area is added to the image bounding box during serial 

276 overscan correction. 

277 overscanBBox: `lsst.geom.Box2I` 

278 Bounding box for the overscan data. 

279 isTransposed: `bool` 

280 If true, then the data will be transposed before fitting 

281 the overscan. 

282 

283 Returns 

284 ------- 

285 results : `lsst.pipe.base.Struct` 

286 ``ampOverscanModel`` 

287 Overscan model broadcast to the full image size. 

288 (`lsst.afw.image.Exposure`) 

289 ``overscanOverscanModel`` 

290 Overscan model broadcast to the full overscan image 

291 size. (`lsst.afw.image.Exposure`) 

292 ``overscanImage`` 

293 Overscan image with the overscan fit subtracted. 

294 (`lsst.afw.image.Exposure`) 

295 ``overscanValue`` 

296 Overscan model. (`float` or `np.array`) 

297 ``overscanMean`` 

298 Mean value of the overscan fit. (`float`) 

299 ``overscanMedian`` 

300 Median value of the overscan fit. (`float`) 

301 ``overscanSigma`` 

302 Standard deviation of the overscan fit. (`float`) 

303 ``overscanMeanResidual`` 

304 Mean value of the overscan region after overscan 

305 subtraction. (`float`) 

306 ``overscanMedianResidual`` 

307 Median value of the overscan region after overscan 

308 subtraction. (`float`) 

309 ``overscanSigmaResidual`` 

310 Standard deviation of the overscan region after 

311 overscan subtraction. (`float`) 

312 """ 

313 overscanBox = self.trimOverscan(exposure, amp, overscanBBox, 

314 self.config.leadingColumnsToSkip, 

315 self.config.trailingColumnsToSkip, 

316 transpose=isTransposed) 

317 overscanImage = exposure[overscanBox].getMaskedImage() 

318 overscanArray = overscanImage.image.array 

319 

320 # Mask pixels. 

321 maskVal = overscanImage.mask.getPlaneBitMask(self.config.maskPlanes) 

322 overscanMask = ~((overscanImage.mask.array & maskVal) == 0) 

323 

324 median = np.ma.median(np.ma.masked_where(overscanMask, overscanArray)) 

325 bad = np.where(np.abs(overscanArray - median) > self.config.maxDeviation) 

326 overscanImage.mask.array[bad] = overscanImage.mask.getPlaneBitMask("SAT") 

327 

328 # Do overscan fit. 

329 # CZW: Handle transposed correctly. 

330 overscanResults = self.fitOverscan(overscanImage, isTransposed=isTransposed) 

331 

332 # Correct image region (and possibly parallel-overscan region). 

333 ampImage = exposure[imageBBox] 

334 ampOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue, 

335 ampImage.image.array, 

336 transpose=isTransposed) 

337 ampImage.image.array -= ampOverscanModel 

338 

339 # Correct overscan region (and possibly doubly-overscaned 

340 # region). 

341 overscanImage = exposure[overscanBBox] 

342 # CZW: Transposed? 

343 overscanOverscanModel = self.broadcastFitToImage(overscanResults.overscanValue, 

344 overscanImage.image.array) 

345 self.debugView(overscanImage, overscanResults.overscanValue, amp, isTransposed=isTransposed) 

346 overscanImage.image.array -= overscanOverscanModel 

347 

348 # Find residual fit statistics. 

349 stats = afwMath.makeStatistics(overscanImage.getMaskedImage(), 

350 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP, self.statControl) 

351 residualMean = stats.getValue(afwMath.MEAN) 

352 residualMedian = stats.getValue(afwMath.MEDIAN) 

353 residualSigma = stats.getValue(afwMath.STDEVCLIP) 

354 

355 return pipeBase.Struct(ampOverscanModel=ampOverscanModel, 

356 overscanOverscanModel=overscanOverscanModel, 

357 overscanImage=overscanImage, 

358 overscanValue=overscanResults.overscanValue, 

359 

360 overscanMean=overscanResults.overscanMean, 

361 overscanMedian=overscanResults.overscanMedian, 

362 overscanSigma=overscanResults.overscanSigma, 

363 overscanMeanResidual=residualMean, 

364 overscanMedianResidual=residualMedian, 

365 overscanSigmaResidual=residualSigma 

366 ) 

367 

368 def broadcastFitToImage(self, overscanValue, imageArray, transpose=False): 

369 """Broadcast 0 or 1 dimension fit to appropriate shape. 

370 

371 Parameters 

372 ---------- 

373 overscanValue : `numpy.ndarray`, (Nrows, ) or scalar 

374 Overscan fit to broadcast. 

375 imageArray : `numpy.ndarray`, (Nrows, Ncols) 

376 Image array that we want to match. 

377 transpose : `bool`, optional 

378 Switch order to broadcast along the other axis. 

379 

380 Returns 

381 ------- 

382 overscanModel : `numpy.ndarray`, (Nrows, Ncols) or scalar 

383 Expanded overscan fit. 

384 

385 Raises 

386 ------ 

387 RuntimeError 

388 Raised if no axis has the appropriate dimension. 

389 """ 

390 if isinstance(overscanValue, np.ndarray): 

391 overscanModel = np.zeros_like(imageArray) 

392 

393 if transpose is False: 

394 if imageArray.shape[0] == overscanValue.shape[0]: 

395 overscanModel[:, :] = overscanValue[:, np.newaxis] 

396 elif imageArray.shape[1] == overscanValue.shape[0]: 

397 overscanModel[:, :] = overscanValue[np.newaxis, :] 

398 elif imageArray.shape[0] == overscanValue.shape[1]: 

399 overscanModel[:, :] = overscanValue[np.newaxis, :] 

400 else: 

401 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to " 

402 f"match {imageArray.shape}") 

403 else: 

404 if imageArray.shape[1] == overscanValue.shape[0]: 

405 overscanModel[:, :] = overscanValue[np.newaxis, :] 

406 elif imageArray.shape[0] == overscanValue.shape[0]: 

407 overscanModel[:, :] = overscanValue[:, np.newaxis] 

408 elif imageArray.shape[1] == overscanValue.shape[1]: 

409 overscanModel[:, :] = overscanValue[:, np.newaxis] 

410 else: 

411 raise RuntimeError(f"Could not broadcast {overscanValue.shape} to " 

412 f"match {imageArray.shape}") 

413 else: 

414 overscanModel = overscanValue 

415 

416 return overscanModel 

417 

418 def trimOverscan(self, exposure, amp, bbox, skipLeading, skipTrailing, transpose=False): 

419 """Trim overscan region to remove edges. 

420 

421 Parameters 

422 ---------- 

423 exposure : `lsst.afw.image.Exposure` 

424 Exposure containing data. 

425 amp : `lsst.afw.cameraGeom.Amplifier` 

426 Amplifier containing geometry information. 

427 bbox : `lsst.geom.Box2I` 

428 Bounding box of the overscan region. 

429 skipLeading : `int` 

430 Number of leading (towards data region) rows/columns to skip. 

431 skipTrailing : `int` 

432 Number of trailing (away from data region) rows/columns to skip. 

433 transpose : `bool`, optional 

434 Operate on the transposed array. 

435 

436 Returns 

437 ------- 

438 overscanArray : `numpy.array`, (N, M) 

439 Data array to fit. 

440 overscanMask : `numpy.array`, (N, M) 

441 Data mask. 

442 """ 

443 dx0, dy0, dx1, dy1 = (0, 0, 0, 0) 

444 dataBBox = amp.getRawDataBBox() 

445 if transpose: 

446 if dataBBox.getBeginY() < bbox.getBeginY(): 

447 dy0 += skipLeading 

448 dy1 -= skipTrailing 

449 else: 

450 dy0 += skipTrailing 

451 dy1 -= skipLeading 

452 else: 

453 if dataBBox.getBeginX() < bbox.getBeginX(): 

454 dx0 += skipLeading 

455 dx1 -= skipTrailing 

456 else: 

457 dx0 += skipTrailing 

458 dx1 -= skipLeading 

459 

460 overscanBBox = geom.Box2I(bbox.getBegin() + geom.Extent2I(dx0, dy0), 

461 geom.Extent2I(bbox.getWidth() - dx0 + dx1, 

462 bbox.getHeight() - dy0 + dy1)) 

463 return overscanBBox 

464 

465 def fitOverscan(self, overscanImage, isTransposed=False): 

466 if self.config.fitType in ('MEAN', 'MEANCLIP', 'MEDIAN'): 

467 # Transposition has no effect here. 

468 overscanResult = self.measureConstantOverscan(overscanImage) 

469 overscanValue = overscanResult.overscanValue 

470 overscanMean = overscanValue 

471 overscanMedian = overscanValue 

472 overscanSigma = 0.0 

473 elif self.config.fitType in ('MEDIAN_PER_ROW', 'POLY', 'CHEB', 'LEG', 

474 'NATURAL_SPLINE', 'CUBIC_SPLINE', 'AKIMA_SPLINE'): 

475 # Force transposes as needed 

476 overscanResult = self.measureVectorOverscan(overscanImage, isTransposed) 

477 overscanValue = overscanResult.overscanValue 

478 

479 stats = afwMath.makeStatistics(overscanResult.overscanValue, 

480 afwMath.MEAN | afwMath.MEDIAN | afwMath.STDEVCLIP, 

481 self.statControl) 

482 overscanMean = stats.getValue(afwMath.MEAN) 

483 overscanMedian = stats.getValue(afwMath.MEDIAN) 

484 overscanSigma = stats.getValue(afwMath.STDEVCLIP) 

485 else: 

486 raise ValueError('%s : %s an invalid overscan type' % 

487 ("overscanCorrection", self.config.fitType)) 

488 

489 return pipeBase.Struct(overscanValue=overscanValue, 

490 overscanMean=overscanMean, 

491 overscanMedian=overscanMedian, 

492 overscanSigma=overscanSigma, 

493 ) 

494 

495 @staticmethod 

496 def integerConvert(image): 

497 """Return an integer version of the input image. 

498 

499 Parameters 

500 ---------- 

501 image : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage` 

502 Image to convert to integers. 

503 

504 Returns 

505 ------- 

506 outI : `numpy.ndarray`, `lsst.afw.image.Image` or `MaskedImage` 

507 The integer converted image. 

508 

509 Raises 

510 ------ 

511 RuntimeError 

512 Raised if the input image could not be converted. 

513 """ 

514 if hasattr(image, "image"): 

515 # Is a maskedImage: 

516 imageI = image.image.convertI() 

517 outI = afwImage.MaskedImageI(imageI, image.mask, image.variance) 

518 elif hasattr(image, "convertI"): 

519 # Is an Image: 

520 outI = image.convertI() 

521 elif hasattr(image, "astype"): 

522 # Is a numpy array: 

523 outI = image.astype(int) 

524 else: 

525 raise RuntimeError("Could not convert this to integers: %s %s %s", 

526 image, type(image), dir(image)) 

527 return outI 

528 

529 def maskParallelOverscan(self, exposure, detector): 

530 """Mask the union of high values on all amplifiers in the parallel 

531 overscan. 

532 

533 This operates on the image in-place. 

534 

535 Parameters 

536 ---------- 

537 exposure : `lsst.afw.image.Exposure` 

538 An untrimmed raw exposure. 

539 detector : `lsst.afw.cameraGeom.Detector` 

540 The detetor to use for amplifier geometry. 

541 """ 

542 parallelMask = None 

543 

544 for amp in detector: 

545 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(), 

546 amp.getRawParallelOverscanBBox(), 

547 afwImage.PARENT) 

548 makeThresholdMask( 

549 maskedImage=dataView, 

550 threshold=self.config.parallelOverscanMaskThreshold, 

551 growFootprints=self.config.parallelOverscanMaskGrowSize, 

552 maskName="BAD" 

553 ) 

554 if parallelMask is None: 

555 parallelMask = dataView.mask.array 

556 else: 

557 parallelMask |= dataView.mask.array 

558 for amp in detector: 

559 dataView = afwImage.MaskedImageF(exposure.getMaskedImage(), 

560 amp.getRawParallelOverscanBBox(), 

561 afwImage.PARENT) 

562 dataView.mask.array |= parallelMask 

563 

564 # Constant methods 

565 def measureConstantOverscan(self, image): 

566 """Measure a constant overscan value. 

567 

568 Parameters 

569 ---------- 

570 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage` 

571 Image data to measure the overscan from. 

572 

573 Returns 

574 ------- 

575 results : `lsst.pipe.base.Struct` 

576 Overscan result with entries: 

577 - ``overscanValue``: Overscan value to subtract (`float`) 

578 - ``isTransposed``: Orientation of the overscan (`bool`) 

579 """ 

580 fitType = afwMath.stringToStatisticsProperty(self.config.fitType) 

581 overscanValue = afwMath.makeStatistics(image, fitType, self.statControl).getValue() 

582 

583 return pipeBase.Struct(overscanValue=overscanValue, 

584 isTransposed=False) 

585 

586 # Vector correction utilities 

587 def getImageArray(self, image): 

588 """Extract the numpy array from the input image. 

589 

590 Parameters 

591 ---------- 

592 image : `lsst.afw.image.Image` or `lsst.afw.image.MaskedImage` 

593 Image data to pull array from. 

594 

595 calcImage : `numpy.ndarray` 

596 Image data array for numpy operating. 

597 """ 

598 if hasattr(image, "getImage"): 

599 calcImage = image.getImage().getArray() 

600 calcImage = np.ma.masked_where(image.getMask().getArray() & self.statControl.getAndMask(), 

601 calcImage) 

602 else: 

603 calcImage = image.getArray() 

604 return calcImage 

605 

606 def maskOutliers(self, imageArray): 

607 """Mask outliers in a row of overscan data from a robust sigma 

608 clipping procedure. 

609 

610 Parameters 

611 ---------- 

612 imageArray : `numpy.ndarray` 

613 Image to filter along numpy axis=1. 

614 

615 Returns 

616 ------- 

617 maskedArray : `numpy.ma.masked_array` 

618 Masked image marking outliers. 

619 """ 

620 lq, median, uq = np.percentile(np.ma.getdata(imageArray), 

621 [25.0, 50.0, 75.0], axis=1) 

622 axisMedians = median 

623 axisStdev = 0.74*(uq - lq) # robust stdev 

624 

625 # Replace pixels that have excessively large stdev values 

626 # with the median of stdev values. A large stdev likely 

627 # indicates a bleed is spilling into the overscan. 

628 axisStdev = np.where(axisStdev > 2.0 * np.median(axisStdev), 

629 np.median(axisStdev), axisStdev) 

630 

631 # Mask pixels that are N-sigma away from their array medians. 

632 diff = np.abs(imageArray - axisMedians[:, np.newaxis]) 

633 masked = np.ma.masked_where(diff > self.statControl.getNumSigmaClip() 

634 * axisStdev[:, np.newaxis], imageArray) 

635 

636 return masked 

637 

638 def fillMaskedPixels(self, overscanVector): 

639 """Fill masked/NaN pixels in the overscan. 

640 

641 Parameters 

642 ---------- 

643 overscanVector : `np.array` or `np.ma.masked_array` 

644 Overscan vector to fill. 

645 

646 Returns 

647 ------- 

648 overscanVector : `np.ma.masked_array` 

649 Filled vector. 

650 

651 Notes 

652 ----- 

653 Each maskSlice is a section of overscan with contiguous masks. 

654 Ideally this adds 5 pixels from the left and right of that 

655 mask slice, and takes the median of those values to fill the 

656 slice. If this isn't possible, the median of all non-masked 

657 values is used. The mask is removed for the pixels filled. 

658 """ 

659 workingCopy = overscanVector 

660 if not isinstance(overscanVector, np.ma.MaskedArray): 

661 workingCopy = np.ma.masked_array(overscanVector, 

662 mask=~np.isfinite(overscanVector)) 

663 

664 defaultValue = np.median(workingCopy.data[~workingCopy.mask]) 

665 for maskSlice in np.ma.clump_masked(workingCopy): 

666 neighborhood = [] 

667 if maskSlice.start > 5: 

668 neighborhood.extend(workingCopy[maskSlice.start - 5:maskSlice.start].data) 

669 if maskSlice.stop < workingCopy.size - 5: 

670 neighborhood.extend(workingCopy[maskSlice.stop:maskSlice.stop+5].data) 

671 if len(neighborhood) > 0: 

672 workingCopy.data[maskSlice] = np.nanmedian(neighborhood) 

673 workingCopy.mask[maskSlice] = False 

674 else: 

675 workingCopy.data[maskSlice] = defaultValue 

676 workingCopy.mask[maskSlice] = False 

677 return workingCopy 

678 

679 def collapseArray(self, maskedArray, fillMasked=True): 

680 """Collapse overscan array (and mask) to a 1-D vector of values. 

681 

682 Parameters 

683 ---------- 

684 maskedArray : `numpy.ma.masked_array` 

685 Masked array of input overscan data. 

686 fillMasked : `bool`, optional 

687 If true, fill any pixels that are masked with a median of 

688 neighbors. 

689 

690 Returns 

691 ------- 

692 collapsed : `numpy.ma.masked_array` 

693 Single dimensional overscan data, combined with the mean. 

694 

695 """ 

696 collapsed = np.mean(maskedArray, axis=1) 

697 if collapsed.mask.sum() > 0 and fillMasked: 

698 collapsed = self.fillMaskedPixels(collapsed) 

699 

700 return collapsed 

701 

702 def collapseArrayMedian(self, maskedArray): 

703 """Collapse overscan array (and mask) to a 1-D vector of using the 

704 correct integer median of row-values. 

705 

706 Parameters 

707 ---------- 

708 maskedArray : `numpy.ma.masked_array` 

709 Masked array of input overscan data. 

710 

711 Returns 

712 ------- 

713 collapsed : `numpy.ma.masked_array` 

714 Single dimensional overscan data, combined with the afwMath median. 

715 """ 

716 integerMI = self.integerConvert(maskedArray) 

717 

718 collapsed = [] 

719 fitType = afwMath.stringToStatisticsProperty('MEDIAN') 

720 for row in integerMI: 

721 newRow = row.compressed() 

722 if len(newRow) > 0: 

723 rowMedian = afwMath.makeStatistics(newRow, fitType, self.statControl).getValue() 

724 else: 

725 rowMedian = np.nan 

726 collapsed.append(rowMedian) 

727 

728 return np.array(collapsed) 

729 

730 def splineFit(self, indices, collapsed, numBins): 

731 """Wrapper function to match spline fit API to polynomial fit API. 

732 

733 Parameters 

734 ---------- 

735 indices : `numpy.ndarray` 

736 Locations to evaluate the spline. 

737 collapsed : `numpy.ndarray` 

738 Collapsed overscan values corresponding to the spline 

739 evaluation points. 

740 numBins : `int` 

741 Number of bins to use in constructing the spline. 

742 

743 Returns 

744 ------- 

745 interp : `lsst.afw.math.Interpolate` 

746 Interpolation object for later evaluation. 

747 """ 

748 if not np.ma.is_masked(collapsed): 

749 collapsed.mask = np.array(len(collapsed)*[np.ma.nomask]) 

750 

751 numPerBin, binEdges = np.histogram(indices, bins=numBins, 

752 weights=1 - collapsed.mask.astype(int)) 

753 with np.errstate(invalid="ignore"): 

754 values = np.histogram(indices, bins=numBins, 

755 weights=collapsed.data*~collapsed.mask)[0]/numPerBin 

756 binCenters = np.histogram(indices, bins=numBins, 

757 weights=indices*~collapsed.mask)[0]/numPerBin 

758 

759 if len(binCenters[numPerBin > 0]) < 5: 

760 self.log.warn("Cannot do spline fitting for overscan: %s valid points.", 

761 len(binCenters[numPerBin > 0])) 

762 # Return a scalar value if we have one, otherwise 

763 # return zero. This amplifier is hopefully already 

764 # masked. 

765 if len(values[numPerBin > 0]) != 0: 

766 return float(values[numPerBin > 0][0]) 

767 else: 

768 return 0.0 

769 

770 interp = afwMath.makeInterpolate(binCenters.astype(float)[numPerBin > 0], 

771 values.astype(float)[numPerBin > 0], 

772 afwMath.stringToInterpStyle(self.config.fitType)) 

773 return interp 

774 

775 @staticmethod 

776 def splineEval(indices, interp): 

777 """Wrapper function to match spline evaluation API to polynomial fit 

778 API. 

779 

780 Parameters 

781 ---------- 

782 indices : `numpy.ndarray` 

783 Locations to evaluate the spline. 

784 interp : `lsst.afw.math.interpolate` 

785 Interpolation object to use. 

786 

787 Returns 

788 ------- 

789 values : `numpy.ndarray` 

790 Evaluated spline values at each index. 

791 """ 

792 

793 return interp.interpolate(indices.astype(float)) 

794 

795 @staticmethod 

796 def maskExtrapolated(collapsed): 

797 """Create mask if edges are extrapolated. 

798 

799 Parameters 

800 ---------- 

801 collapsed : `numpy.ma.masked_array` 

802 Masked array to check the edges of. 

803 

804 Returns 

805 ------- 

806 maskArray : `numpy.ndarray` 

807 Boolean numpy array of pixels to mask. 

808 """ 

809 maskArray = np.full_like(collapsed, False, dtype=bool) 

810 if np.ma.is_masked(collapsed): 

811 num = len(collapsed) 

812 for low in range(num): 

813 if not collapsed.mask[low]: 

814 break 

815 if low > 0: 

816 maskArray[:low] = True 

817 for high in range(1, num): 

818 if not collapsed.mask[-high]: 

819 break 

820 if high > 1: 

821 maskArray[-high:] = True 

822 return maskArray 

823 

824 def measureVectorOverscan(self, image, isTransposed=False): 

825 """Calculate the 1-d vector overscan from the input overscan image. 

826 

827 Parameters 

828 ---------- 

829 image : `lsst.afw.image.MaskedImage` 

830 Image containing the overscan data. 

831 isTransposed : `bool` 

832 If true, the image has been transposed. 

833 

834 Returns 

835 ------- 

836 results : `lsst.pipe.base.Struct` 

837 Overscan result with entries: 

838 

839 ``overscanValue`` 

840 Overscan value to subtract (`float`) 

841 ``maskArray`` 

842 List of rows that should be masked as ``SUSPECT`` when the 

843 overscan solution is applied. (`list` [ `bool` ]) 

844 ``isTransposed`` 

845 Indicates if the overscan data was transposed during 

846 calcuation, noting along which axis the overscan should be 

847 subtracted. (`bool`) 

848 """ 

849 calcImage = self.getImageArray(image) 

850 

851 # operate on numpy-arrays from here 

852 if isTransposed: 

853 calcImage = np.transpose(calcImage) 

854 masked = self.maskOutliers(calcImage) 

855 

856 if self.config.fitType == 'MEDIAN_PER_ROW': 

857 mi = afwImage.MaskedImageI(image.getBBox()) 

858 masked = masked.astype(int) 

859 if isTransposed: 

860 masked = masked.transpose() 

861 

862 mi.image.array[:, :] = masked.data[:, :] 

863 if bool(masked.mask.shape): 

864 mi.mask.array[:, :] = masked.mask[:, :] 

865 

866 overscanVector = fitOverscanImage(mi, self.config.maskPlanes, isTransposed) 

867 overscanVector = self.fillMaskedPixels(overscanVector) 

868 maskArray = self.maskExtrapolated(overscanVector) 

869 else: 

870 collapsed = self.collapseArray(masked) 

871 

872 num = len(collapsed) 

873 indices = 2.0*np.arange(num)/float(num) - 1.0 

874 

875 poly = np.polynomial 

876 fitter, evaler = { 

877 'POLY': (poly.polynomial.polyfit, poly.polynomial.polyval), 

878 'CHEB': (poly.chebyshev.chebfit, poly.chebyshev.chebval), 

879 'LEG': (poly.legendre.legfit, poly.legendre.legval), 

880 'NATURAL_SPLINE': (self.splineFit, self.splineEval), 

881 'CUBIC_SPLINE': (self.splineFit, self.splineEval), 

882 'AKIMA_SPLINE': (self.splineFit, self.splineEval) 

883 }[self.config.fitType] 

884 

885 # These are the polynomial coefficients, or an 

886 # interpolation object. 

887 coeffs = fitter(indices, collapsed, self.config.order) 

888 

889 if isinstance(coeffs, float): 

890 self.log.warn("Using fallback value %f due to fitter failure. Amplifier will be masked.", 

891 coeffs) 

892 overscanVector = np.full_like(indices, coeffs) 

893 maskArray = np.full_like(collapsed, True, dtype=bool) 

894 else: 

895 # Otherwise we can just use things as normal. 

896 overscanVector = evaler(indices, coeffs) 

897 maskArray = self.maskExtrapolated(collapsed) 

898 

899 return pipeBase.Struct(overscanValue=np.array(overscanVector), 

900 maskArray=maskArray, 

901 isTransposed=isTransposed) 

902 

903 def debugView(self, image, model, amp=None, isTransposed=True): 

904 """Debug display for the final overscan solution. 

905 

906 Parameters 

907 ---------- 

908 image : `lsst.afw.image.Image` 

909 Input image the overscan solution was determined from. 

910 model : `numpy.ndarray` or `float` 

911 Overscan model determined for the image. 

912 amp : `lsst.afw.cameraGeom.Amplifier`, optional 

913 Amplifier to extract diagnostic information. 

914 isTransposed : `bool`, optional 

915 Does the data need to be transposed before display? 

916 """ 

917 import lsstDebug 

918 if not lsstDebug.Info(__name__).display: 

919 return 

920 if not self.allowDebug: 

921 return 

922 

923 calcImage = self.getImageArray(image) 

924 # CZW: Check that this is ok 

925 if isTransposed: 

926 calcImage = np.transpose(calcImage) 

927 masked = self.maskOutliers(calcImage) 

928 collapsed = self.collapseArray(masked, fillMasked=False) 

929 

930 num = len(collapsed) 

931 indices = 2.0 * np.arange(num)/float(num) - 1.0 

932 indices = np.arange(num) 

933 

934 if np.ma.is_masked(collapsed): 

935 collapsedMask = collapsed.mask 

936 else: 

937 collapsedMask = np.array(num*[np.ma.nomask]) 

938 

939 import matplotlib.pyplot as plot 

940 figure = plot.figure(1) 

941 figure.clear() 

942 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8)) 

943 axes.plot(indices[~collapsedMask], collapsed[~collapsedMask], 'k+') 

944 if collapsedMask.sum() > 0: 

945 axes.plot(indices[collapsedMask], collapsed.data[collapsedMask], 'b+') 

946 if isinstance(model, np.ndarray): 

947 plotModel = model 

948 else: 

949 plotModel = np.zeros_like(indices) 

950 plotModel += model 

951 

952 axes.plot(indices, plotModel, 'r-') 

953 plot.xlabel("position along overscan region") 

954 plot.ylabel("pixel value/fit value") 

955 if amp: 

956 plot.title(f"{amp.getName()} DataX: " 

957 f"[{amp.getRawDataBBox().getBeginX()}:{amp.getRawBBox().getEndX()}]" 

958 f"OscanX: [{amp.getRawHorizontalOverscanBBox().getBeginX()}:" 

959 f"{amp.getRawHorizontalOverscanBBox().getEndX()}] {self.config.fitType}") 

960 else: 

961 plot.title("No amp supplied.") 

962 figure.show() 

963 prompt = "Press Enter or c to continue [chp]..." 

964 while True: 

965 ans = input(prompt).lower() 

966 if ans in ("", " ", "c",): 

967 break 

968 elif ans in ("p", ): 

969 import pdb 

970 pdb.set_trace() 

971 elif ans in ('x', ): 

972 self.allowDebug = False 

973 break 

974 elif ans in ("h", ): 

975 print("[h]elp [c]ontinue [p]db e[x]itDebug") 

976 plot.close()