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 cp_pipe 

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 <http://www.gnu.org/licenses/>. 

21import itertools 

22import numpy as np 

23from scipy.stats import norm 

24 

25from collections import defaultdict 

26 

27import lsst.pipe.base as pipeBase 

28import lsst.pipe.base.connectionTypes as cT 

29 

30from lsstDebug import getDebugFrame 

31from lsst.afw.detection import FootprintSet, Threshold 

32from lsst.afw.display import getDisplay 

33from lsst.pex.config import Config, Field, ListField, ConfigurableField 

34from lsst.ip.isr import CrosstalkCalib, IsrProvenance 

35from lsst.pipe.tasks.getRepositoryData import DataRefListRunner 

36 

37from ._lookupStaticCalibration import lookupStaticCalibration 

38 

39__all__ = ["CrosstalkExtractConfig", "CrosstalkExtractTask", 

40 "CrosstalkSolveTask", "CrosstalkSolveConfig", 

41 "MeasureCrosstalkConfig", "MeasureCrosstalkTask"] 

42 

43 

44class CrosstalkExtractConnections(pipeBase.PipelineTaskConnections, 

45 dimensions=("instrument", "exposure", "detector")): 

46 inputExp = cT.Input( 

47 name="crosstalkInputs", 

48 doc="Input post-ISR processed exposure to measure crosstalk from.", 

49 storageClass="Exposure", 

50 dimensions=("instrument", "exposure", "detector"), 

51 multiple=False, 

52 ) 

53 # TODO: Depends on DM-21904. 

54 sourceExp = cT.Input( 

55 name="crosstalkSource", 

56 doc="Post-ISR exposure to measure for inter-chip crosstalk onto inputExp.", 

57 storageClass="Exposure", 

58 dimensions=("instrument", "exposure", "detector"), 

59 multiple=True, 

60 deferLoad=True, 

61 # lookupFunction=None, 

62 ) 

63 

64 outputRatios = cT.Output( 

65 name="crosstalkRatios", 

66 doc="Extracted crosstalk pixel ratios.", 

67 storageClass="StructuredDataDict", 

68 dimensions=("instrument", "exposure", "detector"), 

69 ) 

70 outputFluxes = cT.Output( 

71 name="crosstalkFluxes", 

72 doc="Source pixel fluxes used in ratios.", 

73 storageClass="StructuredDataDict", 

74 dimensions=("instrument", "exposure", "detector"), 

75 ) 

76 

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

78 super().__init__(config=config) 

79 # Discard sourceExp until DM-21904 allows full interchip 

80 # measurements. 

81 self.inputs.discard("sourceExp") 

82 

83 

84class CrosstalkExtractConfig(pipeBase.PipelineTaskConfig, 

85 pipelineConnections=CrosstalkExtractConnections): 

86 """Configuration for the measurement of pixel ratios. 

87 """ 

88 doMeasureInterchip = Field( 

89 dtype=bool, 

90 default=False, 

91 doc="Measure inter-chip crosstalk as well?", 

92 ) 

93 threshold = Field( 

94 dtype=float, 

95 default=30000, 

96 doc="Minimum level of source pixels for which to measure crosstalk." 

97 ) 

98 ignoreSaturatedPixels = Field( 

99 dtype=bool, 

100 default=True, 

101 doc="Should saturated pixels be ignored?" 

102 ) 

103 badMask = ListField( 

104 dtype=str, 

105 default=["BAD", "INTRP"], 

106 doc="Mask planes to ignore when identifying source pixels." 

107 ) 

108 isTrimmed = Field( 

109 dtype=bool, 

110 default=True, 

111 doc="Is the input exposure trimmed?" 

112 ) 

113 

114 def validate(self): 

115 super().validate() 

116 

117 # Ensure the handling of the SAT mask plane is consistent 

118 # with the ignoreSaturatedPixels value. 

119 if self.ignoreSaturatedPixels: 

120 if 'SAT' not in self.badMask: 

121 self.badMask.append('SAT') 

122 else: 

123 if 'SAT' in self.badMask: 

124 self.badMask = [mask for mask in self.badMask if mask != 'SAT'] 

