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# 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 

22import math 

23import numpy 

24 

25import lsst.geom 

26import lsst.afw.image as afwImage 

27import lsst.afw.math as afwMath 

28import lsst.pex.config as pexConfig 

29import lsst.pipe.base as pipeBase 

30import lsst.pipe.base.connectionTypes as cT 

31 

32from contextlib import contextmanager 

33from lsstDebug import getDebugFrame 

34 

35from lsst.afw.cameraGeom import (PIXELS, FOCAL_PLANE, NullLinearityType, 

36 ReadoutCorner) 

37from lsst.afw.display import getDisplay 

38from lsst.afw.geom import Polygon 

39from lsst.daf.persistence import ButlerDataRef 

40from lsst.daf.persistence.butler import NoResults 

41from lsst.meas.algorithms.detection import SourceDetectionTask 

42 

43from . import isrFunctions 

44from . import isrQa 

45from . import linearize 

46from .defects import Defects 

47 

48from .assembleCcdTask import AssembleCcdTask 

49from .crosstalk import CrosstalkTask, CrosstalkCalib 

50from .fringe import FringeTask 

51from .isr import maskNans 

52from .masking import MaskingTask 

53from .overscan import OverscanCorrectionTask 

54from .straylight import StrayLightTask 

55from .vignette import VignetteTask 

56from lsst.daf.butler import DimensionGraph 

57 

58 

59__all__ = ["IsrTask", "IsrTaskConfig", "RunIsrTask", "RunIsrConfig"] 

60 

61 

62def crosstalkSourceLookup(datasetType, registry, quantumDataId, collections): 

63 """Lookup function to identify crosstalkSource entries. 

64 

65 This should return an empty list under most circumstances. Only 

66 when inter-chip crosstalk has been identified should this be 

67 populated. 

68 

69 This will be unused until DM-25348 resolves the quantum graph 

70 generation issue. 

71 

72 Parameters 

73 ---------- 

74 datasetType : `str` 

75 Dataset to lookup. 

76 registry : `lsst.daf.butler.Registry` 

77 Butler registry to query. 

78 quantumDataId : `lsst.daf.butler.ExpandedDataCoordinate` 

79 Data id to transform to identify crosstalkSources. The 

80 ``detector`` entry will be stripped. 

81 collections : `lsst.daf.butler.CollectionSearch` 

82 Collections to search through. 

83 

84 Returns 

85 ------- 

86 results : `list` [`lsst.daf.butler.DatasetRef`] 

87 List of datasets that match the query that will be used as 

88 crosstalkSources. 

89 """ 

90 newDataId = quantumDataId.subset(DimensionGraph(registry.dimensions, names=["instrument", "exposure"])) 

91 results = list(registry.queryDatasets(datasetType, 

92 collections=collections, 

93 dataId=newDataId, 

94 findFirst=True, 

95 ).expanded()) 

96 return results 

97 

98 

99class IsrTaskConnections(pipeBase.PipelineTaskConnections, 

100 dimensions={"instrument", "exposure", "detector"}, 

101 defaultTemplates={}): 

102 ccdExposure = cT.Input( 

103 name="raw", 

104 doc="Input exposure to process.", 

105 storageClass="Exposure", 

106 dimensions=["instrument", "exposure", "detector"], 

107 ) 

108 camera = cT.PrerequisiteInput( 

109 name="camera", 

110 storageClass="Camera", 

111 doc="Input camera to construct complete exposures.", 

112 dimensions=["instrument"], 

113 isCalibration=True, 

114 ) 

115 

116 crosstalk = cT.PrerequisiteInput( 

117 name="crosstalk", 

118 doc="Input crosstalk object", 

119 storageClass="CrosstalkCalib", 

120 dimensions=["instrument", "detector"], 

121 isCalibration=True, 

122 ) 

123 # TODO: DM-25348. This does not work yet to correctly load 

124 # possible crosstalk sources. 

125 crosstalkSources = cT.PrerequisiteInput( 

126 name="isrOverscanCorrected", 

127 doc="Overscan corrected input images.", 

128 storageClass="Exposure", 

129 dimensions=["instrument", "exposure", "detector"], 

130 deferLoad=True, 

131 multiple=True, 

132 lookupFunction=crosstalkSourceLookup, 

133 ) 

134 bias = cT.PrerequisiteInput( 

135 name="bias", 

136 doc="Input bias calibration.", 

137 storageClass="ExposureF", 

138 dimensions=["instrument", "detector"], 

139 isCalibration=True, 

140 ) 

141 dark = cT.PrerequisiteInput( 

142 name='dark', 

143 doc="Input dark calibration.", 

144 storageClass="ExposureF", 

145 dimensions=["instrument", "detector"], 

146 isCalibration=True, 

147 ) 

148 flat = cT.PrerequisiteInput( 

149 name="flat", 

150 doc="Input flat calibration.", 

151 storageClass="ExposureF", 

152 dimensions=["instrument", "physical_filter", "detector"], 

153 isCalibration=True, 

154 ) 

155 fringes = cT.PrerequisiteInput( 

156 name="fringe", 

157 doc="Input fringe calibration.", 

158 storageClass="ExposureF", 

159 dimensions=["instrument", "physical_filter", "detector"], 

160 isCalibration=True, 

161 ) 

162 strayLightData = cT.PrerequisiteInput( 

163 name='yBackground', 

164 doc="Input stray light calibration.", 

165 storageClass="StrayLightData", 

166 dimensions=["instrument", "physical_filter", "detector"], 

167 isCalibration=True, 

168 ) 

169 bfKernel = cT.PrerequisiteInput( 

170 name='bfKernel', 

171 doc="Input brighter-fatter kernel.", 

172 storageClass="NumpyArray", 

173 dimensions=["instrument"], 

174 isCalibration=True, 

175 ) 

176 newBFKernel = cT.PrerequisiteInput( 

177 name='brighterFatterKernel', 

178 doc="Newer complete kernel + gain solutions.", 

179 storageClass="BrighterFatterKernel", 

180 dimensions=["instrument", "detector"], 

181 isCalibration=True, 

182 ) 

183 defects = cT.PrerequisiteInput( 

184 name='defects', 

185 doc="Input defect tables.", 

186 storageClass="Defects", 

187 dimensions=["instrument", "detector"], 

188 isCalibration=True, 

189 ) 

190 linearizer = cT.PrerequisiteInput( 

191 name='linearizer', 

192 storageClass="Linearizer", 

193 doc="Linearity correction calibration.", 

194 dimensions=["instrument", "detector"], 

195 isCalibration=True, 

196 ) 

197 opticsTransmission = cT.PrerequisiteInput( 

198 name="transmission_optics", 

199 storageClass="TransmissionCurve", 

200 doc="Transmission curve due to the optics.", 

201 dimensions=["instrument"], 

202 isCalibration=True, 

203 ) 

204 filterTransmission = cT.PrerequisiteInput( 

205 name="transmission_filter", 

206 storageClass="TransmissionCurve", 

207 doc="Transmission curve due to the filter.", 

208 dimensions=["instrument", "physical_filter"], 

209 isCalibration=True, 

210 ) 

211 sensorTransmission = cT.PrerequisiteInput( 

212 name="transmission_sensor", 

213 storageClass="TransmissionCurve", 

214 doc="Transmission curve due to the sensor.", 

215 dimensions=["instrument", "detector"], 

216 isCalibration=True, 

217 ) 

218 atmosphereTransmission = cT.PrerequisiteInput( 

219 name="transmission_atmosphere", 

220 storageClass="TransmissionCurve", 

221 doc="Transmission curve due to the atmosphere.", 

222 dimensions=["instrument"], 

223 isCalibration=True, 

224 ) 

225 illumMaskedImage = cT.PrerequisiteInput( 

226 name="illum", 

227 doc="Input illumination correction.", 

228 storageClass="MaskedImageF", 

229 dimensions=["instrument", "physical_filter", "detector"], 

230 isCalibration=True, 

231 ) 

232 

233 outputExposure = cT.Output( 

234 name='postISRCCD', 

235 doc="Output ISR processed exposure.", 

236 storageClass="Exposure", 

237 dimensions=["instrument", "exposure", "detector"], 

238 ) 

239 preInterpExposure = cT.Output( 

240 name='preInterpISRCCD', 

241 doc="Output ISR processed exposure, with pixels left uninterpolated.", 

242 storageClass="ExposureF", 

243 dimensions=["instrument", "exposure", "detector"], 

244 ) 

245 outputOssThumbnail = cT.Output( 

246 name="OssThumb", 

247 doc="Output Overscan-subtracted thumbnail image.", 

248 storageClass="Thumbnail", 

249 dimensions=["instrument", "exposure", "detector"], 

250 ) 

251 outputFlattenedThumbnail = cT.Output( 

252 name="FlattenedThumb", 

253 doc="Output flat-corrected thumbnail image.", 

254 storageClass="Thumbnail", 

255 dimensions=["instrument", "exposure", "detector"], 

256 ) 

257 

258 def __init__(self, *, config=None): 

259 super().__init__(config=config) 

260 

261 if config.doBias is not True: 

262 self.prerequisiteInputs.discard("bias") 

263 if config.doLinearize is not True: 

264 self.prerequisiteInputs.discard("linearizer") 

265 if config.doCrosstalk is not True: 

266 self.inputs.discard("crosstalkSources") 

267 self.prerequisiteInputs.discard("crosstalk") 

268 if config.doBrighterFatter is not True: 

269 self.prerequisiteInputs.discard("bfKernel") 

270 self.prerequisiteInputs.discard("newBFKernel") 

271 if config.doDefect is not True: 

272 self.prerequisiteInputs.discard("defects") 

273 if config.doDark is not True: 

274 self.prerequisiteInputs.discard("dark") 

275 if config.doFlat is not True: 

276 self.prerequisiteInputs.discard("flat") 

277 if config.doAttachTransmissionCurve is not True: 

278 self.prerequisiteInputs.discard("opticsTransmission") 

279 self.prerequisiteInputs.discard("filterTransmission") 

280 self.prerequisiteInputs.discard("sensorTransmission") 

281 self.prerequisiteInputs.discard("atmosphereTransmission") 

282 if config.doUseOpticsTransmission is not True: 

283 self.prerequisiteInputs.discard("opticsTransmission") 

284 if config.doUseFilterTransmission is not True: 

285 self.prerequisiteInputs.discard("filterTransmission") 

286 if config.doUseSensorTransmission is not True: 

287 self.prerequisiteInputs.discard("sensorTransmission") 

288 if config.doUseAtmosphereTransmission is not True: 

289 self.prerequisiteInputs.discard("atmosphereTransmission") 

290 if config.doIlluminationCorrection is not True: 

291 self.prerequisiteInputs.discard("illumMaskedImage") 

292 

293 if config.doWrite is not True: 

294 self.outputs.discard("outputExposure") 

295 self.outputs.discard("preInterpExposure") 

296 self.outputs.discard("outputFlattenedThumbnail") 

297 self.outputs.discard("outputOssThumbnail") 

298 if config.doSaveInterpPixels is not True: 

299 self.outputs.discard("preInterpExposure") 

300 if config.qa.doThumbnailOss is not True: 

301 self.outputs.discard("outputOssThumbnail") 

302 if config.qa.doThumbnailFlattened is not True: 

303 self.outputs.discard("outputFlattenedThumbnail") 

304 

305 

306class IsrTaskConfig(pipeBase.PipelineTaskConfig, 

307 pipelineConnections=IsrTaskConnections): 

308 """Configuration parameters for IsrTask. 

309 

310 Items are grouped in the order in which they are executed by the task. 

311 """ 

312 datasetType = pexConfig.Field( 

313 dtype=str, 

314 doc="Dataset type for input data; users will typically leave this alone, " 

315 "but camera-specific ISR tasks will override it", 

316 default="raw", 

317 ) 

318 

319 fallbackFilterName = pexConfig.Field( 

320 dtype=str, 

321 doc="Fallback default filter name for calibrations.", 

322 optional=True 

323 ) 

324 useFallbackDate = pexConfig.Field( 

325 dtype=bool, 

326 doc="Pass observation date when using fallback filter.", 

327 default=False, 

328 ) 

329 expectWcs = pexConfig.Field( 

330 dtype=bool, 

331 default=True, 

332 doc="Expect input science images to have a WCS (set False for e.g. spectrographs)." 

333 ) 

334 fwhm = pexConfig.Field( 

335 dtype=float, 

336 doc="FWHM of PSF in arcseconds.", 

337 default=1.0, 

338 ) 

339 qa = pexConfig.ConfigField( 

340 dtype=isrQa.IsrQaConfig, 

341 doc="QA related configuration options.", 

342 ) 

343 

344 # Image conversion configuration 

345 doConvertIntToFloat = pexConfig.Field( 

346 dtype=bool, 

347 doc="Convert integer raw images to floating point values?", 

348 default=True, 

349 ) 

350 

351 # Saturated pixel handling. 

352 doSaturation = pexConfig.Field( 

353 dtype=bool, 

354 doc="Mask saturated pixels? NB: this is totally independent of the" 

355 " interpolation option - this is ONLY setting the bits in the mask." 

356 " To have them interpolated make sure doSaturationInterpolation=True", 

357 default=True, 

358 ) 

359 saturatedMaskName = pexConfig.Field( 

360 dtype=str, 

361 doc="Name of mask plane to use in saturation detection and interpolation", 

362 default="SAT", 

363 ) 

364 saturation = pexConfig.Field( 

365 dtype=float, 

366 doc="The saturation level to use if no Detector is present in the Exposure (ignored if NaN)", 

367 default=float("NaN"), 

368 ) 

369 growSaturationFootprintSize = pexConfig.Field( 

370 dtype=int, 

371 doc="Number of pixels by which to grow the saturation footprints", 

372 default=1, 

373 ) 

374 

375 # Suspect pixel handling. 

376 doSuspect = pexConfig.Field( 

377 dtype=bool, 

378 doc="Mask suspect pixels?", 

379 default=False, 

380 ) 

381 suspectMaskName = pexConfig.Field( 

382 dtype=str, 

383 doc="Name of mask plane to use for suspect pixels", 

384 default="SUSPECT", 

385 ) 

386 numEdgeSuspect = pexConfig.Field( 

387 dtype=int, 

388 doc="Number of edge pixels to be flagged as untrustworthy.", 

389 default=0, 

390 ) 