125 

126 

127class CrosstalkExtractTask(pipeBase.PipelineTask, 

128 pipeBase.CmdLineTask): 

129 """Task to measure pixel ratios to find crosstalk. 

130 """ 

131 ConfigClass = CrosstalkExtractConfig 

132 _DefaultName = 'cpCrosstalkExtract' 

133 

134 def run(self, inputExp, sourceExps=[]): 

135 """Measure pixel ratios between amplifiers in inputExp. 

136 

137 Extract crosstalk ratios between different amplifiers. 

138 

139 For pixels above ``config.threshold``, we calculate the ratio 

140 between each background-subtracted target amp and the source 

141 amp. We return a list of ratios for each pixel for each 

142 target/source combination, as nested dictionary containing the 

143 ratio. 

144 

145 Parameters 

146 ---------- 

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

148 Input exposure to measure pixel ratios on. 

149 sourceExp : `list` [`lsst.afw.image.Exposure`], optional 

150 List of chips to use as sources to measure inter-chip 

151 crosstalk. 

152 

153 Returns 

154 ------- 

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

156 The results struct containing: 

157 

158 ``outputRatios`` : `dict` [`dict` [`dict` [`dict` [`list`]]]] 

159 A catalog of ratio lists. The dictionaries are 

160 indexed such that: 

161 outputRatios[targetChip][sourceChip][targetAmp][sourceAmp] 

162 contains the ratio list for that combination. 

163 ``outputFluxes`` : `dict` [`dict` [`list`]] 

164 A catalog of flux lists. The dictionaries are 

165 indexed such that: 

166 outputFluxes[sourceChip][sourceAmp] 

167 contains the flux list used in the outputRatios. 

168 

169 Notes 

170 ----- 

171 The lsstDebug.Info() method can be rewritten for __name__ = 

172 `lsst.cp.pipe.measureCrosstalk`, and supports the parameters: 

173 

174 debug.display['extract'] : `bool` 

175 Display the exposure under consideration, with the pixels used 

176 for crosstalk measurement indicated by the DETECTED mask plane. 

177 debug.display['pixels'] : `bool` 

178 Display a plot of the ratio calculated for each pixel used in this 

179 exposure, split by amplifier pairs. The median value is listed 

180 for reference. 

181 """ 

182 outputRatios = defaultdict(lambda: defaultdict(dict)) 

183 outputFluxes = defaultdict(lambda: defaultdict(dict)) 

184 

185 threshold = self.config.threshold 

186 badPixels = list(self.config.badMask) 

187 

188 targetDetector = inputExp.getDetector() 

189 targetChip = targetDetector.getName() 

190 

191 # Always look at the target chip first, then go to any other supplied exposures. 

192 sourceExtractExps = [inputExp] 

193 sourceExtractExps.extend(sourceExps) 

194 

195 self.log.info("Measuring full detector background for target: %s", targetChip) 

196 targetIm = inputExp.getMaskedImage() 

197 FootprintSet(targetIm, Threshold(threshold), "DETECTED") 

198 detected = targetIm.getMask().getPlaneBitMask("DETECTED") 

199 bg = CrosstalkCalib.calculateBackground(targetIm, badPixels + ["DETECTED"]) 

200 

201 self.debugView('extract', inputExp) 

202 

203 for sourceExp in sourceExtractExps: 

204 sourceDetector = sourceExp.getDetector() 

205 sourceChip = sourceDetector.getName() 

206 sourceIm = sourceExp.getMaskedImage() 

207 bad = sourceIm.getMask().getPlaneBitMask(badPixels) 

208 self.log.info("Measuring crosstalk from source: %s", sourceChip) 

209 

210 if sourceExp != inputExp: 

211 FootprintSet(sourceIm, Threshold(threshold), "DETECTED") 

212 detected = sourceIm.getMask().getPlaneBitMask("DETECTED") 

213 

214 # The dictionary of amp-to-amp ratios for this pair of source->target detectors. 

215 ratioDict = defaultdict(lambda: defaultdict(list)) 

216 extractedCount = 0 

217 