391 edgeMaskLevel = pexConfig.ChoiceField( 

392 dtype=str, 

393 doc="Mask edge pixels in which coordinate frame: DETECTOR or AMP?", 

394 default="DETECTOR", 

395 allowed={ 

396 'DETECTOR': 'Mask only the edges of the full detector.', 

397 'AMP': 'Mask edges of each amplifier.', 

398 }, 

399 ) 

400 

401 # Initial masking options. 

402 doSetBadRegions = pexConfig.Field( 

403 dtype=bool, 

404 doc="Should we set the level of all BAD patches of the chip to the chip's average value?", 

405 default=True, 

406 ) 

407 badStatistic = pexConfig.ChoiceField( 

408 dtype=str, 

409 doc="How to estimate the average value for BAD regions.", 

410 default='MEANCLIP', 

411 allowed={ 

412 "MEANCLIP": "Correct using the (clipped) mean of good data", 

413 "MEDIAN": "Correct using the median of the good data", 

414 }, 

415 ) 

416 

417 # Overscan subtraction configuration. 

418 doOverscan = pexConfig.Field( 

419 dtype=bool, 

420 doc="Do overscan subtraction?", 

421 default=True, 

422 ) 

423 overscan = pexConfig.ConfigurableField( 

424 target=OverscanCorrectionTask, 

425 doc="Overscan subtraction task for image segments.", 

426 ) 

427 

428 overscanFitType = pexConfig.ChoiceField( 

429 dtype=str, 

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

431 default='MEDIAN', 

432 allowed={ 

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

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

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

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

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

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

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

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

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

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

443 }, 

444 deprecated=("Please configure overscan via the OverscanCorrectionConfig interface." 

445 " This option will no longer be used, and will be removed after v20.") 

446 ) 

447 overscanOrder = pexConfig.Field( 

448 dtype=int, 

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

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

451 default=1, 

452 deprecated=("Please configure overscan via the OverscanCorrectionConfig interface." 

453 " This option will no longer be used, and will be removed after v20.") 

454 ) 

455 overscanNumSigmaClip = pexConfig.Field( 

456 dtype=float, 

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

458 default=3.0, 

459 deprecated=("Please configure overscan via the OverscanCorrectionConfig interface." 

460 " This option will no longer be used, and will be removed after v20.") 

461 ) 

462 overscanIsInt = pexConfig.Field( 

463 dtype=bool, 

464 doc="Treat overscan as an integer image for purposes of overscan.FitType=MEDIAN" 

465 " and overscan.FitType=MEDIAN_PER_ROW.", 

466 default=True, 

467 deprecated=("Please configure overscan via the OverscanCorrectionConfig interface." 

468 " This option will no longer be used, and will be removed after v20.") 

469 ) 

470 # These options do not get deprecated, as they define how we slice up the image data. 

471 overscanNumLeadingColumnsToSkip = pexConfig.Field( 

472 dtype=int, 

473 doc="Number of columns to skip in overscan, i.e. those closest to amplifier", 

474 default=0, 

475 ) 

476 overscanNumTrailingColumnsToSkip = pexConfig.Field( 

477 dtype=int, 

478 doc="Number of columns to skip in overscan, i.e. those farthest from amplifier", 

479 default=0, 

480 ) 

481 overscanMaxDev = pexConfig.Field( 481 ↛ exitline 481 didn't jump to the function exit

482 dtype=float, 

483 doc="Maximum deviation from the median for overscan", 

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

485 ) 

486 overscanBiasJump = pexConfig.Field( 

487 dtype=bool, 

488 doc="Fit the overscan in a piecewise-fashion to correct for bias jumps?", 

489 default=False, 

490 ) 

491 overscanBiasJumpKeyword = pexConfig.Field( 

492 dtype=str, 

493 doc="Header keyword containing information about devices.", 

494 default="NO_SUCH_KEY", 

495 ) 

496 overscanBiasJumpDevices = pexConfig.ListField( 

497 dtype=str, 

498 doc="List of devices that need piecewise overscan correction.", 

499 default=(), 

500 ) 

501 overscanBiasJumpLocation = pexConfig.Field( 

502 dtype=int, 

503 doc="Location of bias jump along y-axis.", 

504 default=0, 

505 ) 

506 

507 # Amplifier to CCD assembly configuration 

508 doAssembleCcd = pexConfig.Field( 

509 dtype=bool, 

510 default=True, 

511 doc="Assemble amp-level exposures into a ccd-level exposure?" 

512 ) 

513 assembleCcd = pexConfig.ConfigurableField( 

514 target=AssembleCcdTask, 

515 doc="CCD assembly task", 

516 ) 

517 

518 # General calibration configuration. 

519 doAssembleIsrExposures = pexConfig.Field( 

520 dtype=bool, 

521 default=False, 

522 doc="Assemble amp-level calibration exposures into ccd-level exposure?" 

523 ) 

524 doTrimToMatchCalib = pexConfig.Field( 

525 dtype=bool, 

526 default=False, 

527 doc="Trim raw data to match calibration bounding boxes?" 

528 ) 

529 

530 # Bias subtraction. 

531 doBias = pexConfig.Field( 

532 dtype=bool, 

533 doc="Apply bias frame correction?", 

534 default=True, 

535 ) 

536 biasDataProductName = pexConfig.Field( 

537 dtype=str, 

538 doc="Name of the bias data product", 

539 default="bias", 

540 ) 

541 doBiasBeforeOverscan = pexConfig.Field( 

542 dtype=bool, 

543 doc="Reverse order of overscan and bias correction.", 

544 default=False 

545 ) 

546 

547 # Variance construction 

548 doVariance = pexConfig.Field( 

549 dtype=bool, 

550 doc="Calculate variance?", 

551 default=True 

552 ) 

553 gain = pexConfig.Field( 

554 dtype=float, 

555 doc="The gain to use if no Detector is present in the Exposure (ignored if NaN)", 

556 default=float("NaN"), 

557 ) 

558 readNoise = pexConfig.Field( 

559 dtype=float, 

560 doc="The read noise to use if no Detector is present in the Exposure", 

561 default=0.0, 

562 ) 

563 doEmpiricalReadNoise = pexConfig.Field( 

564 dtype=bool, 

565 default=False, 

566 doc="Calculate empirical read noise instead of value from AmpInfo data?" 

567 ) 

568 

569 # Linearization. 

570 doLinearize = pexConfig.Field( 

571 dtype=bool, 

572 doc="Correct for nonlinearity of the detector's response?", 

573 default=True, 

574 ) 

575 

576 # Crosstalk. 

577 doCrosstalk = pexConfig.Field( 

578 dtype=bool, 

579 doc="Apply intra-CCD crosstalk correction?", 

580 default=False, 

581 ) 

582 doCrosstalkBeforeAssemble = pexConfig.Field( 

583 dtype=bool, 

584 doc="Apply crosstalk correction before CCD assembly, and before trimming?", 

585 default=False, 

586 ) 

587 crosstalk = pexConfig.ConfigurableField( 

588 target=CrosstalkTask, 

589 doc="Intra-CCD crosstalk correction", 

590 ) 

591 

592 # Masking options. 

593 doDefect = pexConfig.Field( 

594 dtype=bool, 

595 doc="Apply correction for CCD defects, e.g. hot pixels?", 

596 default=True, 

597 ) 

598 doNanMasking = pexConfig.Field( 

599 dtype=bool, 

600 doc="Mask non-finite (NAN, inf) pixels?", 

601 default=True, 

602 ) 

603 doWidenSaturationTrails = pexConfig.Field( 

604 dtype=bool, 

605 doc="Widen bleed trails based on their width?", 

606 default=True 

607 ) 

608 

609 # Brighter-Fatter correction. 

610 doBrighterFatter = pexConfig.Field( 

611 dtype=bool, 

612 default=False, 

613 doc="Apply the brighter-fatter correction?" 

614 ) 

615 brighterFatterLevel = pexConfig.ChoiceField( 

616 dtype=str, 

617 default="DETECTOR", 

618 doc="The level at which to correct for brighter-fatter.", 

619 allowed={ 

620 "AMP": "Every amplifier treated separately.", 

621 "DETECTOR": "One kernel per detector", 

622 } 

623 ) 

624 brighterFatterMaxIter = pexConfig.Field( 

625 dtype=int, 

626 default=10, 

627 doc="Maximum number of iterations for the brighter-fatter correction" 

628 ) 

629 brighterFatterThreshold = pexConfig.Field( 

630 dtype=float, 

631 default=1000, 

632 doc="Threshold used to stop iterating the brighter-fatter correction. It is the " 

633 "absolute value of the difference between the current corrected image and the one " 

634 "from the previous iteration summed over all the pixels." 

635 ) 

636 brighterFatterApplyGain = pexConfig.Field( 

637 dtype=bool, 

638 default=True, 

639 doc="Should the gain be applied when applying the brighter-fatter correction?" 

640 ) 

641 brighterFatterMaskListToInterpolate = pexConfig.ListField( 

642 dtype=str, 

643 doc="List of mask planes that should be interpolated over when applying the brighter-fatter " 

644 "correction.", 

645 default=["SAT", "BAD", "NO_DATA", "UNMASKEDNAN"], 

646 ) 

647 brighterFatterMaskGrowSize = pexConfig.Field( 

648 dtype=int, 

649 default=0, 

650 doc="Number of pixels to grow the masks listed in config.brighterFatterMaskListToInterpolate " 

651 "when brighter-fatter correction is applied." 

652 ) 

653 

654 # Dark subtraction. 

655 doDark = pexConfig.Field( 

656 dtype=bool, 

657 doc="Apply dark frame correction?", 

658 default=True, 

659 ) 

660 darkDataProductName = pexConfig.Field( 

661 dtype=str, 

662 doc="Name of the dark data product", 

663 default="dark", 

664 ) 

665 

666 # Camera-specific stray light removal. 

667 doStrayLight = pexConfig.Field( 

668 dtype=bool, 

669 doc="Subtract stray light in the y-band (due to encoder LEDs)?", 

670 default=False, 

671 ) 

672 strayLight = pexConfig.ConfigurableField( 

673 target=StrayLightTask, 

674 doc="y-band stray light correction" 

675 ) 

676 

677 # Flat correction. 

678 doFlat = pexConfig.Field( 

679 dtype=bool, 

680 doc="Apply flat field correction?", 

681 default=True, 

682 ) 

683 flatDataProductName = pexConfig.Field( 

684 dtype=str, 

685 doc="Name of the flat data product", 

686 default="flat", 

687 ) 

688 flatScalingType = pexConfig.ChoiceField( 

689 dtype=str, 

690 doc="The method for scaling the flat on the fly.", 

691 default='USER', 

692 allowed={ 

693 "USER": "Scale by flatUserScale", 

694 "MEAN": "Scale by the inverse of the mean", 

695 "MEDIAN": "Scale by the inverse of the median", 

696 }, 

697 ) 

698 flatUserScale = pexConfig.Field( 

699 dtype=float, 

700 doc="If flatScalingType is 'USER' then scale flat by this amount; ignored otherwise", 

701 default=1.0, 

702 ) 

703 doTweakFlat = pexConfig.Field( 

704 dtype=bool, 

705 doc="Tweak flats to match observed amplifier ratios?", 

706 default=False 

707 ) 

708 

709 # Amplifier normalization based on gains instead of using flats configuration. 

710 doApplyGains = pexConfig.Field( 

711 dtype=bool, 

712 doc="Correct the amplifiers for their gains instead of applying flat correction", 

713 default=False, 

714 ) 

715 normalizeGains = pexConfig.Field( 

716 dtype=bool, 

717 doc="Normalize all the amplifiers in each CCD to have the same median value.", 

718 default=False, 

719 ) 

720 

721 # Fringe correction. 

722 doFringe = pexConfig.Field( 

723 dtype=bool, 

724 doc="Apply fringe correction?", 

725 default=True, 

726 ) 

727 fringe = pexConfig.ConfigurableField( 

728 target=FringeTask, 

729 doc="Fringe subtraction task", 

730 ) 

731 fringeAfterFlat = pexConfig.Field( 

732 dtype=bool, 

733 doc="Do fringe subtraction after flat-fielding?", 

734 default=True, 

735 ) 

736 

737 # Initial CCD-level background statistics options. 

738 doMeasureBackground = pexConfig.Field( 

739 dtype=bool, 

740 doc="Measure the background level on the reduced image?", 

741 default=False, 

742 ) 

743 

744 # Camera-specific masking configuration. 

745 doCameraSpecificMasking = pexConfig.Field( 

746 dtype=bool, 

747 doc="Mask camera-specific bad regions?", 

748 default=False, 

749 ) 

750 masking = pexConfig.ConfigurableField( 

751 target=MaskingTask, 

752 doc="Masking task." 

753 ) 

754 

755 # Interpolation options. 

756 

757 doInterpolate = pexConfig.Field( 

758 dtype=bool, 

759 doc="Interpolate masked pixels?", 

760 default=True, 

761 ) 

762 doSaturationInterpolation = pexConfig.Field( 

763 dtype=bool, 

764 doc="Perform interpolation over pixels masked as saturated?" 

765 " NB: This is independent of doSaturation; if that is False this plane" 

766 " will likely be blank, resulting in a no-op here.", 

767 default=True, 

768 ) 

769 doNanInterpolation = pexConfig.Field( 

770 dtype=bool, 

771 doc="Perform interpolation over pixels masked as NaN?" 

772 " NB: This is independent of doNanMasking; if that is False this plane" 

773 " will likely be blank, resulting in a no-op here.", 

774 default=True, 

775 ) 

776 doNanInterpAfterFlat = pexConfig.Field( 

777 dtype=bool, 

778 doc=("If True, ensure we interpolate NaNs after flat-fielding, even if we " 

779 "also have to interpolate them before flat-fielding."), 

780 default=False, 

781 ) 

782 maskListToInterpolate = pexConfig.ListField( 

783 dtype=str, 

784 doc="List of mask planes that should be interpolated.", 

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

786 ) 

787 doSaveInterpPixels = pexConfig.Field( 

788 dtype=bool, 

789 doc="Save a copy of the pre-interpolated pixel values?", 

790 default=False, 

791 ) 

792 

793 # Default photometric calibration options. 

794 fluxMag0T1 = pexConfig.DictField( 

795 keytype=str, 

796 itemtype=float, 

797 doc="The approximate flux of a zero-magnitude object in a one-second exposure, per filter.", 

798 default=dict((f, pow(10.0, 0.4*m)) for f, m in (("Unknown", 28.0), 

799 )) 

800 ) 