218 for sourceAmp in sourceDetector: 

219 sourceAmpName = sourceAmp.getName() 

220 sourceAmpImage = sourceIm[sourceAmp.getBBox()] 

221 sourceMask = sourceAmpImage.mask.array 

222 select = ((sourceMask & detected > 0) & 

223 (sourceMask & bad == 0) & 

224 np.isfinite(sourceAmpImage.image.array)) 

225 count = np.sum(select) 

226 self.log.debug(" Source amplifier: %s", sourceAmpName) 

227 

228 outputFluxes[sourceChip][sourceAmpName] = sourceAmpImage.image.array[select] 

229 

230 for targetAmp in targetDetector: 

231 # iterate over targetExposure 

232 targetAmpName = targetAmp.getName() 

233 if sourceAmpName == targetAmpName and sourceChip == targetChip: 

234 ratioDict[sourceAmpName][targetAmpName] = [] 

235 continue 

236 self.log.debug(" Target amplifier: %s", targetAmpName) 

237 

238 targetAmpImage = CrosstalkCalib.extractAmp(targetIm.image, 

239 targetAmp, sourceAmp, 

240 isTrimmed=self.config.isTrimmed) 

241 ratios = (targetAmpImage.array[select] - bg)/sourceAmpImage.image.array[select] 

242 ratioDict[targetAmpName][sourceAmpName] = ratios.tolist() 

243 extractedCount += count 

244 

245 self.debugPixels('pixels', 

246 sourceAmpImage.image.array[select], 

247 targetAmpImage.array[select] - bg, 

248 sourceAmpName, targetAmpName) 

249 

250 self.log.info("Extracted %d pixels from %s -> %s (targetBG: %f)", 

251 extractedCount, sourceChip, targetChip, bg) 

252 outputRatios[targetChip][sourceChip] = ratioDict 

253 

254 return pipeBase.Struct( 

255 outputRatios=outputRatios, 

256 outputFluxes=outputFluxes 

257 ) 

258 

259 def debugView(self, stepname, exposure): 

260 """Utility function to examine the image being processed. 

261 

262 Parameters 

263 ---------- 

264 stepname : `str` 

265 State of processing to view. 

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

267 Exposure to view. 

268 """ 

269 frame = getDebugFrame(self._display, stepname) 

270 if frame: 

271 display = getDisplay(frame) 

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

273 display.mtv(exposure) 

274 

275 prompt = "Press Enter to continue: " 

276 while True: 

277 ans = input(prompt).lower() 

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

279 break 

280 

281 def debugPixels(self, stepname, pixelsIn, pixelsOut, sourceName, targetName): 

282 """Utility function to examine the CT ratio pixel values. 

283 

284 Parameters 

285 ---------- 

286 stepname : `str` 

287 State of processing to view. 

288 pixelsIn : `np.ndarray` 

289 Pixel values from the potential crosstalk source. 

290 pixelsOut : `np.ndarray` 

291 Pixel values from the potential crosstalk target. 

292 sourceName : `str` 

293 Source amplifier name 

294 targetName : `str` 

295 Target amplifier name 

296 """ 

297 frame = getDebugFrame(self._display, stepname) 

298 if frame: 

299 import matplotlib.pyplot as plt 

300 figure = plt.figure(1) 

301 figure.clear() 

302 

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

304 axes.plt(pixelsIn, pixelsOut / pixelsIn, 'k+') 

305 plt.xlabel("Source amplifier pixel value") 

306 plt.ylabel("Measured pixel ratio") 

307 plt.title(f"(Source {sourceName} -> Target {targetName}) median ratio: " 

308 f"{(np.median(pixelsOut / pixelsIn))}") 

309 figure.show() 

310 

311 prompt = "Press Enter to continue: " 

312 while True: 

313 ans = input(prompt).lower() 

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

315 break 

316 plt.close() 

317 

318 

319class CrosstalkSolveConnections(pipeBase.PipelineTaskConnections, 

320 dimensions=("instrument", "detector")): 