801 defaultFluxMag0T1 = pexConfig.Field( 

802 dtype=float, 

803 doc="Default value for fluxMag0T1 (for an unrecognized filter).", 

804 default=pow(10.0, 0.4*28.0) 

805 ) 

806 

807 # Vignette correction configuration. 

808 doVignette = pexConfig.Field( 

809 dtype=bool, 

810 doc="Apply vignetting parameters?", 

811 default=False, 

812 ) 

813 vignette = pexConfig.ConfigurableField( 

814 target=VignetteTask, 

815 doc="Vignetting task.", 

816 ) 

817 

818 # Transmission curve configuration. 

819 doAttachTransmissionCurve = pexConfig.Field( 

820 dtype=bool, 

821 default=False, 

822 doc="Construct and attach a wavelength-dependent throughput curve for this CCD image?" 

823 ) 

824 doUseOpticsTransmission = pexConfig.Field( 

825 dtype=bool, 

826 default=True, 

827 doc="Load and use transmission_optics (if doAttachTransmissionCurve is True)?" 

828 ) 

829 doUseFilterTransmission = pexConfig.Field( 

830 dtype=bool, 

831 default=True, 

832 doc="Load and use transmission_filter (if doAttachTransmissionCurve is True)?" 

833 ) 

834 doUseSensorTransmission = pexConfig.Field( 

835 dtype=bool, 

836 default=True, 

837 doc="Load and use transmission_sensor (if doAttachTransmissionCurve is True)?" 

838 ) 

839 doUseAtmosphereTransmission = pexConfig.Field( 

840 dtype=bool, 

841 default=True, 

842 doc="Load and use transmission_atmosphere (if doAttachTransmissionCurve is True)?" 

843 ) 

844 

845 # Illumination correction. 

846 doIlluminationCorrection = pexConfig.Field( 

847 dtype=bool, 

848 default=False, 

849 doc="Perform illumination correction?" 

850 ) 

851 illuminationCorrectionDataProductName = pexConfig.Field( 

852 dtype=str, 

853 doc="Name of the illumination correction data product.", 

854 default="illumcor", 

855 ) 

856 illumScale = pexConfig.Field( 

857 dtype=float, 

858 doc="Scale factor for the illumination correction.", 

859 default=1.0, 

860 ) 

861 illumFilters = pexConfig.ListField( 

862 dtype=str, 

863 default=[], 

864 doc="Only perform illumination correction for these filters." 

865 ) 

866 

867 # Write the outputs to disk. If ISR is run as a subtask, this may not be needed. 

868 doWrite = pexConfig.Field( 

869 dtype=bool, 

870 doc="Persist postISRCCD?", 

871 default=True, 

872 ) 

873 

874 def validate(self): 

875 super().validate() 

876 if self.doFlat and self.doApplyGains: 

877 raise ValueError("You may not specify both doFlat and doApplyGains") 

878 if self.doBiasBeforeOverscan and self.doTrimToMatchCalib: 

879 raise ValueError("You may not specify both doBiasBeforeOverscan and doTrimToMatchCalib") 

880 if self.doSaturationInterpolation and self.saturatedMaskName not in self.maskListToInterpolate: 

881 self.maskListToInterpolate.append(self.saturatedMaskName) 

882 if not self.doSaturationInterpolation and self.saturatedMaskName in self.maskListToInterpolate: 

883 self.maskListToInterpolate.remove(self.saturatedMaskName) 

884 if self.doNanInterpolation and "UNMASKEDNAN" not in self.maskListToInterpolate: 

885 self.maskListToInterpolate.append("UNMASKEDNAN") 

886 

887 

888class IsrTask(pipeBase.PipelineTask, pipeBase.CmdLineTask): 

889 """Apply common instrument signature correction algorithms to a raw frame. 

890 

891 The process for correcting imaging data is very similar from 

892 camera to camera. This task provides a vanilla implementation of 

893 doing these corrections, including the ability to turn certain 

894 corrections off if they are not needed. The inputs to the primary 

895 method, `run()`, are a raw exposure to be corrected and the 

896 calibration data products. The raw input is a single chip sized 

897 mosaic of all amps including overscans and other non-science 

898 pixels. The method `runDataRef()` identifies and defines the 

899 calibration data products, and is intended for use by a 

900 `lsst.pipe.base.cmdLineTask.CmdLineTask` and takes as input only a 

901 `daf.persistence.butlerSubset.ButlerDataRef`. This task may be 

902 subclassed for different camera, although the most camera specific 

903 methods have been split into subtasks that can be redirected 

904 appropriately. 

905 

906 The __init__ method sets up the subtasks for ISR processing, using 

907 the defaults from `lsst.ip.isr`. 

908 

909 Parameters 

910 ---------- 

911 args : `list` 

912 Positional arguments passed to the Task constructor. None used at this time. 

913 kwargs : `dict`, optional 

914 Keyword arguments passed on to the Task constructor. None used at this time. 

915 """ 

916 ConfigClass = IsrTaskConfig 

917 _DefaultName = "isr" 

918 

919 def __init__(self, **kwargs): 

920 super().__init__(**kwargs) 

921 self.makeSubtask("assembleCcd") 

922 self.makeSubtask("crosstalk") 

923 self.makeSubtask("strayLight") 

924 self.makeSubtask("fringe") 

925 self.makeSubtask("masking") 

926 self.makeSubtask("overscan") 

927 self.makeSubtask("vignette") 

928 

929 def runQuantum(self, butlerQC, inputRefs, outputRefs): 

930 inputs = butlerQC.get(inputRefs) 

931 

932 try: 

933 inputs['detectorNum'] = inputRefs.ccdExposure.dataId['detector'] 

934 except Exception as e: 

935 raise ValueError("Failure to find valid detectorNum value for Dataset %s: %s." % 

936 (inputRefs, e)) 

937 

938 inputs['isGen3'] = True 

939 

940 detector = inputs['ccdExposure'].getDetector() 

941 

942 if self.config.doCrosstalk is True: 

943 # Crosstalk sources need to be defined by the pipeline 

944 # yaml if they exist. 

945 if 'crosstalk' in inputs and inputs['crosstalk'] is not None: 

946 if not isinstance(inputs['crosstalk'], CrosstalkCalib): 

947 inputs['crosstalk'] = CrosstalkCalib.fromTable(inputs['crosstalk']) 

948 else: 

949 coeffVector = (self.config.crosstalk.crosstalkValues 

950 if self.config.crosstalk.useConfigCoefficients else None) 

951 crosstalkCalib = CrosstalkCalib().fromDetector(detector, coeffVector=coeffVector) 

952 inputs['crosstalk'] = crosstalkCalib 

953 if inputs['crosstalk'].interChip and len(inputs['crosstalk'].interChip) > 0: 

954 if 'crosstalkSources' not in inputs: 

955 self.log.warn("No crosstalkSources found for chip with interChip terms!") 

956 

957 if self.doLinearize(detector) is True: 

958 if 'linearizer' in inputs: 

959 if isinstance(inputs['linearizer'], dict): 

960 linearizer = linearize.Linearizer(detector=detector, log=self.log) 

961 linearizer.fromYaml(inputs['linearizer']) 

962 self.log.warn("Dictionary linearizers will be deprecated in DM-28741.") 

963 elif isinstance(inputs['linearizer'], numpy.ndarray): 

964 linearizer = linearize.Linearizer(table=inputs.get('linearizer', None), 

965 detector=detector, 

966 log=self.log) 

967 self.log.warn("Bare lookup table linearizers will be deprecated in DM-28741.") 

968 else: 

969 linearizer = inputs['linearizer'] 

970 linearizer.log = self.log 

971 inputs['linearizer'] = linearizer 

972 else: 

973 inputs['linearizer'] = linearize.Linearizer(detector=detector, log=self.log) 

974 self.log.warn("Constructing linearizer from cameraGeom information.") 

975 

976 if self.config.doDefect is True: 

977 if "defects" in inputs and inputs['defects'] is not None: 

978 # defects is loaded as a BaseCatalog with columns x0, y0, width, height. 

979 # masking expects a list of defects defined by their bounding box 

980 if not isinstance(inputs["defects"], Defects): 

981 inputs["defects"] = Defects.fromTable(inputs["defects"]) 

982 

983 # Load the correct style of brighter-fatter kernel, and repack 

984 # the information as a numpy array. 

985 if self.config.doBrighterFatter: 

986 brighterFatterKernel = inputs.pop('newBFKernel', None) 

987 if brighterFatterKernel is None: 

988 brighterFatterKernel = inputs.get('bfKernel', None) 

989 

990 if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray): 

991 # This is a ISR calib kernel 

992 detName = detector.getName() 

993 level = brighterFatterKernel.level 

994 

995 # This is expected to be a dictionary of amp-wise gains. 

996 inputs['bfGains'] = brighterFatterKernel.gain 

997 if self.config.brighterFatterLevel == 'DETECTOR': 

998 if level == 'DETECTOR': 

999 if detName in brighterFatterKernel.detKernels: 

1000 inputs['bfKernel'] = brighterFatterKernel.detKernels[detName] 

1001 else: 

1002 raise RuntimeError("Failed to extract kernel from new-style BF kernel.") 

1003 elif level == 'AMP': 

1004 self.log.warn("Making DETECTOR level kernel from AMP based brighter fatter kernels.") 

1005 brighterFatterKernel.makeDetectorKernelFromAmpwiseKernels(detName) 

1006 inputs['bfKernel'] = brighterFatterKernel.detKernels[detName] 

1007 elif self.config.brighterFatterLevel == 'AMP': 

1008 raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented") 

1009 

1010 if self.config.doFringe is True and self.fringe.checkFilter(inputs['ccdExposure']): 

1011 expId = inputs['ccdExposure'].getInfo().getVisitInfo().getExposureId() 

1012 inputs['fringes'] = self.fringe.loadFringes(inputs['fringes'], 

1013 expId=expId, 

1014 assembler=self.assembleCcd 

1015 if self.config.doAssembleIsrExposures else None) 

1016 else: 

1017 inputs['fringes'] = pipeBase.Struct(fringes=None) 

1018 

1019 if self.config.doStrayLight is True and self.strayLight.checkFilter(inputs['ccdExposure']): 

1020 if 'strayLightData' not in inputs: 

1021 inputs['strayLightData'] = None 

1022 

1023 outputs = self.run(**inputs) 

1024 butlerQC.put(outputs, outputRefs) 

1025 

1026 def readIsrData(self, dataRef, rawExposure): 

1027 """Retrieve necessary frames for instrument signature removal. 

1028 

1029 Pre-fetching all required ISR data products limits the IO 

1030 required by the ISR. Any conflict between the calibration data 

1031 available and that needed for ISR is also detected prior to 

1032 doing processing, allowing it to fail quickly. 

1033 

1034 Parameters 

1035 ---------- 

1036 dataRef : `daf.persistence.butlerSubset.ButlerDataRef` 

1037 Butler reference of the detector data to be processed 

1038 rawExposure : `afw.image.Exposure` 

1039 The raw exposure that will later be corrected with the 

1040 retrieved calibration data; should not be modified in this 

1041 method. 

1042 

1043 Returns 

1044 ------- 

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

1046 Result struct with components (which may be `None`): 

1047 - ``bias``: bias calibration frame (`afw.image.Exposure`) 

1048 - ``linearizer``: functor for linearization (`ip.isr.linearize.LinearizeBase`) 

1049 - ``crosstalkSources``: list of possible crosstalk sources (`list`) 

1050 - ``dark``: dark calibration frame (`afw.image.Exposure`) 

1051 - ``flat``: flat calibration frame (`afw.image.Exposure`) 

1052 - ``bfKernel``: Brighter-Fatter kernel (`numpy.ndarray`) 

1053 - ``defects``: list of defects (`lsst.ip.isr.Defects`) 

1054 - ``fringes``: `lsst.pipe.base.Struct` with components: 

1055 - ``fringes``: fringe calibration frame (`afw.image.Exposure`) 

1056 - ``seed``: random seed derived from the ccdExposureId for random 

1057 number generator (`uint32`). 

1058 - ``opticsTransmission``: `lsst.afw.image.TransmissionCurve` 

1059 A ``TransmissionCurve`` that represents the throughput of the optics, 

1060 to be evaluated in focal-plane coordinates. 

1061 - ``filterTransmission`` : `lsst.afw.image.TransmissionCurve` 

1062 A ``TransmissionCurve`` that represents the throughput of the filter 

1063 itself, to be evaluated in focal-plane coordinates. 

1064 - ``sensorTransmission`` : `lsst.afw.image.TransmissionCurve` 

1065 A ``TransmissionCurve`` that represents the throughput of the sensor 

1066 itself, to be evaluated in post-assembly trimmed detector coordinates. 

1067 - ``atmosphereTransmission`` : `lsst.afw.image.TransmissionCurve` 

1068 A ``TransmissionCurve`` that represents the throughput of the 

1069 atmosphere, assumed to be spatially constant. 

1070 - ``strayLightData`` : `object` 

1071 An opaque object containing calibration information for 

1072 stray-light correction. If `None`, no correction will be 

1073 performed. 

1074 - ``illumMaskedImage`` : illumination correction image (`lsst.afw.image.MaskedImage`) 

1075 

1076 Raises 

1077 ------ 

1078 NotImplementedError : 

1079 Raised if a per-amplifier brighter-fatter kernel is requested by the configuration. 

1080 """ 

1081 try: 

1082 dateObs = rawExposure.getInfo().getVisitInfo().getDate() 

1083 dateObs = dateObs.toPython().isoformat() 

1084 except RuntimeError: 

1085 self.log.warn("Unable to identify dateObs for rawExposure.") 

1086 dateObs = None 

1087 

1088 ccd = rawExposure.getDetector() 

1089 filterLabel = rawExposure.getFilterLabel() 

1090 rawExposure.mask.addMaskPlane("UNMASKEDNAN") # needed to match pre DM-15862 processing. 

1091 biasExposure = (self.getIsrExposure(dataRef, self.config.biasDataProductName) 

1092 if self.config.doBias else None) 

1093 # immediate=True required for functors and linearizers are functors; see ticket DM-6515 

1094 linearizer = (dataRef.get("linearizer", immediate=True) 

1095 if self.doLinearize(ccd) else None) 

1096 if linearizer is not None and not isinstance(linearizer, numpy.ndarray): 

1097 linearizer.log = self.log 

1098 if isinstance(linearizer, numpy.ndarray): 

1099 linearizer = linearize.Linearizer(table=linearizer, detector=ccd) 

1100 

1101 crosstalkCalib = None 

1102 if self.config.doCrosstalk: 

1103 try: 

1104 crosstalkCalib = dataRef.get("crosstalk", immediate=True) 

1105 except NoResults: 

1106 coeffVector = (self.config.crosstalk.crosstalkValues 

1107 if self.config.crosstalk.useConfigCoefficients else None) 

1108 crosstalkCalib = CrosstalkCalib().fromDetector(ccd, coeffVector=coeffVector) 

1109 crosstalkSources = (self.crosstalk.prepCrosstalk(dataRef, crosstalkCalib) 

1110 if self.config.doCrosstalk else None) 

1111 

1112 darkExposure = (self.getIsrExposure(dataRef, self.config.darkDataProductName) 

1113 if self.config.doDark else None) 

1114 flatExposure = (self.getIsrExposure(dataRef, self.config.flatDataProductName, 

1115 dateObs=dateObs) 

1116 if self.config.doFlat else None) 

1117 

1118 brighterFatterKernel = None 

1119 brighterFatterGains = None 

1120 if self.config.doBrighterFatter is True: 

1121 try: 

1122 # Use the new-style cp_pipe version of the kernel if it exists 

1123 # If using a new-style kernel, always use the self-consistent 

1124 # gains, i.e. the ones inside the kernel object itself 

1125 brighterFatterKernel = dataRef.get("brighterFatterKernel") 

1126 brighterFatterGains = brighterFatterKernel.gain 

1127 self.log.info("New style brighter-fatter kernel (brighterFatterKernel) loaded") 

1128 except NoResults: 

1129 try: # Fall back to the old-style numpy-ndarray style kernel if necessary. 

1130 brighterFatterKernel = dataRef.get("bfKernel") 

1131 self.log.info("Old style brighter-fatter kernel (np.array) loaded") 

1132 except NoResults: 

1133 brighterFatterKernel = None 

1134 if brighterFatterKernel is not None and not isinstance(brighterFatterKernel, numpy.ndarray): 

1135 # If the kernel is not an ndarray, it's the cp_pipe version 

1136 # so extract the kernel for this detector, or raise an error 

1137 if self.config.brighterFatterLevel == 'DETECTOR': 

1138 if brighterFatterKernel.detectorKernel: 

1139 brighterFatterKernel = brighterFatterKernel.detectorKernel[ccd.getId()] 

1140 else: 

1141 raise RuntimeError("Failed to extract kernel from new-style BF kernel.") 

1142 else: 

1143 # TODO DM-15631 for implementing this 

1144 raise NotImplementedError("Per-amplifier brighter-fatter correction not implemented") 

1145 

1146 defectList = (dataRef.get("defects") 

1147 if self.config.doDefect else None) 

1148 fringeStruct = (self.fringe.readFringes(dataRef, assembler=self.assembleCcd 

1149 if self.config.doAssembleIsrExposures else None) 

1150 if self.config.doFringe and self.fringe.checkFilter(rawExposure) 

1151 else pipeBase.Struct(fringes=None)) 

1152 

1153 if self.config.doAttachTransmissionCurve: 

1154 opticsTransmission = (dataRef.get("transmission_optics") 

1155 if self.config.doUseOpticsTransmission else None) 

1156 filterTransmission = (dataRef.get("transmission_filter") 

1157 if self.config.doUseFilterTransmission else None) 

1158 sensorTransmission = (dataRef.get("transmission_sensor") 

1159 if self.config.doUseSensorTransmission else None) 

1160 atmosphereTransmission = (dataRef.get("transmission_atmosphere") 

1161 if self.config.doUseAtmosphereTransmission else None) 

1162 else: 

1163 opticsTransmission = None 

1164 filterTransmission = None 

1165 sensorTransmission = None 

1166 atmosphereTransmission = None 

1167 

1168 if self.config.doStrayLight: 

1169 strayLightData = self.strayLight.readIsrData(dataRef, rawExposure) 

1170 else: 

1171 strayLightData = None 

1172 

1173 illumMaskedImage = (self.getIsrExposure(dataRef, 

1174 self.config.illuminationCorrectionDataProductName).getMaskedImage() 

1175 if (self.config.doIlluminationCorrection 

1176 and filterLabel in self.config.illumFilters) 

1177 else None) 

1178 

1179 # Struct should include only kwargs to run() 

1180 return pipeBase.Struct(bias=biasExposure, 

1181 linearizer=linearizer, 

1182 crosstalk=crosstalkCalib, 

1183 crosstalkSources=crosstalkSources, 

1184 dark=darkExposure, 

1185 flat=flatExposure, 

1186 bfKernel=brighterFatterKernel, 

1187 bfGains=brighterFatterGains, 

1188 defects=defectList, 

1189 fringes=fringeStruct, 

1190 opticsTransmission=opticsTransmission, 

1191 filterTransmission=filterTransmission, 

1192 sensorTransmission=sensorTransmission, 

1193 atmosphereTransmission=atmosphereTransmission, 

1194 strayLightData=strayLightData, 

1195 illumMaskedImage=illumMaskedImage 

1196 ) 

1197 

1198 @pipeBase.timeMethod 

1199 def run(self, ccdExposure, camera=None, bias=None, linearizer=None, 

1200 crosstalk=None, crosstalkSources=None, 

1201 dark=None, flat=None, bfKernel=None, bfGains=None, defects=None, 

1202 fringes=pipeBase.Struct(fringes=None), opticsTransmission=None, filterTransmission=None, 

1203 sensorTransmission=None, atmosphereTransmission=None, 

1204 detectorNum=None, strayLightData=None, illumMaskedImage=None, 

1205 isGen3=False, 

1206 ): 

1207 """Perform instrument signature removal on an exposure. 

1208 

1209 Steps included in the ISR processing, in order performed, are: 

1210 - saturation and suspect pixel masking 

1211 - overscan subtraction 

1212 - CCD assembly of individual amplifiers 

1213 - bias subtraction 

1214 - variance image construction 

1215 - linearization of non-linear response 

1216 - crosstalk masking 

1217 - brighter-fatter correction 

1218 - dark subtraction 

1219 - fringe correction 

1220 - stray light subtraction 

1221 - flat correction 

1222 - masking of known defects and camera specific features 

1223 - vignette calculation 

1224 - appending transmission curve and distortion model 

1225 

1226 Parameters 

1227 ---------- 

1228 ccdExposure : `lsst.afw.image.Exposure` 

1229 The raw exposure that is to be run through ISR. The 

1230 exposure is modified by this method. 

1231 camera : `lsst.afw.cameraGeom.Camera`, optional 

1232 The camera geometry for this exposure. Required if ``isGen3`` is 

1233 `True` and one or more of ``ccdExposure``, ``bias``, ``dark``, or 

1234 ``flat`` does not have an associated detector. 

1235 bias : `lsst.afw.image.Exposure`, optional 

1236 Bias calibration frame. 

1237 linearizer : `lsst.ip.isr.linearize.LinearizeBase`, optional 

1238 Functor for linearization. 

1239 crosstalk : `lsst.ip.isr.crosstalk.CrosstalkCalib`, optional 

1240 Calibration for crosstalk. 

1241 crosstalkSources : `list`, optional 

1242 List of possible crosstalk sources. 

1243 dark : `lsst.afw.image.Exposure`, optional 

1244 Dark calibration frame. 

1245 flat : `lsst.afw.image.Exposure`, optional 

1246 Flat calibration frame. 

1247 bfKernel : `numpy.ndarray`, optional 

1248 Brighter-fatter kernel. 

1249 bfGains : `dict` of `float`, optional 

1250 Gains used to override the detector's nominal gains for the 

1251 brighter-fatter correction. A dict keyed by amplifier name for 

1252 the detector in question. 

1253 defects : `lsst.ip.isr.Defects`, optional 

1254 List of defects. 

1255 fringes : `lsst.pipe.base.Struct`, optional 

1256 Struct containing the fringe correction data, with 

1257 elements: 

1258 - ``fringes``: fringe calibration frame (`afw.image.Exposure`) 

1259 - ``seed``: random seed derived from the ccdExposureId for random 

1260 number generator (`uint32`) 

1261 opticsTransmission: `lsst.afw.image.TransmissionCurve`, optional 

1262 A ``TransmissionCurve`` that represents the throughput of the optics, 

1263 to be evaluated in focal-plane coordinates. 

1264 filterTransmission : `lsst.afw.image.TransmissionCurve` 

1265 A ``TransmissionCurve`` that represents the throughput of the filter 

1266 itself, to be evaluated in focal-plane coordinates. 

1267 sensorTransmission : `lsst.afw.image.TransmissionCurve` 

1268 A ``TransmissionCurve`` that represents the throughput of the sensor 

1269 itself, to be evaluated in post-assembly trimmed detector coordinates. 

1270 atmosphereTransmission : `lsst.afw.image.TransmissionCurve` 

1271 A ``TransmissionCurve`` that represents the throughput of the 

1272 atmosphere, assumed to be spatially constant. 

1273 detectorNum : `int`, optional 

1274 The integer number for the detector to process. 

1275 isGen3 : bool, optional 

1276 Flag this call to run() as using the Gen3 butler environment. 

1277 strayLightData : `object`, optional 

1278 Opaque object containing calibration information for stray-light 

1279 correction. If `None`, no correction will be performed. 

1280 illumMaskedImage : `lsst.afw.image.MaskedImage`, optional 

1281 Illumination correction image. 

1282 

1283 Returns 

1284 ------- 

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

1286 Result struct with component: 

1287 - ``exposure`` : `afw.image.Exposure` 

1288 The fully ISR corrected exposure. 

1289 - ``outputExposure`` : `afw.image.Exposure` 

1290 An alias for `exposure` 

1291 - ``ossThumb`` : `numpy.ndarray` 

1292 Thumbnail image of the exposure after overscan subtraction. 

1293 - ``flattenedThumb`` : `numpy.ndarray` 

1294 Thumbnail image of the exposure after flat-field correction. 

1295 

1296 Raises 

1297 ------ 

1298 RuntimeError 

1299 Raised if a configuration option is set to True, but the 

1300 required calibration data has not been specified. 

1301 

1302 Notes 

1303 ----- 

1304 The current processed exposure can be viewed by setting the 

1305 appropriate lsstDebug entries in the `debug.display` 

1306 dictionary. The names of these entries correspond to some of 

1307 the IsrTaskConfig Boolean options, with the value denoting the 

1308 frame to use. The exposure is shown inside the matching 

1309 option check and after the processing of that step has 

1310 finished. The steps with debug points are: 

1311 

1312 doAssembleCcd 

1313 doBias 

1314 doCrosstalk 

1315 doBrighterFatter 

1316 doDark 

1317 doFringe 

1318 doStrayLight 

1319 doFlat 

1320 

1321 In addition, setting the "postISRCCD" entry displays the 

1322 exposure after all ISR processing has finished. 

1323 

1324 """ 

1325 

1326 if isGen3 is True: 

1327 # Gen3 currently cannot automatically do configuration overrides. 

1328 # DM-15257 looks to discuss this issue. 

1329 # Configure input exposures; 

1330 if detectorNum is None: 

1331 raise RuntimeError("Must supply the detectorNum if running as Gen3.") 

1332 

1333 ccdExposure = self.ensureExposure(ccdExposure, camera, detectorNum) 

1334 bias = self.ensureExposure(bias, camera, detectorNum) 

1335 dark = self.ensureExposure(dark, camera, detectorNum) 

1336 flat = self.ensureExposure(flat, camera, detectorNum) 

1337 else: 

1338 if isinstance(ccdExposure, ButlerDataRef): 

1339 return self.runDataRef(ccdExposure) 

1340 

1341 ccd = ccdExposure.getDetector() 

1342 filterLabel = ccdExposure.getFilterLabel() 

1343 

1344 if not ccd: 

1345 assert not self.config.doAssembleCcd, "You need a Detector to run assembleCcd." 

1346 ccd = [FakeAmp(ccdExposure, self.config)] 

1347 

1348 # Validate Input 

1349 if self.config.doBias and bias is None: 

1350 raise RuntimeError("Must supply a bias exposure if config.doBias=True.") 

1351 if self.doLinearize(ccd) and linearizer is None: 

1352 raise RuntimeError("Must supply a linearizer if config.doLinearize=True for this detector.") 

1353 if self.config.doBrighterFatter and bfKernel is None: 

1354 raise RuntimeError("Must supply a kernel if config.doBrighterFatter=True.") 

1355 if self.config.doDark and dark is None: 

1356 raise RuntimeError("Must supply a dark exposure if config.doDark=True.") 

1357 if self.config.doFlat and flat is None: 

1358 raise RuntimeError("Must supply a flat exposure if config.doFlat=True.") 

1359 if self.config.doDefect and defects is None: 

1360 raise RuntimeError("Must supply defects if config.doDefect=True.") 

1361 if (self.config.doFringe and filterLabel in self.fringe.config.filters 

1362 and fringes.fringes is None): 

1363 # The `fringes` object needs to be a pipeBase.Struct, as 

1364 # we use it as a `dict` for the parameters of 

1365 # `FringeTask.run()`. The `fringes.fringes` `list` may 

1366 # not be `None` if `doFringe=True`. Otherwise, raise. 

1367 raise RuntimeError("Must supply fringe exposure as a pipeBase.Struct.") 

1368 if (self.config.doIlluminationCorrection and filterLabel in self.config.illumFilters 

1369 and illumMaskedImage is None): 

1370 raise RuntimeError("Must supply an illumcor if config.doIlluminationCorrection=True.") 

1371 

1372 # Begin ISR processing. 

1373 if self.config.doConvertIntToFloat: 

1374 self.log.info("Converting exposure to floating point values.") 

1375 ccdExposure = self.convertIntToFloat(ccdExposure) 

1376 

1377 if self.config.doBias and self.config.doBiasBeforeOverscan: 

1378 self.log.info("Applying bias correction.") 

1379 isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(), 

1380 trimToFit=self.config.doTrimToMatchCalib) 

1381 self.debugView(ccdExposure, "doBias") 

1382 

1383 # Amplifier level processing. 

1384 overscans = [] 

1385 for amp in ccd: 

1386 # if ccdExposure is one amp, check for coverage to prevent performing ops multiple times 

1387 if ccdExposure.getBBox().contains(amp.getBBox()): 

1388 # Check for fully masked bad amplifiers, and generate masks for SUSPECT and SATURATED values. 

1389 badAmp = self.maskAmplifier(ccdExposure, amp, defects) 

1390 

1391 if self.config.doOverscan and not badAmp: 