321 inputRatios = cT.Input( 

322 name="crosstalkRatios", 

323 doc="Ratios measured for an input exposure.", 

324 storageClass="StructuredDataDict", 

325 dimensions=("instrument", "exposure", "detector"), 

326 multiple=True, 

327 ) 

328 inputFluxes = cT.Input( 

329 name="crosstalkFluxes", 

330 doc="Fluxes of CT source pixels, for nonlinear fits.", 

331 storageClass="StructuredDataDict", 

332 dimensions=("instrument", "exposure", "detector"), 

333 multiple=True, 

334 ) 

335 camera = cT.PrerequisiteInput( 

336 name="camera", 

337 doc="Camera the input data comes from.", 

338 storageClass="Camera", 

339 dimensions=("instrument",), 

340 isCalibration=True, 

341 lookupFunction=lookupStaticCalibration, 

342 ) 

343 

344 outputCrosstalk = cT.Output( 

345 name="crosstalk", 

346 doc="Output proposed crosstalk calibration.", 

347 storageClass="CrosstalkCalib", 

348 dimensions=("instrument", "detector"), 

349 multiple=False, 

350 isCalibration=True, 

351 ) 

352 

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

354 super().__init__(config=config) 

355 

356 if config.fluxOrder == 0: 

357 self.inputs.discard("inputFluxes") 

358 

359 

360class CrosstalkSolveConfig(pipeBase.PipelineTaskConfig, 

361 pipelineConnections=CrosstalkSolveConnections): 

362 """Configuration for the solving of crosstalk from pixel ratios. 

363 """ 

364 rejIter = Field( 

365 dtype=int, 

366 default=3, 

367 doc="Number of rejection iterations for final coefficient calculation.", 

368 ) 

369 rejSigma = Field( 

370 dtype=float, 

371 default=2.0, 

372 doc="Rejection threshold (sigma) for final coefficient calculation.", 

373 ) 

374 fluxOrder = Field( 

375 dtype=int, 

376 default=0, 

377 doc="Polynomial order in source flux to fit crosstalk.", 

378 ) 

379 doFiltering = Field( 

380 dtype=bool, 

381 default=False, 

382 doc="Filter generated crosstalk to remove marginal measurements.", 

383 ) 

384 

385 

386class CrosstalkSolveTask(pipeBase.PipelineTask, 

387 pipeBase.CmdLineTask): 

388 """Task to solve crosstalk from pixel ratios. 

389 """ 

390 ConfigClass = CrosstalkSolveConfig 

391 _DefaultName = 'cpCrosstalkSolve' 

392 

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

394 """Ensure that the input and output dimensions are passed along. 

395 

396 Parameters 

397 ---------- 

398 butlerQC : `lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext` 

399 Butler to operate on. 

400 inputRefs : `lsst.pipe.base.connections.InputQuantizedConnection` 

401 Input data refs to load. 

402 ouptutRefs : `lsst.pipe.base.connections.OutputQuantizedConnection` 

403 Output data refs to persist. 

404 """ 

405 inputs = butlerQC.get(inputRefs) 

406 

407 # Use the dimensions to set calib/provenance information. 

408 inputs['inputDims'] = [exp.dataId.byName() for exp in inputRefs.inputRatios] 

409 inputs['outputDims'] = outputRefs.outputCrosstalk.dataId.byName() 

410 

411 outputs = self.run(**inputs) 

412 butlerQC.put(outputs, outputRefs) 

413 

414 def run(self, inputRatios, inputFluxes=None, camera=None, inputDims=None, outputDims=None): 