1392 # Overscan correction on amp-by-amp basis. 

1393 overscanResults = self.overscanCorrection(ccdExposure, amp) 

1394 self.log.debug("Corrected overscan for amplifier %s.", amp.getName()) 

1395 if overscanResults is not None and \ 

1396 self.config.qa is not None and self.config.qa.saveStats is True: 

1397 if isinstance(overscanResults.overscanFit, float): 

1398 qaMedian = overscanResults.overscanFit 

1399 qaStdev = float("NaN") 

1400 else: 

1401 qaStats = afwMath.makeStatistics(overscanResults.overscanFit, 

1402 afwMath.MEDIAN | afwMath.STDEVCLIP) 

1403 qaMedian = qaStats.getValue(afwMath.MEDIAN) 

1404 qaStdev = qaStats.getValue(afwMath.STDEVCLIP) 

1405 

1406 self.metadata.set(f"FIT MEDIAN {amp.getName()}", qaMedian) 

1407 self.metadata.set(f"FIT STDEV {amp.getName()}", qaStdev) 

1408 self.log.debug(" Overscan stats for amplifer %s: %f +/- %f", 

1409 amp.getName(), qaMedian, qaStdev) 

1410 

1411 # Residuals after overscan correction 

1412 qaStatsAfter = afwMath.makeStatistics(overscanResults.overscanImage, 

1413 afwMath.MEDIAN | afwMath.STDEVCLIP) 

1414 qaMedianAfter = qaStatsAfter.getValue(afwMath.MEDIAN) 

1415 qaStdevAfter = qaStatsAfter.getValue(afwMath.STDEVCLIP) 

1416 

1417 self.metadata.set(f"RESIDUAL MEDIAN {amp.getName()}", qaMedianAfter) 

1418 self.metadata.set(f"RESIDUAL STDEV {amp.getName()}", qaStdevAfter) 

1419 self.log.debug(" Overscan stats for amplifer %s after correction: %f +/- %f", 

1420 amp.getName(), qaMedianAfter, qaStdevAfter) 

1421 

1422 ccdExposure.getMetadata().set('OVERSCAN', "Overscan corrected") 

1423 else: 

1424 if badAmp: 

1425 self.log.warn("Amplifier %s is bad.", amp.getName()) 

1426 overscanResults = None 

1427 

1428 overscans.append(overscanResults if overscanResults is not None else None) 

1429 else: 

1430 self.log.info("Skipped OSCAN for %s.", amp.getName()) 

1431 

1432 if self.config.doCrosstalk and self.config.doCrosstalkBeforeAssemble: 

1433 self.log.info("Applying crosstalk correction.") 

1434 self.crosstalk.run(ccdExposure, crosstalk=crosstalk, 

1435 crosstalkSources=crosstalkSources, camera=camera) 

1436 self.debugView(ccdExposure, "doCrosstalk") 

1437 

1438 if self.config.doAssembleCcd: 

1439 self.log.info("Assembling CCD from amplifiers.") 

1440 ccdExposure = self.assembleCcd.assembleCcd(ccdExposure) 

1441 

1442 if self.config.expectWcs and not ccdExposure.getWcs(): 

1443 self.log.warn("No WCS found in input exposure.") 

1444 self.debugView(ccdExposure, "doAssembleCcd") 

1445 

1446 ossThumb = None 

1447 if self.config.qa.doThumbnailOss: 

1448 ossThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa) 

1449 

1450 if self.config.doBias and not self.config.doBiasBeforeOverscan: 

1451 self.log.info("Applying bias correction.") 

1452 isrFunctions.biasCorrection(ccdExposure.getMaskedImage(), bias.getMaskedImage(), 

1453 trimToFit=self.config.doTrimToMatchCalib) 

1454 self.debugView(ccdExposure, "doBias") 

1455 

1456 if self.config.doVariance: 

1457 for amp, overscanResults in zip(ccd, overscans): 

1458 if ccdExposure.getBBox().contains(amp.getBBox()): 

1459 self.log.debug("Constructing variance map for amplifer %s.", amp.getName()) 

1460 ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox()) 

1461 if overscanResults is not None: 

1462 self.updateVariance(ampExposure, amp, 

1463 overscanImage=overscanResults.overscanImage) 

1464 else: 

1465 self.updateVariance(ampExposure, amp, 

1466 overscanImage=None) 

1467 if self.config.qa is not None and self.config.qa.saveStats is True: 

1468 qaStats = afwMath.makeStatistics(ampExposure.getVariance(), 

1469 afwMath.MEDIAN | afwMath.STDEVCLIP) 

1470 self.metadata.set(f"ISR VARIANCE {amp.getName()} MEDIAN", 

1471 qaStats.getValue(afwMath.MEDIAN)) 

1472 self.metadata.set(f"ISR VARIANCE {amp.getName()} STDEV", 

1473 qaStats.getValue(afwMath.STDEVCLIP)) 

1474 self.log.debug(" Variance stats for amplifer %s: %f +/- %f.", 

1475 amp.getName(), qaStats.getValue(afwMath.MEDIAN), 

1476 qaStats.getValue(afwMath.STDEVCLIP)) 

1477 

1478 if self.doLinearize(ccd): 

1479 self.log.info("Applying linearizer.") 

1480 linearizer.applyLinearity(image=ccdExposure.getMaskedImage().getImage(), 

1481 detector=ccd, log=self.log) 

1482 

1483 if self.config.doCrosstalk and not self.config.doCrosstalkBeforeAssemble: 

1484 self.log.info("Applying crosstalk correction.") 

1485 self.crosstalk.run(ccdExposure, crosstalk=crosstalk, 

1486 crosstalkSources=crosstalkSources, isTrimmed=True) 

1487 self.debugView(ccdExposure, "doCrosstalk") 

1488 

1489 # Masking block. Optionally mask known defects, NAN/inf pixels, widen trails, and do 

1490 # anything else the camera needs. Saturated and suspect pixels have already been masked. 

1491 if self.config.doDefect: 

1492 self.log.info("Masking defects.") 

1493 self.maskDefect(ccdExposure, defects) 

1494 

1495 if self.config.numEdgeSuspect > 0: 

1496 self.log.info("Masking edges as SUSPECT.") 

1497 self.maskEdges(ccdExposure, numEdgePixels=self.config.numEdgeSuspect, 

1498 maskPlane="SUSPECT", level=self.config.edgeMaskLevel) 

1499 

1500 if self.config.doNanMasking: 

1501 self.log.info("Masking non-finite (NAN, inf) value pixels.") 

1502 self.maskNan(ccdExposure) 

1503 

1504 if self.config.doWidenSaturationTrails: 

1505 self.log.info("Widening saturation trails.") 

1506 isrFunctions.widenSaturationTrails(ccdExposure.getMaskedImage().getMask()) 

1507 

1508 if self.config.doCameraSpecificMasking: 

1509 self.log.info("Masking regions for camera specific reasons.") 

1510 self.masking.run(ccdExposure) 

1511 

1512 if self.config.doBrighterFatter: 

1513 # We need to apply flats and darks before we can interpolate, and we 

1514 # need to interpolate before we do B-F, but we do B-F without the 

1515 # flats and darks applied so we can work in units of electrons or holes. 

1516 # This context manager applies and then removes the darks and flats. 

1517 # 

1518 # We also do not want to interpolate values here, so operate on temporary 

1519 # images so we can apply only the BF-correction and roll back the 

1520 # interpolation. 

1521 interpExp = ccdExposure.clone() 

1522 with self.flatContext(interpExp, flat, dark): 

1523 isrFunctions.interpolateFromMask( 

1524 maskedImage=interpExp.getMaskedImage(), 

1525 fwhm=self.config.fwhm, 

1526 growSaturatedFootprints=self.config.growSaturationFootprintSize, 

1527 maskNameList=list(self.config.brighterFatterMaskListToInterpolate) 

1528 ) 

1529 bfExp = interpExp.clone() 

1530 

1531 self.log.info("Applying brighter-fatter correction using kernel type %s / gains %s.", 

1532 type(bfKernel), type(bfGains)) 

1533 bfResults = isrFunctions.brighterFatterCorrection(bfExp, bfKernel, 

1534 self.config.brighterFatterMaxIter, 

1535 self.config.brighterFatterThreshold, 

1536 self.config.brighterFatterApplyGain, 

1537 bfGains) 

1538 if bfResults[1] == self.config.brighterFatterMaxIter: 

1539 self.log.warn("Brighter-fatter correction did not converge, final difference %f.", 

1540 bfResults[0]) 

1541 else: 

1542 self.log.info("Finished brighter-fatter correction in %d iterations.", 

1543 bfResults[1]) 

1544 image = ccdExposure.getMaskedImage().getImage() 

1545 bfCorr = bfExp.getMaskedImage().getImage() 

1546 bfCorr -= interpExp.getMaskedImage().getImage() 

1547 image += bfCorr 

1548 

1549 # Applying the brighter-fatter correction applies a 

1550 # convolution to the science image. At the edges this 

1551 # convolution may not have sufficient valid pixels to 

1552 # produce a valid correction. Mark pixels within the size 

1553 # of the brighter-fatter kernel as EDGE to warn of this 

1554 # fact. 

1555 self.log.info("Ensuring image edges are masked as EDGE to the brighter-fatter kernel size.") 

1556 self.maskEdges(ccdExposure, numEdgePixels=numpy.max(bfKernel.shape) // 2, 

1557 maskPlane="EDGE") 

1558 

1559 if self.config.brighterFatterMaskGrowSize > 0: 

1560 self.log.info("Growing masks to account for brighter-fatter kernel convolution.") 

1561 for maskPlane in self.config.brighterFatterMaskListToInterpolate: 

1562 isrFunctions.growMasks(ccdExposure.getMask(), 

1563 radius=self.config.brighterFatterMaskGrowSize, 

1564 maskNameList=maskPlane, 

1565 maskValue=maskPlane) 

1566 

1567 self.debugView(ccdExposure, "doBrighterFatter") 

1568 

1569 if self.config.doDark: 

1570 self.log.info("Applying dark correction.") 

1571 self.darkCorrection(ccdExposure, dark) 

1572 self.debugView(ccdExposure, "doDark") 

1573 

1574 if self.config.doFringe and not self.config.fringeAfterFlat: 

1575 self.log.info("Applying fringe correction before flat.") 

1576 self.fringe.run(ccdExposure, **fringes.getDict()) 

1577 self.debugView(ccdExposure, "doFringe") 

1578 

1579 if self.config.doStrayLight and self.strayLight.check(ccdExposure): 

1580 self.log.info("Checking strayLight correction.") 

1581 self.strayLight.run(ccdExposure, strayLightData) 

1582 self.debugView(ccdExposure, "doStrayLight") 

1583 

1584 if self.config.doFlat: 

1585 self.log.info("Applying flat correction.") 

1586 self.flatCorrection(ccdExposure, flat) 

1587 self.debugView(ccdExposure, "doFlat") 

1588 

1589 if self.config.doApplyGains: 

1590 self.log.info("Applying gain correction instead of flat.") 

1591 isrFunctions.applyGains(ccdExposure, self.config.normalizeGains) 

1592 

1593 if self.config.doFringe and self.config.fringeAfterFlat: 

1594 self.log.info("Applying fringe correction after flat.") 

1595 self.fringe.run(ccdExposure, **fringes.getDict()) 

1596 

1597 if self.config.doVignette: 

1598 self.log.info("Constructing Vignette polygon.") 

1599 self.vignettePolygon = self.vignette.run(ccdExposure) 

1600 

1601 if self.config.vignette.doWriteVignettePolygon: 

1602 self.setValidPolygonIntersect(ccdExposure, self.vignettePolygon) 

1603 

1604 if self.config.doAttachTransmissionCurve: 

1605 self.log.info("Adding transmission curves.") 

1606 isrFunctions.attachTransmissionCurve(ccdExposure, opticsTransmission=opticsTransmission, 

1607 filterTransmission=filterTransmission, 

1608 sensorTransmission=sensorTransmission, 

1609 atmosphereTransmission=atmosphereTransmission) 

1610 

1611 flattenedThumb = None 

1612 if self.config.qa.doThumbnailFlattened: 

1613 flattenedThumb = isrQa.makeThumbnail(ccdExposure, isrQaConfig=self.config.qa) 

1614 

1615 if self.config.doIlluminationCorrection and filterLabel in self.config.illumFilters: 

1616 self.log.info("Performing illumination correction.") 

1617 isrFunctions.illuminationCorrection(ccdExposure.getMaskedImage(), 

1618 illumMaskedImage, illumScale=self.config.illumScale, 

1619 trimToFit=self.config.doTrimToMatchCalib) 

1620 

1621 preInterpExp = None 

1622 if self.config.doSaveInterpPixels: 

1623 preInterpExp = ccdExposure.clone() 

1624 

1625 # Reset and interpolate bad pixels. 

1626 # 

1627 # Large contiguous bad regions (which should have the BAD mask 

1628 # bit set) should have their values set to the image median. 

1629 # This group should include defects and bad amplifiers. As the 

1630 # area covered by these defects are large, there's little 

1631 # reason to expect that interpolation would provide a more 

1632 # useful value. 

1633 # 

1634 # Smaller defects can be safely interpolated after the larger 

1635 # regions have had their pixel values reset. This ensures 

1636 # that the remaining defects adjacent to bad amplifiers (as an 

1637 # example) do not attempt to interpolate extreme values. 

1638 if self.config.doSetBadRegions: 

1639 badPixelCount, badPixelValue = isrFunctions.setBadRegions(ccdExposure) 

1640 if badPixelCount > 0: 

1641 self.log.info("Set %d BAD pixels to %f.", badPixelCount, badPixelValue) 

1642 

1643 if self.config.doInterpolate: 

1644 self.log.info("Interpolating masked pixels.") 

1645 isrFunctions.interpolateFromMask( 

1646 maskedImage=ccdExposure.getMaskedImage(), 

1647 fwhm=self.config.fwhm, 

1648 growSaturatedFootprints=self.config.growSaturationFootprintSize, 

1649 maskNameList=list(self.config.maskListToInterpolate) 

1650 ) 

1651 

1652 self.roughZeroPoint(ccdExposure) 

1653 

1654 if self.config.doMeasureBackground: 

1655 self.log.info("Measuring background level.") 

1656 self.measureBackground(ccdExposure, self.config.qa) 

1657 

1658 if self.config.qa is not None and self.config.qa.saveStats is True: 

1659 for amp in ccd: 