415 """Combine ratios to produce crosstalk coefficients. 

416 

417 Parameters 

418 ---------- 

419 inputRatios : `list` [`dict` [`dict` [`dict` [`dict` [`list`]]]]] 

420 A list of nested dictionaries of ratios indexed by target 

421 and source chip, then by target and source amplifier. 

422 inputFluxes : `list` [`dict` [`dict` [`list`]]] 

423 A list of nested dictionaries of source pixel fluxes, indexed 

424 by source chip and amplifier. 

425 camera : `lsst.afw.cameraGeom.Camera` 

426 Input camera. 

427 inputDims : `list` [`lsst.daf.butler.DataCoordinate`] 

428 DataIds to use to construct provenance. 

429 outputDims : `list` [`lsst.daf.butler.DataCoordinate`] 

430 DataIds to use to populate the output calibration. 

431 

432 Returns 

433 ------- 

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

435 The results struct containing: 

436 

437 ``outputCrosstalk`` : `lsst.ip.isr.CrosstalkCalib` 

438 Final crosstalk calibration. 

439 ``outputProvenance`` : `lsst.ip.isr.IsrProvenance` 

440 Provenance data for the new calibration. 

441 

442 Raises 

443 ------ 

444 RuntimeError 

445 Raised if the input data contains multiple target detectors. 

446 

447 Notes 

448 ----- 

449 The lsstDebug.Info() method can be rewritten for __name__ = 

450 `lsst.ip.isr.measureCrosstalk`, and supports the parameters: 

451 

452 debug.display['reduce'] : `bool` 

453 Display a histogram of the combined ratio measurements for 

454 a pair of source/target amplifiers from all input 

455 exposures/detectors. 

456 

457 """ 

458 if outputDims: 

459 calibChip = outputDims['detector'] 

460 instrument = outputDims['instrument'] 

461 else: 

462 # calibChip needs to be set manually in Gen2. 

463 calibChip = None 

464 instrument = None 

465 

466 self.log.info("Combining measurements from %d ratios and %d fluxes", 

467 len(inputRatios), len(inputFluxes) if inputFluxes else 0) 

468 

469 if inputFluxes is None: 

470 inputFluxes = [None for exp in inputRatios] 

471 

472 combinedRatios = defaultdict(lambda: defaultdict(list)) 

473 combinedFluxes = defaultdict(lambda: defaultdict(list)) 

474 for ratioDict, fluxDict in zip(inputRatios, inputFluxes): 

475 for targetChip in ratioDict: 

476 if calibChip and targetChip != calibChip: 

477 raise RuntimeError("Received multiple target chips!") 

478 

479 sourceChip = targetChip 

480 if sourceChip in ratioDict[targetChip]: 

481 ratios = ratioDict[targetChip][sourceChip] 

482 

483 for targetAmp in ratios: 

484 for sourceAmp in ratios[targetAmp]: 

485 combinedRatios[targetAmp][sourceAmp].extend(ratios[targetAmp][sourceAmp]) 

486 if fluxDict: 

487 combinedFluxes[targetAmp][sourceAmp].extend(fluxDict[sourceChip][sourceAmp]) 

488 # TODO: DM-21904 

489 # Iterating over all other entries in ratioDict[targetChip] will yield 

490 # inter-chip terms. 

491 

492 for targetAmp in combinedRatios: 

493 for sourceAmp in combinedRatios[targetAmp]: 

494 self.log.info("Read %d pixels for %s -> %s", 

495 len(combinedRatios[targetAmp][sourceAmp]), 

496 targetAmp, sourceAmp) 

497 if len(combinedRatios[targetAmp][sourceAmp]) > 1: 

498 self.debugRatios('reduce', combinedRatios, targetAmp, sourceAmp) 

499 

500 if self.config.fluxOrder == 0: 

501 self.log.info("Fitting crosstalk coefficients.") 

502 calib = self.measureCrosstalkCoefficients(combinedRatios, 

503 self.config.rejIter, self.config.rejSigma) 

504 else: 

505 raise NotImplementedError("Non-linear crosstalk terms are not yet supported.") 

506 

507 self.log.info("Number of valid coefficients: %d", np.sum(calib.coeffValid)) 

508 

509 if self.config.doFiltering: 

510 # This step will apply the calculated validity values to 

511 # censor poorly measured coefficients. 

512 self.log.info("Filtering measured crosstalk to remove invalid solutions.") 

513 calib = self.filterCrosstalkCalib(calib) 

514 

515 # Populate the remainder of the calibration information. 

516 calib.hasCrosstalk = True 

517 calib.interChip = {} 

518 calib._detectorName = calibChip 

519 if camera: 

520 for chip in camera: 

521 if chip.getName() == calibChip: 

522 calib._detectorSerial = chip.getSerial() 