1660 ampExposure = ccdExposure.Factory(ccdExposure, amp.getBBox()) 

1661 qaStats = afwMath.makeStatistics(ampExposure.getImage(), 

1662 afwMath.MEDIAN | afwMath.STDEVCLIP) 

1663 self.metadata.set("ISR BACKGROUND {} MEDIAN".format(amp.getName()), 

1664 qaStats.getValue(afwMath.MEDIAN)) 

1665 self.metadata.set("ISR BACKGROUND {} STDEV".format(amp.getName()), 

1666 qaStats.getValue(afwMath.STDEVCLIP)) 

1667 self.log.debug(" Background stats for amplifer %s: %f +/- %f", 

1668 amp.getName(), qaStats.getValue(afwMath.MEDIAN), 

1669 qaStats.getValue(afwMath.STDEVCLIP)) 

1670 

1671 self.debugView(ccdExposure, "postISRCCD") 

1672 

1673 return pipeBase.Struct( 

1674 exposure=ccdExposure, 

1675 ossThumb=ossThumb, 

1676 flattenedThumb=flattenedThumb, 

1677 

1678 preInterpolatedExposure=preInterpExp, 

1679 outputExposure=ccdExposure, 

1680 outputOssThumbnail=ossThumb, 

1681 outputFlattenedThumbnail=flattenedThumb, 

1682 ) 

1683 

1684 @pipeBase.timeMethod 

1685 def runDataRef(self, sensorRef): 

1686 """Perform instrument signature removal on a ButlerDataRef of a Sensor. 

1687 

1688 This method contains the `CmdLineTask` interface to the ISR 

1689 processing. All IO is handled here, freeing the `run()` method 

1690 to manage only pixel-level calculations. The steps performed 

1691 are: 

1692 - Read in necessary detrending/isr/calibration data. 

1693 - Process raw exposure in `run()`. 

1694 - Persist the ISR-corrected exposure as "postISRCCD" if 

1695 config.doWrite=True. 

1696 

1697 Parameters 

1698 ---------- 

1699 sensorRef : `daf.persistence.butlerSubset.ButlerDataRef` 

1700 DataRef of the detector data to be processed 

1701 

1702 Returns 

1703 ------- 

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

1705 Result struct with component: 

1706 - ``exposure`` : `afw.image.Exposure` 

1707 The fully ISR corrected exposure. 

1708 

1709 Raises 

1710 ------ 

1711 RuntimeError 

1712 Raised if a configuration option is set to True, but the 

1713 required calibration data does not exist. 

1714 

1715 """ 

1716 self.log.info("Performing ISR on sensor %s.", sensorRef.dataId) 

1717 

1718 ccdExposure = sensorRef.get(self.config.datasetType) 

1719 

1720 camera = sensorRef.get("camera") 

1721 isrData = self.readIsrData(sensorRef, ccdExposure) 

1722 

1723 result = self.run(ccdExposure, camera=camera, **isrData.getDict()) 

1724 

1725 if self.config.doWrite: 

1726 sensorRef.put(result.exposure, "postISRCCD") 

1727 if result.preInterpolatedExposure is not None: 

1728 sensorRef.put(result.preInterpolatedExposure, "postISRCCD_uninterpolated") 

1729 if result.ossThumb is not None: 

1730 isrQa.writeThumbnail(sensorRef, result.ossThumb, "ossThumb") 

1731 if result.flattenedThumb is not None: 

1732 isrQa.writeThumbnail(sensorRef, result.flattenedThumb, "flattenedThumb") 

1733 

1734 return result 

1735 

1736 def getIsrExposure(self, dataRef, datasetType, dateObs=None, immediate=True): 

1737 """Retrieve a calibration dataset for removing instrument signature. 

1738 

1739 Parameters 

1740 ---------- 

1741 

1742 dataRef : `daf.persistence.butlerSubset.ButlerDataRef` 

1743 DataRef of the detector data to find calibration datasets 

1744 for. 

1745 datasetType : `str` 

1746 Type of dataset to retrieve (e.g. 'bias', 'flat', etc). 

1747 dateObs : `str`, optional 

1748 Date of the observation. Used to correct butler failures 

1749 when using fallback filters. 

1750 immediate : `Bool` 

1751 If True, disable butler proxies to enable error handling 

1752 within this routine. 

1753 

1754 Returns 

1755 ------- 

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

1757 Requested calibration frame. 

1758 

1759 Raises 

1760 ------ 

1761 RuntimeError 

1762 Raised if no matching calibration frame can be found. 

1763 """ 

1764 try: 

1765 exp = dataRef.get(datasetType, immediate=immediate) 

1766 except Exception as exc1: 

1767 if not self.config.fallbackFilterName: 

1768 raise RuntimeError("Unable to retrieve %s for %s: %s." % (datasetType, dataRef.dataId, exc1)) 

1769 try: 

1770 if self.config.useFallbackDate and dateObs: 

1771 exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName, 

1772 dateObs=dateObs, immediate=immediate) 

1773 else: 

1774 exp = dataRef.get(datasetType, filter=self.config.fallbackFilterName, immediate=immediate) 

1775 except Exception as exc2: 

1776 raise RuntimeError("Unable to retrieve %s for %s, even with fallback filter %s: %s AND %s." % 

1777 (datasetType, dataRef.dataId, self.config.fallbackFilterName, exc1, exc2)) 

1778 self.log.warn("Using fallback calibration from filter %s.", self.config.fallbackFilterName) 

1779 

1780 if self.config.doAssembleIsrExposures: 

1781 exp = self.assembleCcd.assembleCcd(exp) 

1782 return exp 

1783 

1784 def ensureExposure(self, inputExp, camera, detectorNum): 

1785 """Ensure that the data returned by Butler is a fully constructed exposure. 

1786 

1787 ISR requires exposure-level image data for historical reasons, so if we did 

1788 not recieve that from Butler, construct it from what we have, modifying the 

1789 input in place. 

1790 

1791 Parameters 

1792 ---------- 

1793 inputExp : `lsst.afw.image.Exposure`, `lsst.afw.image.DecoratedImageU`, or 

1794 `lsst.afw.image.ImageF` 

1795 The input data structure obtained from Butler. 

1796 camera : `lsst.afw.cameraGeom.camera` 

1797 The camera associated with the image. Used to find the appropriate 

1798 detector. 

1799 detectorNum : `int` 

1800 The detector this exposure should match. 

1801 

1802 Returns 

1803 ------- 

1804 inputExp : `lsst.afw.image.Exposure` 

1805 The re-constructed exposure, with appropriate detector parameters. 

1806 

1807 Raises 

1808 ------ 

1809 TypeError 

1810 Raised if the input data cannot be used to construct an exposure. 

1811 """ 

1812 if isinstance(inputExp, afwImage.DecoratedImageU): 

1813 inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp)) 

1814 elif isinstance(inputExp, afwImage.ImageF): 

1815 inputExp = afwImage.makeExposure(afwImage.makeMaskedImage(inputExp)) 

1816 elif isinstance(inputExp, afwImage.MaskedImageF): 

1817 inputExp = afwImage.makeExposure(inputExp) 

1818 elif isinstance(inputExp, afwImage.Exposure): 

1819 pass 

1820 elif inputExp is None: 

1821 # Assume this will be caught by the setup if it is a problem. 

1822 return inputExp 

1823 else: 

1824 raise TypeError("Input Exposure is not known type in isrTask.ensureExposure: %s." % 

1825 (type(inputExp), )) 

1826 

1827 if inputExp.getDetector() is None: 

1828 inputExp.setDetector(camera[detectorNum]) 

1829 

1830 return inputExp 

1831 

1832 def convertIntToFloat(self, exposure): 

1833 """Convert exposure image from uint16 to float. 

1834 

1835 If the exposure does not need to be converted, the input is 

1836 immediately returned. For exposures that are converted to use 

1837 floating point pixels, the variance is set to unity and the 

1838 mask to zero. 

1839 

1840 Parameters 

1841 ---------- 

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

1843 The raw exposure to be converted. 

1844 

1845 Returns 

1846 ------- 

1847 newexposure : `lsst.afw.image.Exposure` 

1848 The input ``exposure``, converted to floating point pixels. 

1849 

1850 Raises 

1851 ------ 

1852 RuntimeError 

1853 Raised if the exposure type cannot be converted to float. 

1854 

1855 """ 

1856 if isinstance(exposure, afwImage.ExposureF): 

1857 # Nothing to be done 

1858 self.log.debug("Exposure already of type float.") 

1859 return exposure 

1860 if not hasattr(exposure, "convertF"): 

1861 raise RuntimeError("Unable to convert exposure (%s) to float." % type(exposure)) 

1862 

1863 newexposure = exposure.convertF() 

1864 newexposure.variance[:] = 1 

1865 newexposure.mask[:] = 0x0 

1866 

1867 return newexposure 

1868 

1869 def maskAmplifier(self, ccdExposure, amp, defects): 

1870 """Identify bad amplifiers, saturated and suspect pixels. 

1871 

1872 Parameters 

1873 ---------- 

1874 ccdExposure : `lsst.afw.image.Exposure` 

1875 Input exposure to be masked. 

1876 amp : `lsst.afw.table.AmpInfoCatalog` 

1877 Catalog of parameters defining the amplifier on this 

1878 exposure to mask. 

1879 defects : `lsst.ip.isr.Defects` 

1880 List of defects. Used to determine if the entire 

1881 amplifier is bad. 

1882 

1883 Returns 

1884 ------- 

1885 badAmp : `Bool` 

1886 If this is true, the entire amplifier area is covered by 

1887 defects and unusable. 

1888 

1889 """ 

1890 maskedImage = ccdExposure.getMaskedImage() 

1891 

1892 badAmp = False 

1893 

1894 # Check if entire amp region is defined as a defect (need to use amp.getBBox() for correct 

1895 # comparison with current defects definition. 

1896 if defects is not None: 

1897 badAmp = bool(sum([v.getBBox().contains(amp.getBBox()) for v in defects])) 

1898 

1899 # In the case of a bad amp, we will set mask to "BAD" (here use amp.getRawBBox() for correct 

1900 # association with pixels in current ccdExposure). 

1901 if badAmp: 

1902 dataView = afwImage.MaskedImageF(maskedImage, amp.getRawBBox(), 

1903 afwImage.PARENT) 

1904 maskView = dataView.getMask() 

1905 maskView |= maskView.getPlaneBitMask("BAD") 

1906 del maskView 

1907 return badAmp 

1908 

1909 # Mask remaining defects after assembleCcd() to allow for defects that cross amplifier boundaries. 

1910 # Saturation and suspect pixels can be masked now, though. 

1911 limits = dict() 

1912 if self.config.doSaturation and not badAmp: 

1913 limits.update({self.config.saturatedMaskName: amp.getSaturation()}) 

1914 if self.config.doSuspect and not badAmp: 

1915 limits.update({self.config.suspectMaskName: amp.getSuspectLevel()}) 

1916 if math.isfinite(self.config.saturation): 

1917 limits.update({self.config.saturatedMaskName: self.config.saturation}) 

1918 

1919 for maskName, maskThreshold in limits.items(): 

1920 if not math.isnan(maskThreshold): 

1921 dataView = maskedImage.Factory(maskedImage, amp.getRawBBox()) 

1922 isrFunctions.makeThresholdMask( 

1923 maskedImage=dataView, 

1924 threshold=maskThreshold, 

1925 growFootprints=0, 

1926 maskName=maskName 

1927 ) 

1928 

1929 # Determine if we've fully masked this amplifier with SUSPECT and SAT pixels. 

1930 maskView = afwImage.Mask(maskedImage.getMask(), amp.getRawDataBBox(), 

1931 afwImage.PARENT) 

1932 maskVal = maskView.getPlaneBitMask([self.config.saturatedMaskName, 

1933 self.config.suspectMaskName]) 

1934 if numpy.all(maskView.getArray() & maskVal > 0): 

1935 badAmp = True 

1936 maskView |= maskView.getPlaneBitMask("BAD") 

1937 

1938 return badAmp 

1939 

1940 def overscanCorrection(self, ccdExposure, amp): 

1941 """Apply overscan correction in place. 

1942 

1943 This method does initial pixel rejection of the overscan 

1944 region. The overscan can also be optionally segmented to 

1945 allow for discontinuous overscan responses to be fit 

1946 separately. The actual overscan subtraction is performed by 

1947 the `lsst.ip.isr.isrFunctions.overscanCorrection` function, 

1948 which is called here after the amplifier is preprocessed. 

1949 

1950 Parameters 

1951 ---------- 

1952 ccdExposure : `lsst.afw.image.Exposure` 

1953 Exposure to have overscan correction performed. 

1954 amp : `lsst.afw.cameraGeom.Amplifer` 

1955 The amplifier to consider while correcting the overscan. 

1956 

1957 Returns 

1958 ------- 

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

1960 Result struct with components: 

1961 - ``imageFit`` : scalar or `lsst.afw.image.Image` 

1962 Value or fit subtracted from the amplifier image data. 

1963 - ``overscanFit`` : scalar or `lsst.afw.image.Image` 

1964 Value or fit subtracted from the overscan image data. 

1965 - ``overscanImage`` : `lsst.afw.image.Image` 

1966 Image of the overscan region with the overscan 

1967 correction applied. This quantity is used to estimate 

1968 the amplifier read noise empirically. 

1969 

1970 Raises 

1971 ------ 

1972 RuntimeError 

1973 Raised if the ``amp`` does not contain raw pixel information. 

1974 

1975 See Also 

1976 -------- 

1977 lsst.ip.isr.isrFunctions.overscanCorrection 

1978 """ 

1979 if amp.getRawHorizontalOverscanBBox().isEmpty(): 

1980 self.log.info("ISR_OSCAN: No overscan region. Not performing overscan correction.") 

1981 return None 

1982 

1983 statControl = afwMath.StatisticsControl() 

1984 statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT")) 

1985 

1986 # Determine the bounding boxes 

1987 dataBBox = amp.getRawDataBBox() 

1988 oscanBBox = amp.getRawHorizontalOverscanBBox() 

1989 dx0 = 0 

1990 dx1 = 0 

1991 

1992 prescanBBox = amp.getRawPrescanBBox() 

1993 if (oscanBBox.getBeginX() > prescanBBox.getBeginX()): # amp is at the right 

1994 dx0 += self.config.overscanNumLeadingColumnsToSkip 

1995 dx1 -= self.config.overscanNumTrailingColumnsToSkip 

1996 else: 