523 calib._instrument = instrument 

524 calib.updateMetadata() 

525 

526 # Make an IsrProvenance(). 

527 provenance = IsrProvenance(calibType="CROSSTALK") 

528 provenance._detectorName = calibChip 

529 if inputDims: 

530 provenance.fromDataIds(inputDims) 

531 provenance._instrument = instrument 

532 provenance.updateMetadata() 

533 

534 return pipeBase.Struct( 

535 outputCrosstalk=calib, 

536 outputProvenance=provenance, 

537 ) 

538 

539 def measureCrosstalkCoefficients(self, ratios, rejIter, rejSigma): 

540 """Measure crosstalk coefficients from the ratios. 

541 

542 Given a list of ratios for each target/source amp combination, 

543 we measure a sigma clipped mean and error. 

544 

545 The coefficient errors returned are the standard deviation of 

546 the final set of clipped input ratios. 

547 

548 Parameters 

549 ---------- 

550 ratios : `dict` of `dict` of `numpy.ndarray` 

551 Catalog of arrays of ratios. 

552 rejIter : `int` 

553 Number of rejection iterations. 

554 rejSigma : `float` 

555 Rejection threshold (sigma). 

556 

557 Returns 

558 ------- 

559 calib : `lsst.ip.isr.CrosstalkCalib` 

560 The output crosstalk calibration. 

561 

562 Notes 

563 ----- 

564 The lsstDebug.Info() method can be rewritten for __name__ = 

565 `lsst.ip.isr.measureCrosstalk`, and supports the parameters: 

566 

567 debug.display['measure'] : `bool` 

568 Display the CDF of the combined ratio measurements for 

569 a pair of source/target amplifiers from the final set of 

570 clipped input ratios. 

571 """ 

572 calib = CrosstalkCalib(nAmp=len(ratios)) 

573 

574 # Calibration stores coefficients as a numpy ndarray. 

575 ordering = list(ratios.keys()) 

576 for ii, jj in itertools.product(range(calib.nAmp), range(calib.nAmp)): 

577 if ii == jj: 

578 values = [0.0] 

579 else: 

580 values = np.array(ratios[ordering[ii]][ordering[jj]]) 

581 values = values[np.abs(values) < 1.0] # Discard unreasonable values 

582 

583 calib.coeffNum[ii][jj] = len(values) 

584 

585 if len(values) == 0: 

586 self.log.warn("No values for matrix element %d,%d" % (ii, jj)) 

587 calib.coeffs[ii][jj] = np.nan 

588 calib.coeffErr[ii][jj] = np.nan 

589 calib.coeffValid[ii][jj] = False 

590 else: 

591 if ii != jj: 

592 for rej in range(rejIter): 

593 lo, med, hi = np.percentile(values, [25.0, 50.0, 75.0]) 

594 sigma = 0.741*(hi - lo) 

595 good = np.abs(values - med) < rejSigma*sigma 

596 if good.sum() == len(good): 

597 break 

598 values = values[good] 

599 

600 calib.coeffs[ii][jj] = np.mean(values) 

601 if calib.coeffNum[ii][jj] == 1: 

602 calib.coeffErr[ii][jj] = np.nan 

603 else: 

604 correctionFactor = self.sigmaClipCorrection(rejSigma) 

605 calib.coeffErr[ii][jj] = np.std(values) * correctionFactor 

606 calib.coeffValid[ii][jj] = (np.abs(calib.coeffs[ii][jj]) > 

607 calib.coeffErr[ii][jj] / np.sqrt(calib.coeffNum[ii][jj])) 

608 

609 if calib.coeffNum[ii][jj] > 1: 

610 self.debugRatios('measure', ratios, ordering[ii], ordering[jj], 

611 calib.coeffs[ii][jj], calib.coeffValid[ii][jj]) 

612 

613 return calib 

614 

615 @staticmethod 

616 def sigmaClipCorrection(nSigClip): 