1997 dx0 += self.config.overscanNumTrailingColumnsToSkip 

1998 dx1 -= self.config.overscanNumLeadingColumnsToSkip 

1999 

2000 # Determine if we need to work on subregions of the amplifier and overscan. 

2001 imageBBoxes = [] 

2002 overscanBBoxes = [] 

2003 

2004 if ((self.config.overscanBiasJump 

2005 and self.config.overscanBiasJumpLocation) 

2006 and (ccdExposure.getMetadata().exists(self.config.overscanBiasJumpKeyword) 

2007 and ccdExposure.getMetadata().getScalar(self.config.overscanBiasJumpKeyword) in 

2008 self.config.overscanBiasJumpDevices)): 

2009 if amp.getReadoutCorner() in (ReadoutCorner.LL, ReadoutCorner.LR): 

2010 yLower = self.config.overscanBiasJumpLocation 

2011 yUpper = dataBBox.getHeight() - yLower 

2012 else: 

2013 yUpper = self.config.overscanBiasJumpLocation 

2014 yLower = dataBBox.getHeight() - yUpper 

2015 

2016 imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(), 

2017 lsst.geom.Extent2I(dataBBox.getWidth(), yLower))) 

2018 overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0), 

2019 lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1, 

2020 yLower))) 

2021 

2022 imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin() + lsst.geom.Extent2I(0, yLower), 

2023 lsst.geom.Extent2I(dataBBox.getWidth(), yUpper))) 

2024 overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, yLower), 

2025 lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1, 

2026 yUpper))) 

2027 else: 

2028 imageBBoxes.append(lsst.geom.Box2I(dataBBox.getBegin(), 

2029 lsst.geom.Extent2I(dataBBox.getWidth(), dataBBox.getHeight()))) 

2030 overscanBBoxes.append(lsst.geom.Box2I(oscanBBox.getBegin() + lsst.geom.Extent2I(dx0, 0), 

2031 lsst.geom.Extent2I(oscanBBox.getWidth() - dx0 + dx1, 

2032 oscanBBox.getHeight()))) 

2033 

2034 # Perform overscan correction on subregions, ensuring saturated pixels are masked. 

2035 for imageBBox, overscanBBox in zip(imageBBoxes, overscanBBoxes): 

2036 ampImage = ccdExposure.maskedImage[imageBBox] 

2037 overscanImage = ccdExposure.maskedImage[overscanBBox] 

2038 

2039 overscanArray = overscanImage.image.array 

2040 median = numpy.ma.median(numpy.ma.masked_where(overscanImage.mask.array, overscanArray)) 

2041 bad = numpy.where(numpy.abs(overscanArray - median) > self.config.overscanMaxDev) 

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

2043 

2044 statControl = afwMath.StatisticsControl() 

2045 statControl.setAndMask(ccdExposure.mask.getPlaneBitMask("SAT")) 

2046 

2047 overscanResults = self.overscan.run(ampImage.getImage(), overscanImage, amp) 

2048 

2049 # Measure average overscan levels and record them in the metadata. 

2050 levelStat = afwMath.MEDIAN 

2051 sigmaStat = afwMath.STDEVCLIP 

2052 

2053 sctrl = afwMath.StatisticsControl(self.config.qa.flatness.clipSigma, 

2054 self.config.qa.flatness.nIter) 

2055 metadata = ccdExposure.getMetadata() 

2056 ampNum = amp.getName() 

2057 # if self.config.overscanFitType in ("MEDIAN", "MEAN", "MEANCLIP"): 

2058 if isinstance(overscanResults.overscanFit, float): 

2059 metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, overscanResults.overscanFit) 

2060 metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, 0.0) 

2061 else: 

2062 stats = afwMath.makeStatistics(overscanResults.overscanFit, levelStat | sigmaStat, sctrl) 

2063 metadata.set("ISR_OSCAN_LEVEL%s" % ampNum, stats.getValue(levelStat)) 

2064 metadata.set("ISR_OSCAN_SIGMA%s" % ampNum, stats.getValue(sigmaStat)) 

2065 

2066 return overscanResults 

2067 

2068 def updateVariance(self, ampExposure, amp, overscanImage=None): 

2069 """Set the variance plane using the amplifier gain and read noise 

2070 

2071 The read noise is calculated from the ``overscanImage`` if the 

2072 ``doEmpiricalReadNoise`` option is set in the configuration; otherwise 

2073 the value from the amplifier data is used. 

2074 

2075 Parameters 

2076 ---------- 

2077 ampExposure : `lsst.afw.image.Exposure` 

2078 Exposure to process. 

2079 amp : `lsst.afw.table.AmpInfoRecord` or `FakeAmp` 

2080 Amplifier detector data. 

2081 overscanImage : `lsst.afw.image.MaskedImage`, optional. 

2082 Image of overscan, required only for empirical read noise. 

2083 

2084 See also 

2085 -------- 

2086 lsst.ip.isr.isrFunctions.updateVariance 

2087 """ 

2088 maskPlanes = [self.config.saturatedMaskName, self.config.suspectMaskName] 

2089 gain = amp.getGain() 

2090 

2091 if math.isnan(gain): 

2092 gain = 1.0 

2093 self.log.warn("Gain set to NAN! Updating to 1.0 to generate Poisson variance.") 

2094 elif gain <= 0: 

2095 patchedGain = 1.0 

2096 self.log.warn("Gain for amp %s == %g <= 0; setting to %f.", 

2097 amp.getName(), gain, patchedGain) 

2098 gain = patchedGain 

2099 

2100 if self.config.doEmpiricalReadNoise and overscanImage is None: 

2101 self.log.info("Overscan is none for EmpiricalReadNoise.") 

2102 

2103 if self.config.doEmpiricalReadNoise and overscanImage is not None: 

2104 stats = afwMath.StatisticsControl() 

2105 stats.setAndMask(overscanImage.mask.getPlaneBitMask(maskPlanes)) 

2106 readNoise = afwMath.makeStatistics(overscanImage, afwMath.STDEVCLIP, stats).getValue() 

2107 self.log.info("Calculated empirical read noise for amp %s: %f.", 

2108 amp.getName(), readNoise) 

2109 else: 

2110 readNoise = amp.getReadNoise() 

2111 

2112 isrFunctions.updateVariance( 

2113 maskedImage=ampExposure.getMaskedImage(), 

2114 gain=gain, 

2115 readNoise=readNoise, 

2116 ) 

2117 

2118 def darkCorrection(self, exposure, darkExposure, invert=False): 

2119 """Apply dark correction in place. 

2120 

2121 Parameters 

2122 ---------- 

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

2124 Exposure to process. 

2125 darkExposure : `lsst.afw.image.Exposure` 

2126 Dark exposure of the same size as ``exposure``. 

2127 invert : `Bool`, optional 

2128 If True, re-add the dark to an already corrected image. 

2129 

2130 Raises 

2131 ------ 

2132 RuntimeError 

2133 Raised if either ``exposure`` or ``darkExposure`` do not 

2134 have their dark time defined. 

2135 

2136 See Also 

2137 -------- 

2138 lsst.ip.isr.isrFunctions.darkCorrection 

2139 """ 

2140 expScale = exposure.getInfo().getVisitInfo().getDarkTime() 

2141 if math.isnan(expScale): 

2142 raise RuntimeError("Exposure darktime is NAN.") 

2143 if darkExposure.getInfo().getVisitInfo() is not None \ 

2144 and not math.isnan(darkExposure.getInfo().getVisitInfo().getDarkTime()): 

2145 darkScale = darkExposure.getInfo().getVisitInfo().getDarkTime() 

2146 else: 

2147 # DM-17444: darkExposure.getInfo.getVisitInfo() is None 

2148 # so getDarkTime() does not exist. 

2149 self.log.warn("darkExposure.getInfo().getVisitInfo() does not exist. Using darkScale = 1.0.") 

2150 darkScale = 1.0 

2151 

2152 isrFunctions.darkCorrection( 

2153 maskedImage=exposure.getMaskedImage(), 

2154 darkMaskedImage=darkExposure.getMaskedImage(), 

2155 expScale=expScale, 

2156 darkScale=darkScale, 

2157 invert=invert, 

2158 trimToFit=self.config.doTrimToMatchCalib 

2159 ) 

2160 

2161 def doLinearize(self, detector): 

2162 """Check if linearization is needed for the detector cameraGeom. 

2163 

2164 Checks config.doLinearize and the linearity type of the first 

2165 amplifier. 

2166 

2167 Parameters 

2168 ---------- 

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

2170 Detector to get linearity type from. 

2171 

2172 Returns 

2173 ------- 

2174 doLinearize : `Bool` 

2175 If True, linearization should be performed. 

2176 """ 

2177 return self.config.doLinearize and \ 

2178 detector.getAmplifiers()[0].getLinearityType() != NullLinearityType 

2179 

2180 def flatCorrection(self, exposure, flatExposure, invert=False): 

2181 """Apply flat correction in place. 

2182 

2183 Parameters 

2184 ---------- 

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

2186 Exposure to process. 

2187 flatExposure : `lsst.afw.image.Exposure` 

2188 Flat exposure of the same size as ``exposure``. 

2189 invert : `Bool`, optional 

2190 If True, unflatten an already flattened image. 

2191 

2192 See Also 

2193 -------- 

2194 lsst.ip.isr.isrFunctions.flatCorrection 

2195 """ 

2196 isrFunctions.flatCorrection( 

2197 maskedImage=exposure.getMaskedImage(), 

2198 flatMaskedImage=flatExposure.getMaskedImage(), 

2199 scalingType=self.config.flatScalingType, 

2200 userScale=self.config.flatUserScale, 

2201 invert=invert, 

2202 trimToFit=self.config.doTrimToMatchCalib 

2203 ) 

2204 

2205 def saturationDetection(self, exposure, amp): 

2206 """Detect saturated pixels and mask them using mask plane config.saturatedMaskName, in place. 

2207 

2208 Parameters 

2209 ---------- 

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

2211 Exposure to process. Only the amplifier DataSec is processed. 

2212 amp : `lsst.afw.table.AmpInfoCatalog` 

2213 Amplifier detector data. 

2214 

2215 See Also 

2216 -------- 

2217 lsst.ip.isr.isrFunctions.makeThresholdMask 

2218 """ 

2219 if not math.isnan(amp.getSaturation()): 

2220 maskedImage = exposure.getMaskedImage() 

2221 dataView = maskedImage.Factory(maskedImage, amp.getRawBBox()) 

2222 isrFunctions.makeThresholdMask( 

2223 maskedImage=dataView, 

2224 threshold=amp.getSaturation(), 

2225 growFootprints=0, 

2226 maskName=self.config.saturatedMaskName, 

2227 ) 

2228 

2229 def saturationInterpolation(self, exposure): 

2230 """Interpolate over saturated pixels, in place. 

2231 

2232 This method should be called after `saturationDetection`, to 

2233 ensure that the saturated pixels have been identified in the 

2234 SAT mask. It should also be called after `assembleCcd`, since 

2235 saturated regions may cross amplifier boundaries. 

2236 

2237 Parameters 

2238 ---------- 

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

2240 Exposure to process. 

2241 

2242 See Also 

2243 -------- 

2244 lsst.ip.isr.isrTask.saturationDetection 

2245 lsst.ip.isr.isrFunctions.interpolateFromMask 

2246 """ 

2247 isrFunctions.interpolateFromMask( 

2248 maskedImage=exposure.getMaskedImage(), 

2249 fwhm=self.config.fwhm, 

2250 growSaturatedFootprints=self.config.growSaturationFootprintSize, 

2251 maskNameList=list(self.config.saturatedMaskName), 

2252 ) 

2253 

2254 def suspectDetection(self, exposure, amp): 

2255 """Detect suspect pixels and mask them using mask plane config.suspectMaskName, in place. 

2256 

2257 Parameters 

2258 ---------- 

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

2260 Exposure to process. Only the amplifier DataSec is processed. 

2261 amp : `lsst.afw.table.AmpInfoCatalog` 

2262 Amplifier detector data. 

2263 

2264 See Also 

2265 -------- 

2266 lsst.ip.isr.isrFunctions.makeThresholdMask 

2267 

2268 Notes 

2269 ----- 

2270 Suspect pixels are pixels whose value is greater than amp.getSuspectLevel(). 

2271 This is intended to indicate pixels that may be affected by unknown systematics; 

2272 for example if non-linearity corrections above a certain level are unstable 

2273 then that would be a useful value for suspectLevel. A value of `nan` indicates 

2274 that no such level exists and no pixels are to be masked as suspicious. 

2275 """ 

2276 suspectLevel = amp.getSuspectLevel() 

2277 if math.isnan(suspectLevel): 

2278 return 

2279 

2280 maskedImage = exposure.getMaskedImage() 

2281 dataView = maskedImage.Factory(maskedImage, amp.getRawBBox()) 

2282 isrFunctions.makeThresholdMask( 

2283 maskedImage=dataView, 

2284 threshold=suspectLevel, 

2285 growFootprints=0, 

2286 maskName=self.config.suspectMaskName, 

2287 ) 

2288 

2289 def maskDefect(self, exposure, defectBaseList): 

2290 """Mask defects using mask plane "BAD", in place. 

2291 

2292 Parameters 

2293 ---------- 

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

2295 Exposure to process. 

2296 defectBaseList : `lsst.ip.isr.Defects` or `list` of 

2297 `lsst.afw.image.DefectBase`. 

2298 List of defects to mask. 

2299 

2300 Notes 

2301 ----- 

2302 Call this after CCD assembly, since defects may cross amplifier boundaries. 

2303 """ 

2304 maskedImage = exposure.getMaskedImage() 

2305 if not isinstance(defectBaseList, Defects): 

2306 # Promotes DefectBase to Defect 

2307 defectList = Defects(defectBaseList) 

2308 else: 

2309 defectList = defectBaseList 

2310 defectList.maskPixels(maskedImage, maskName="BAD") 

2311 

2312 def maskEdges(self, exposure, numEdgePixels=0, maskPlane="SUSPECT", level='DETECTOR'): 

2313 """Mask edge pixels with applicable mask plane. 

2314 

2315 Parameters 

2316 ---------- 

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

2318 Exposure to process. 

2319 numEdgePixels : `int`, optional 

2320 Number of edge pixels to mask. 

2321 maskPlane : `str`, optional 

2322 Mask plane name to use. 

2323 level : `str`, optional 

2324 Level at which to mask edges. 

2325 """ 

2326 maskedImage = exposure.getMaskedImage() 