617 """Correct measured sigma to account for clipping. 

618 

619 If we clip our input data and then measure sigma, then the 

620 measured sigma is smaller than the true value because real 

621 points beyond the clip threshold have been removed. This is a 

622 small (1.5% at nSigClip=3) effect when nSigClip >~ 3, but the 

623 default parameters for measure crosstalk use nSigClip=2.0. 

624 This causes the measured sigma to be about 15% smaller than 

625 real. This formula corrects the issue, for the symmetric case 

626 (upper clip threshold equal to lower clip threshold). 

627 

628 Parameters 

629 ---------- 

630 nSigClip : `float` 

631 Number of sigma the measurement was clipped by. 

632 

633 Returns 

634 ------- 

635 scaleFactor : `float` 

636 Scale factor to increase the measured sigma by. 

637 

638 """ 

639 varFactor = 1.0 + (2 * nSigClip * norm.pdf(nSigClip)) / (norm.cdf(nSigClip) - norm.cdf(-nSigClip)) 

640 return 1.0 / np.sqrt(varFactor) 

641 

642 @staticmethod 

643 def filterCrosstalkCalib(inCalib): 

644 """Apply valid constraints to the measured values. 

645 

646 Any measured coefficient that is determined to be invalid is 

647 set to zero, and has the error set to nan. The validation is 

648 determined by checking that the measured coefficient is larger 

649 than the calculated standard error of the mean. 

650 

651 Parameters 

652 ---------- 

653 inCalib : `lsst.ip.isr.CrosstalkCalib` 

654 Input calibration to filter. 

655 

656 Returns 

657 ------- 

658 outCalib : `lsst.ip.isr.CrosstalkCalib` 

659 Filtered calibration. 

660 """ 

661 outCalib = CrosstalkCalib() 

662 outCalib.numAmps = inCalib.numAmps 

663 

664 outCalib.coeffs = inCalib.coeffs 

665 outCalib.coeffs[~inCalib.coeffValid] = 0.0 

666 

667 outCalib.coeffErr = inCalib.coeffErr 

668 outCalib.coeffErr[~inCalib.coeffValid] = np.nan 

669 

670 outCalib.coeffNum = inCalib.coeffNum 

671 outCalib.coeffValid = inCalib.coeffValid 

672 

673 return outCalib 

674 

675 def debugRatios(self, stepname, ratios, i, j, coeff=0.0, valid=False): 

676 """Utility function to examine the final CT ratio set. 

677 

678 Parameters 

679 ---------- 

680 stepname : `str` 

681 State of processing to view. 

682 ratios : `dict` of `dict` of `np.ndarray` 

683 Array of measured CT ratios, indexed by source/victim 

684 amplifier. 

685 i : `str` 

686 Index of the source amplifier. 

687 j : `str` 

688 Index of the target amplifier. 

689 coeff : `float`, optional 

690 Coefficient calculated to plot along with the simple mean. 

691 valid : `bool`, optional 

692 Validity to be added to the plot title. 

693 """ 

694 frame = getDebugFrame(self._display, stepname) 

695 if frame: 

696 if i == j or ratios is None or len(ratios) < 1: 

697 pass 

698 

699 ratioList = ratios[i][j] 

700 if ratioList is None or len(ratioList) < 1: 

701 pass 

702 

703 mean = np.mean(ratioList) 

704 std = np.std(ratioList) 

705 import matplotlib.pyplot as plt 

706 figure = plt.figure(1) 

707 figure.clear() 

708 plt.hist(x=ratioList, bins=len(ratioList), 

709 cumulative=True, color='b', density=True, histtype='step') 

710 plt.xlabel("Measured pixel ratio") 

711 plt.ylabel(f"CDF: n={len(ratioList)}") 

712 plt.xlim(np.percentile(ratioList, [1.0, 99])) 

713 plt.axvline(x=mean, color="k") 

714 plt.axvline(x=coeff, color='g') 

715 plt.axvline(x=(std / np.sqrt(len(ratioList))), color='r') 

716 plt.axvline(x=-(std / np.sqrt(len(ratioList))), color='r') 

717 plt.title(f"(Source {i} -> Target {j}) mean: {mean:.2g} coeff: {coeff:.2g} valid: {valid}") 

718 figure.show() 

719 

720 prompt = "Press Enter to continue: " 

721 while True: 

722 ans = input(prompt).lower() 

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

724 break 

725 elif ans in ("pdb", "p",): 

726 import pdb 

727 pdb.set_trace() 

728 plt.close() 

729 

730 

731class MeasureCrosstalkConfig(Config): 

732 extract = ConfigurableField( 

733 target=CrosstalkExtractTask, 

734 doc="Task to measure pixel ratios.", 

735 ) 

736 solver = ConfigurableField( 

737 target=CrosstalkSolveTask, 

738 doc="Task to convert ratio lists to crosstalk coefficients.", 

739 ) 

740 

741 

742class MeasureCrosstalkTask(pipeBase.CmdLineTask): 

743 """Measure intra-detector crosstalk. 

744 

745 See also 

746 -------- 

747 lsst.ip.isr.crosstalk.CrosstalkCalib 

748 lsst.cp.pipe.measureCrosstalk.CrosstalkExtractTask 

749 lsst.cp.pipe.measureCrosstalk.CrosstalkSolveTask 

750 

751 Notes 

752 ----- 

753 The crosstalk this method measures assumes that when a bright 

754 pixel is found in one detector amplifier, all other detector 

755 amplifiers may see a signal change in the same pixel location 

756 (relative to the readout amplifier) as these other pixels are read 

757 out at the same time. 

758 

759 After processing each input exposure through a limited set of ISR 

760 stages, bright unmasked pixels above the threshold are identified. 

761 The potential CT signal is found by taking the ratio of the 

762 appropriate background-subtracted pixel value on the other 

763 amplifiers to the input value on the source amplifier. If the 

764 source amplifier has a large number of bright pixels as well, the 

765 background level may be elevated, leading to poor ratio 

766 measurements. 

767 

768 The set of ratios found between each pair of amplifiers across all 

769 input exposures is then gathered to produce the final CT 

770 coefficients. The sigma-clipped mean and sigma are returned from 

771 these sets of ratios, with the coefficient to supply to the ISR 

772 CrosstalkTask() being the multiplicative inverse of these values. 

773 

774 This Task simply calls the pipetask versions of the measure 

775 crosstalk code. 

776 """ 

777 ConfigClass = MeasureCrosstalkConfig 

778 _DefaultName = "measureCrosstalk" 

779 

780 # Let's use this instead of messing with parseAndRun. 

781 RunnerClass = DataRefListRunner 

782 

783 def __init__(self, **kwargs): 

784 super().__init__(**kwargs) 

785 self.makeSubtask("extract") 

786 self.makeSubtask("solver") 

787 

788 def runDataRef(self, dataRefList): 

789 """Run extract task on each of inputs in the dataRef list, then pass 

790 that to the solver task. 

791 

792 Parameters 

793 ---------- 

794 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`] 

795 Data references for exposures for detectors to process. 

796 

797 Returns 

798 ------- 

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

800 The results struct containing: 

801 

802 ``outputCrosstalk`` : `lsst.ip.isr.CrosstalkCalib` 

803 Final crosstalk calibration. 

804 ``outputProvenance`` : `lsst.ip.isr.IsrProvenance` 

805 Provenance data for the new calibration. 

806 

807 Raises 

808 ------ 

809 RuntimeError 

810 Raised if multiple target detectors are supplied. 

811 """ 

812 dataRef = dataRefList[0] 

813 camera = dataRef.get("camera") 

814 

815 ratios = [] 

816 activeChip = None 

817 for dataRef in dataRefList: 

818 exposure = dataRef.get("postISRCCD") 

819 if activeChip: 

820 if exposure.getDetector().getName() != activeChip: 

821 raise RuntimeError("Too many input detectors supplied!") 

822 else: 

823 activeChip = exposure.getDetector().getName() 

824 

825 self.extract.debugView("extract", exposure) 

826 result = self.extract.run(exposure) 

827 ratios.append(result.outputRatios) 

828 

829 finalResults = self.solver.run(ratios, camera=camera) 

830 dataRef.put(finalResults.outputCrosstalk, "crosstalk") 

831 

832 return finalResults