2327 maskBitMask = maskedImage.getMask().getPlaneBitMask(maskPlane) 

2328 

2329 if numEdgePixels > 0: 

2330 if level == 'DETECTOR': 

2331 boxes = [maskedImage.getBBox()] 

2332 elif level == 'AMP': 

2333 boxes = [amp.getBBox() for amp in exposure.getDetector()] 

2334 

2335 for box in boxes: 

2336 # This makes a bbox numEdgeSuspect pixels smaller than the image on each side 

2337 subImage = maskedImage[box] 

2338 box.grow(-numEdgePixels) 

2339 # Mask pixels outside box 

2340 SourceDetectionTask.setEdgeBits( 

2341 subImage, 

2342 box, 

2343 maskBitMask) 

2344 

2345 def maskAndInterpolateDefects(self, exposure, defectBaseList): 

2346 """Mask and interpolate defects using mask plane "BAD", in place. 

2347 

2348 Parameters 

2349 ---------- 

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

2351 Exposure to process. 

2352 defectBaseList : `lsst.ip.isr.Defects` or `list` of 

2353 `lsst.afw.image.DefectBase`. 

2354 List of defects to mask and interpolate. 

2355 

2356 See Also 

2357 -------- 

2358 lsst.ip.isr.isrTask.maskDefect 

2359 """ 

2360 self.maskDefect(exposure, defectBaseList) 

2361 self.maskEdges(exposure, numEdgePixels=self.config.numEdgeSuspect, 

2362 maskPlane="SUSPECT", level=self.config.edgeMaskLevel) 

2363 isrFunctions.interpolateFromMask( 

2364 maskedImage=exposure.getMaskedImage(), 

2365 fwhm=self.config.fwhm, 

2366 growSaturatedFootprints=0, 

2367 maskNameList=["BAD"], 

2368 ) 

2369 

2370 def maskNan(self, exposure): 

2371 """Mask NaNs using mask plane "UNMASKEDNAN", in place. 

2372 

2373 Parameters 

2374 ---------- 

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

2376 Exposure to process. 

2377 

2378 Notes 

2379 ----- 

2380 We mask over all non-finite values (NaN, inf), including those 

2381 that are masked with other bits (because those may or may not be 

2382 interpolated over later, and we want to remove all NaN/infs). 

2383 Despite this behaviour, the "UNMASKEDNAN" mask plane is used to 

2384 preserve the historical name. 

2385 """ 

2386 maskedImage = exposure.getMaskedImage() 

2387 

2388 # Find and mask NaNs 

2389 maskedImage.getMask().addMaskPlane("UNMASKEDNAN") 

2390 maskVal = maskedImage.getMask().getPlaneBitMask("UNMASKEDNAN") 

2391 numNans = maskNans(maskedImage, maskVal) 

2392 self.metadata.set("NUMNANS", numNans) 

2393 if numNans > 0: 

2394 self.log.warn("There were %d unmasked NaNs.", numNans) 

2395 

2396 def maskAndInterpolateNan(self, exposure): 

2397 """"Mask and interpolate NaN/infs using mask plane "UNMASKEDNAN", 

2398 in place. 

2399 

2400 Parameters 

2401 ---------- 

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

2403 Exposure to process. 

2404 

2405 See Also 

2406 -------- 

2407 lsst.ip.isr.isrTask.maskNan 

2408 """ 

2409 self.maskNan(exposure) 

2410 isrFunctions.interpolateFromMask( 

2411 maskedImage=exposure.getMaskedImage(), 

2412 fwhm=self.config.fwhm, 

2413 growSaturatedFootprints=0, 

2414 maskNameList=["UNMASKEDNAN"], 

2415 ) 

2416 

2417 def measureBackground(self, exposure, IsrQaConfig=None): 

2418 """Measure the image background in subgrids, for quality control purposes. 

2419 

2420 Parameters 

2421 ---------- 

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

2423 Exposure to process. 

2424 IsrQaConfig : `lsst.ip.isr.isrQa.IsrQaConfig` 

2425 Configuration object containing parameters on which background 

2426 statistics and subgrids to use. 

2427 """ 

2428 if IsrQaConfig is not None: 

2429 statsControl = afwMath.StatisticsControl(IsrQaConfig.flatness.clipSigma, 

2430 IsrQaConfig.flatness.nIter) 

2431 maskVal = exposure.getMaskedImage().getMask().getPlaneBitMask(["BAD", "SAT", "DETECTED"]) 

2432 statsControl.setAndMask(maskVal) 

2433 maskedImage = exposure.getMaskedImage() 

2434 stats = afwMath.makeStatistics(maskedImage, afwMath.MEDIAN | afwMath.STDEVCLIP, statsControl) 

2435 skyLevel = stats.getValue(afwMath.MEDIAN) 

2436 skySigma = stats.getValue(afwMath.STDEVCLIP) 

2437 self.log.info("Flattened sky level: %f +/- %f.", skyLevel, skySigma) 

2438 metadata = exposure.getMetadata() 

2439 metadata.set('SKYLEVEL', skyLevel) 

2440 metadata.set('SKYSIGMA', skySigma) 

2441 

2442 # calcluating flatlevel over the subgrids 

2443 stat = afwMath.MEANCLIP if IsrQaConfig.flatness.doClip else afwMath.MEAN 

2444 meshXHalf = int(IsrQaConfig.flatness.meshX/2.) 

2445 meshYHalf = int(IsrQaConfig.flatness.meshY/2.) 

2446 nX = int((exposure.getWidth() + meshXHalf) / IsrQaConfig.flatness.meshX) 

2447 nY = int((exposure.getHeight() + meshYHalf) / IsrQaConfig.flatness.meshY) 

2448 skyLevels = numpy.zeros((nX, nY)) 

2449 

2450 for j in range(nY): 

2451 yc = meshYHalf + j * IsrQaConfig.flatness.meshY 

2452 for i in range(nX): 

2453 xc = meshXHalf + i * IsrQaConfig.flatness.meshX 

2454 

2455 xLLC = xc - meshXHalf 

2456 yLLC = yc - meshYHalf 

2457 xURC = xc + meshXHalf - 1 

2458 yURC = yc + meshYHalf - 1 

2459 

2460 bbox = lsst.geom.Box2I(lsst.geom.Point2I(xLLC, yLLC), lsst.geom.Point2I(xURC, yURC)) 

2461 miMesh = maskedImage.Factory(exposure.getMaskedImage(), bbox, afwImage.LOCAL) 

2462 

2463 skyLevels[i, j] = afwMath.makeStatistics(miMesh, stat, statsControl).getValue() 

2464 

2465 good = numpy.where(numpy.isfinite(skyLevels)) 

2466 skyMedian = numpy.median(skyLevels[good]) 

2467 flatness = (skyLevels[good] - skyMedian) / skyMedian 

2468 flatness_rms = numpy.std(flatness) 

2469 flatness_pp = flatness.max() - flatness.min() if len(flatness) > 0 else numpy.nan 

2470 

2471 self.log.info("Measuring sky levels in %dx%d grids: %f.", nX, nY, skyMedian) 

2472 self.log.info("Sky flatness in %dx%d grids - pp: %f rms: %f.", 

2473 nX, nY, flatness_pp, flatness_rms) 

2474 

2475 metadata.set('FLATNESS_PP', float(flatness_pp)) 

2476 metadata.set('FLATNESS_RMS', float(flatness_rms)) 

2477 metadata.set('FLATNESS_NGRIDS', '%dx%d' % (nX, nY)) 

2478 metadata.set('FLATNESS_MESHX', IsrQaConfig.flatness.meshX) 

2479 metadata.set('FLATNESS_MESHY', IsrQaConfig.flatness.meshY) 

2480 

2481 def roughZeroPoint(self, exposure): 

2482 """Set an approximate magnitude zero point for the exposure. 

2483 

2484 Parameters 

2485 ---------- 

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

2487 Exposure to process. 

2488 """ 

2489 filterLabel = exposure.getFilterLabel() 

2490 if filterLabel in self.config.fluxMag0T1: 

2491 fluxMag0 = self.config.fluxMag0T1[filterLabel] 

2492 else: 

2493 self.log.warn("No rough magnitude zero point set for filter %s.", filterLabel) 

2494 fluxMag0 = self.config.defaultFluxMag0T1 

2495 

2496 expTime = exposure.getInfo().getVisitInfo().getExposureTime() 

2497 if not expTime > 0: # handle NaN as well as <= 0 

2498 self.log.warn("Non-positive exposure time; skipping rough zero point.") 

2499 return 

2500 

2501 self.log.info("Setting rough magnitude zero point: %f", 2.5*math.log10(fluxMag0*expTime)) 

2502 exposure.setPhotoCalib(afwImage.makePhotoCalibFromCalibZeroPoint(fluxMag0*expTime, 0.0)) 

2503 

2504 def setValidPolygonIntersect(self, ccdExposure, fpPolygon): 

2505 """Set the valid polygon as the intersection of fpPolygon and the ccd corners. 

2506 

2507 Parameters 

2508 ---------- 

2509 ccdExposure : `lsst.afw.image.Exposure` 

2510 Exposure to process. 

2511 fpPolygon : `lsst.afw.geom.Polygon` 

2512 Polygon in focal plane coordinates. 

2513 """ 

2514 # Get ccd corners in focal plane coordinates 

2515 ccd = ccdExposure.getDetector() 

2516 fpCorners = ccd.getCorners(FOCAL_PLANE) 

2517 ccdPolygon = Polygon(fpCorners) 

2518 

2519 # Get intersection of ccd corners with fpPolygon 

2520 intersect = ccdPolygon.intersectionSingle(fpPolygon) 

2521 

2522 # Transform back to pixel positions and build new polygon 

2523 ccdPoints = ccd.transform(intersect, FOCAL_PLANE, PIXELS) 

2524 validPolygon = Polygon(ccdPoints) 

2525 ccdExposure.getInfo().setValidPolygon(validPolygon) 

2526 

2527 @contextmanager 

2528 def flatContext(self, exp, flat, dark=None): 

2529 """Context manager that applies and removes flats and darks, 

2530 if the task is configured to apply them. 

2531 

2532 Parameters 

2533 ---------- 

2534 exp : `lsst.afw.image.Exposure` 

2535 Exposure to process. 

2536 flat : `lsst.afw.image.Exposure` 

2537 Flat exposure the same size as ``exp``. 

2538 dark : `lsst.afw.image.Exposure`, optional 

2539 Dark exposure the same size as ``exp``. 

2540 

2541 Yields 

2542 ------ 

2543 exp : `lsst.afw.image.Exposure` 

2544 The flat and dark corrected exposure. 

2545 """ 

2546 if self.config.doDark and dark is not None: 

2547 self.darkCorrection(exp, dark) 

2548 if self.config.doFlat: 

2549 self.flatCorrection(exp, flat) 

2550 try: 

2551 yield exp 

2552 finally: 

2553 if self.config.doFlat: 

2554 self.flatCorrection(exp, flat, invert=True) 

2555 if self.config.doDark and dark is not None: 

2556 self.darkCorrection(exp, dark, invert=True) 

2557 

2558 def debugView(self, exposure, stepname): 

2559 """Utility function to examine ISR exposure at different stages. 

2560 

2561 Parameters 

2562 ---------- 

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

2564 Exposure to view. 

2565 stepname : `str` 

2566 State of processing to view. 

2567 """ 

2568 frame = getDebugFrame(self._display, stepname) 

2569 if frame: 

2570 display = getDisplay(frame) 

2571 display.scale('asinh', 'zscale') 

2572 display.mtv(exposure) 

2573 prompt = "Press Enter to continue [c]... " 

2574 while True: 

2575 ans = input(prompt).lower() 

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

2577 break 

2578 

2579 

2580class FakeAmp(object): 

2581 """A Detector-like object that supports returning gain and saturation level 

2582 

2583 This is used when the input exposure does not have a detector. 

2584 

2585 Parameters 

2586 ---------- 

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

2588 Exposure to generate a fake amplifier for. 

2589 config : `lsst.ip.isr.isrTaskConfig` 

2590 Configuration to apply to the fake amplifier. 

2591 """ 

2592 

2593 def __init__(self, exposure, config): 

2594 self._bbox = exposure.getBBox(afwImage.LOCAL) 

2595 self._RawHorizontalOverscanBBox = lsst.geom.Box2I() 

2596 self._gain = config.gain 

2597 self._readNoise = config.readNoise 

2598 self._saturation = config.saturation 

2599 

2600 def getBBox(self): 

2601 return self._bbox 

2602 

2603 def getRawBBox(self): 

2604 return self._bbox 

2605 

2606 def getRawHorizontalOverscanBBox(self): 

2607 return self._RawHorizontalOverscanBBox 

2608 

2609 def getGain(self): 

2610 return self._gain 

2611 

2612 def getReadNoise(self): 

2613 return self._readNoise 

2614 

2615 def getSaturation(self): 

2616 return self._saturation 

2617 

2618 def getSuspectLevel(self): 

2619 return float("NaN") 

2620 

2621 

2622class RunIsrConfig(pexConfig.Config): 

2623 isr = pexConfig.ConfigurableField(target=IsrTask, doc="Instrument signature removal") 

2624 

2625 

2626class RunIsrTask(pipeBase.CmdLineTask): 

2627 """Task to wrap the default IsrTask to allow it to be retargeted. 

2628 

2629 The standard IsrTask can be called directly from a command line 

2630 program, but doing so removes the ability of the task to be 

2631 retargeted. As most cameras override some set of the IsrTask 

2632 methods, this would remove those data-specific methods in the 

2633 output post-ISR images. This wrapping class fixes the issue, 

2634 allowing identical post-ISR images to be generated by both the 

2635 processCcd and isrTask code. 

2636 """ 

2637 ConfigClass = RunIsrConfig 

2638 _DefaultName = "runIsr" 

2639 

2640 def __init__(self, *args, **kwargs): 

2641 super().__init__(*args, **kwargs) 

2642 self.makeSubtask("isr") 

2643 

2644 def runDataRef(self, dataRef): 

2645 """ 

2646 Parameters 

2647 ---------- 

2648 dataRef : `lsst.daf.persistence.ButlerDataRef` 

2649 data reference of the detector data to be processed 

2650 

2651 Returns 

2652 ------- 

2653 result : `pipeBase.Struct` 

2654 Result struct with component: 

2655 

2656 - exposure : `lsst.afw.image.Exposure` 

2657 Post-ISR processed exposure. 

2658 """ 

2659 return self.isr.runDataRef(dataRef)