Coverage for python/lsst/cp/pipe/measureCrosstalk.py : 16%

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
25from collections import defaultdict
27import lsst.pipe.base as pipeBase
28import lsst.pipe.base.connectionTypes as cT
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
36from lsst.cp.pipe.utils import ddict2dict
38from ._lookupStaticCalibration import lookupStaticCalibration
40__all__ = ["CrosstalkExtractConfig", "CrosstalkExtractTask",
41 "CrosstalkSolveTask", "CrosstalkSolveConfig",
42 "MeasureCrosstalkConfig", "MeasureCrosstalkTask"]
45class CrosstalkExtractConnections(pipeBase.PipelineTaskConnections,
46 dimensions=("instrument", "exposure", "detector")):
47 inputExp = cT.Input(
48 name="crosstalkInputs",
49 doc="Input post-ISR processed exposure to measure crosstalk from.",
50 storageClass="Exposure",
51 dimensions=("instrument", "exposure", "detector"),
52 multiple=False,
53 )
54 # TODO: Depends on DM-21904.
55 sourceExp = cT.Input(
56 name="crosstalkSource",
57 doc="Post-ISR exposure to measure for inter-chip crosstalk onto inputExp.",
58 storageClass="Exposure",
59 dimensions=("instrument", "exposure", "detector"),
60 multiple=True,
61 deferLoad=True,
62 # lookupFunction=None,
63 )
65 outputRatios = cT.Output(
66 name="crosstalkRatios",
67 doc="Extracted crosstalk pixel ratios.",
68 storageClass="StructuredDataDict",
69 dimensions=("instrument", "exposure", "detector"),
70 )
71 outputFluxes = cT.Output(
72 name="crosstalkFluxes",
73 doc="Source pixel fluxes used in ratios.",
74 storageClass="StructuredDataDict",
75 dimensions=("instrument", "exposure", "detector"),
76 )
78 def __init__(self, *, config=None):
79 super().__init__(config=config)
80 # Discard sourceExp until DM-21904 allows full interchip
81 # measurements.
82 self.inputs.discard("sourceExp")
85class CrosstalkExtractConfig(pipeBase.PipelineTaskConfig,
86 pipelineConnections=CrosstalkExtractConnections):
87 """Configuration for the measurement of pixel ratios.
88 """
89 doMeasureInterchip = Field(
90 dtype=bool,
91 default=False,
92 doc="Measure inter-chip crosstalk as well?",
93 )
94 threshold = Field(
95 dtype=float,
96 default=30000,
97 doc="Minimum level of source pixels for which to measure crosstalk."
98 )
99 ignoreSaturatedPixels = Field(
100 dtype=bool,
101 default=True,
102 doc="Should saturated pixels be ignored?"
103 )
104 badMask = ListField(
105 dtype=str,
106 default=["BAD", "INTRP"],
107 doc="Mask planes to ignore when identifying source pixels."
108 )
109 isTrimmed = Field(
110 dtype=bool,
111 default=True,
112 doc="Is the input exposure trimmed?"
113 )
115 def validate(self):
116 super().validate()
118 # Ensure the handling of the SAT mask plane is consistent
119 # with the ignoreSaturatedPixels value.
120 if self.ignoreSaturatedPixels:
121 if 'SAT' not in self.badMask:
122 self.badMask.append('SAT')
123 else:
124 if 'SAT' in self.badMask:
125 self.badMask = [mask for mask in self.badMask if mask != 'SAT']
128class CrosstalkExtractTask(pipeBase.PipelineTask,
129 pipeBase.CmdLineTask):
130 """Task to measure pixel ratios to find crosstalk.
131 """
132 ConfigClass = CrosstalkExtractConfig
133 _DefaultName = 'cpCrosstalkExtract'
135 def run(self, inputExp, sourceExps=[]):
136 """Measure pixel ratios between amplifiers in inputExp.
138 Extract crosstalk ratios between different amplifiers.
140 For pixels above ``config.threshold``, we calculate the ratio
141 between each background-subtracted target amp and the source
142 amp. We return a list of ratios for each pixel for each
143 target/source combination, as nested dictionary containing the
144 ratio.
146 Parameters
147 ----------
148 inputExp : `lsst.afw.image.Exposure`
149 Input exposure to measure pixel ratios on.
150 sourceExp : `list` [`lsst.afw.image.Exposure`], optional
151 List of chips to use as sources to measure inter-chip
152 crosstalk.
154 Returns
155 -------
156 results : `lsst.pipe.base.Struct`
157 The results struct containing:
159 ``outputRatios`` : `dict` [`dict` [`dict` [`dict` [`list`]]]]
160 A catalog of ratio lists. The dictionaries are
161 indexed such that:
162 outputRatios[targetChip][sourceChip][targetAmp][sourceAmp]
163 contains the ratio list for that combination.
164 ``outputFluxes`` : `dict` [`dict` [`list`]]
165 A catalog of flux lists. The dictionaries are
166 indexed such that:
167 outputFluxes[sourceChip][sourceAmp]
168 contains the flux list used in the outputRatios.
170 Notes
171 -----
172 The lsstDebug.Info() method can be rewritten for __name__ =
173 `lsst.cp.pipe.measureCrosstalk`, and supports the parameters:
175 debug.display['extract'] : `bool`
176 Display the exposure under consideration, with the pixels used
177 for crosstalk measurement indicated by the DETECTED mask plane.
178 debug.display['pixels'] : `bool`
179 Display a plot of the ratio calculated for each pixel used in this
180 exposure, split by amplifier pairs. The median value is listed
181 for reference.
182 """
183 outputRatios = defaultdict(lambda: defaultdict(dict))
184 outputFluxes = defaultdict(lambda: defaultdict(dict))
186 threshold = self.config.threshold
187 badPixels = list(self.config.badMask)
189 targetDetector = inputExp.getDetector()
190 targetChip = targetDetector.getName()
192 # Always look at the target chip first, then go to any other supplied exposures.
193 sourceExtractExps = [inputExp]
194 sourceExtractExps.extend(sourceExps)
196 self.log.info("Measuring full detector background for target: %s", targetChip)
197 targetIm = inputExp.getMaskedImage()
198 FootprintSet(targetIm, Threshold(threshold), "DETECTED")
199 detected = targetIm.getMask().getPlaneBitMask("DETECTED")
200 bg = CrosstalkCalib.calculateBackground(targetIm, badPixels + ["DETECTED"])
202 self.debugView('extract', inputExp)
204 for sourceExp in sourceExtractExps:
205 sourceDetector = sourceExp.getDetector()
206 sourceChip = sourceDetector.getName()
207 sourceIm = sourceExp.getMaskedImage()
208 bad = sourceIm.getMask().getPlaneBitMask(badPixels)
209 self.log.info("Measuring crosstalk from source: %s", sourceChip)
211 if sourceExp != inputExp:
212 FootprintSet(sourceIm, Threshold(threshold), "DETECTED")
213 detected = sourceIm.getMask().getPlaneBitMask("DETECTED")
215 # The dictionary of amp-to-amp ratios for this pair of source->target detectors.
216 ratioDict = defaultdict(lambda: defaultdict(list))
217 extractedCount = 0
219 for sourceAmp in sourceDetector:
220 sourceAmpName = sourceAmp.getName()
221 sourceAmpImage = sourceIm[sourceAmp.getBBox()]
222 sourceMask = sourceAmpImage.mask.array
223 select = ((sourceMask & detected > 0) &
224 (sourceMask & bad == 0) &
225 np.isfinite(sourceAmpImage.image.array))
226 count = np.sum(select)
227 self.log.debug(" Source amplifier: %s", sourceAmpName)
229 outputFluxes[sourceChip][sourceAmpName] = sourceAmpImage.image.array[select].tolist()
231 for targetAmp in targetDetector:
232 # iterate over targetExposure
233 targetAmpName = targetAmp.getName()
234 if sourceAmpName == targetAmpName and sourceChip == targetChip:
235 ratioDict[sourceAmpName][targetAmpName] = []
236 continue
237 self.log.debug(" Target amplifier: %s", targetAmpName)
239 targetAmpImage = CrosstalkCalib.extractAmp(targetIm.image,
240 targetAmp, sourceAmp,
241 isTrimmed=self.config.isTrimmed)
242 ratios = (targetAmpImage.array[select] - bg)/sourceAmpImage.image.array[select]
243 ratioDict[targetAmpName][sourceAmpName] = ratios.tolist()
244 extractedCount += count
246 self.debugPixels('pixels',
247 sourceAmpImage.image.array[select],
248 targetAmpImage.array[select] - bg,
249 sourceAmpName, targetAmpName)
251 self.log.info("Extracted %d pixels from %s -> %s (targetBG: %f)",
252 extractedCount, sourceChip, targetChip, bg)
253 outputRatios[targetChip][sourceChip] = ratioDict
255 return pipeBase.Struct(
256 outputRatios=ddict2dict(outputRatios),
257 outputFluxes=ddict2dict(outputFluxes)
258 )
260 def debugView(self, stepname, exposure):
261 """Utility function to examine the image being processed.
263 Parameters
264 ----------
265 stepname : `str`
266 State of processing to view.
267 exposure : `lsst.afw.image.Exposure`
268 Exposure to view.
269 """
270 frame = getDebugFrame(self._display, stepname)
271 if frame:
272 display = getDisplay(frame)
273 display.scale('asinh', 'zscale')
274 display.mtv(exposure)
276 prompt = "Press Enter to continue: "
277 while True:
278 ans = input(prompt).lower()
279 if ans in ("", "c",):
280 break
282 def debugPixels(self, stepname, pixelsIn, pixelsOut, sourceName, targetName):
283 """Utility function to examine the CT ratio pixel values.
285 Parameters
286 ----------
287 stepname : `str`
288 State of processing to view.
289 pixelsIn : `np.ndarray`
290 Pixel values from the potential crosstalk source.
291 pixelsOut : `np.ndarray`
292 Pixel values from the potential crosstalk target.
293 sourceName : `str`
294 Source amplifier name
295 targetName : `str`
296 Target amplifier name
297 """
298 frame = getDebugFrame(self._display, stepname)
299 if frame:
300 import matplotlib.pyplot as plt
301 figure = plt.figure(1)
302 figure.clear()
304 axes = figure.add_axes((0.1, 0.1, 0.8, 0.8))
305 axes.plot(pixelsIn, pixelsOut / pixelsIn, 'k+')
306 plt.xlabel("Source amplifier pixel value")
307 plt.ylabel("Measured pixel ratio")
308 plt.title(f"(Source {sourceName} -> Target {targetName}) median ratio: "
309 f"{(np.median(pixelsOut / pixelsIn))}")
310 figure.show()
312 prompt = "Press Enter to continue: "
313 while True:
314 ans = input(prompt).lower()
315 if ans in ("", "c",):
316 break
317 plt.close()
320class CrosstalkSolveConnections(pipeBase.PipelineTaskConnections,
321 dimensions=("instrument", "detector")):
322 inputRatios = cT.Input(
323 name="crosstalkRatios",
324 doc="Ratios measured for an input exposure.",
325 storageClass="StructuredDataDict",
326 dimensions=("instrument", "exposure", "detector"),
327 multiple=True,
328 )
329 inputFluxes = cT.Input(
330 name="crosstalkFluxes",
331 doc="Fluxes of CT source pixels, for nonlinear fits.",
332 storageClass="StructuredDataDict",
333 dimensions=("instrument", "exposure", "detector"),
334 multiple=True,
335 )
336 camera = cT.PrerequisiteInput(
337 name="camera",
338 doc="Camera the input data comes from.",
339 storageClass="Camera",
340 dimensions=("instrument",),
341 isCalibration=True,
342 lookupFunction=lookupStaticCalibration,
343 )
345 outputCrosstalk = cT.Output(
346 name="crosstalk",
347 doc="Output proposed crosstalk calibration.",
348 storageClass="CrosstalkCalib",
349 dimensions=("instrument", "detector"),
350 multiple=False,
351 isCalibration=True,
352 )
354 def __init__(self, *, config=None):
355 super().__init__(config=config)
357 if config.fluxOrder == 0:
358 self.inputs.discard("inputFluxes")
361class CrosstalkSolveConfig(pipeBase.PipelineTaskConfig,
362 pipelineConnections=CrosstalkSolveConnections):
363 """Configuration for the solving of crosstalk from pixel ratios.
364 """
365 rejIter = Field(
366 dtype=int,
367 default=3,
368 doc="Number of rejection iterations for final coefficient calculation.",
369 )
370 rejSigma = Field(
371 dtype=float,
372 default=2.0,
373 doc="Rejection threshold (sigma) for final coefficient calculation.",
374 )
375 fluxOrder = Field(
376 dtype=int,
377 default=0,
378 doc="Polynomial order in source flux to fit crosstalk.",
379 )
380 doFiltering = Field(
381 dtype=bool,
382 default=False,
383 doc="Filter generated crosstalk to remove marginal measurements.",
384 )
387class CrosstalkSolveTask(pipeBase.PipelineTask,
388 pipeBase.CmdLineTask):
389 """Task to solve crosstalk from pixel ratios.
390 """
391 ConfigClass = CrosstalkSolveConfig
392 _DefaultName = 'cpCrosstalkSolve'
394 def runQuantum(self, butlerQC, inputRefs, outputRefs):
395 """Ensure that the input and output dimensions are passed along.
397 Parameters
398 ----------
399 butlerQC : `lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext`
400 Butler to operate on.
401 inputRefs : `lsst.pipe.base.connections.InputQuantizedConnection`
402 Input data refs to load.
403 ouptutRefs : `lsst.pipe.base.connections.OutputQuantizedConnection`
404 Output data refs to persist.
405 """
406 inputs = butlerQC.get(inputRefs)
408 # Use the dimensions to set calib/provenance information.
409 inputs['inputDims'] = [exp.dataId.byName() for exp in inputRefs.inputRatios]
410 inputs['outputDims'] = outputRefs.outputCrosstalk.dataId.byName()
412 outputs = self.run(**inputs)
413 butlerQC.put(outputs, outputRefs)
415 def run(self, inputRatios, inputFluxes=None, camera=None, inputDims=None, outputDims=None):
416 """Combine ratios to produce crosstalk coefficients.
418 Parameters
419 ----------
420 inputRatios : `list` [`dict` [`dict` [`dict` [`dict` [`list`]]]]]
421 A list of nested dictionaries of ratios indexed by target
422 and source chip, then by target and source amplifier.
423 inputFluxes : `list` [`dict` [`dict` [`list`]]]
424 A list of nested dictionaries of source pixel fluxes, indexed
425 by source chip and amplifier.
426 camera : `lsst.afw.cameraGeom.Camera`
427 Input camera.
428 inputDims : `list` [`lsst.daf.butler.DataCoordinate`]
429 DataIds to use to construct provenance.
430 outputDims : `list` [`lsst.daf.butler.DataCoordinate`]
431 DataIds to use to populate the output calibration.
433 Returns
434 -------
435 results : `lsst.pipe.base.Struct`
436 The results struct containing:
438 ``outputCrosstalk`` : `lsst.ip.isr.CrosstalkCalib`
439 Final crosstalk calibration.
440 ``outputProvenance`` : `lsst.ip.isr.IsrProvenance`
441 Provenance data for the new calibration.
443 Raises
444 ------
445 RuntimeError
446 Raised if the input data contains multiple target detectors.
448 Notes
449 -----
450 The lsstDebug.Info() method can be rewritten for __name__ =
451 `lsst.ip.isr.measureCrosstalk`, and supports the parameters:
453 debug.display['reduce'] : `bool`
454 Display a histogram of the combined ratio measurements for
455 a pair of source/target amplifiers from all input
456 exposures/detectors.
458 """
459 if outputDims:
460 calibChip = outputDims['detector']
461 instrument = outputDims['instrument']
462 else:
463 # calibChip needs to be set manually in Gen2.
464 calibChip = None
465 instrument = None
467 self.log.info("Combining measurements from %d ratios and %d fluxes",
468 len(inputRatios), len(inputFluxes) if inputFluxes else 0)
470 if inputFluxes is None:
471 inputFluxes = [None for exp in inputRatios]
473 combinedRatios = defaultdict(lambda: defaultdict(list))
474 combinedFluxes = defaultdict(lambda: defaultdict(list))
475 for ratioDict, fluxDict in zip(inputRatios, inputFluxes):
476 for targetChip in ratioDict:
477 if calibChip and targetChip != calibChip:
478 raise RuntimeError("Received multiple target chips!")
480 sourceChip = targetChip
481 if sourceChip in ratioDict[targetChip]:
482 ratios = ratioDict[targetChip][sourceChip]
484 for targetAmp in ratios:
485 for sourceAmp in ratios[targetAmp]:
486 combinedRatios[targetAmp][sourceAmp].extend(ratios[targetAmp][sourceAmp])
487 if fluxDict:
488 combinedFluxes[targetAmp][sourceAmp].extend(fluxDict[sourceChip][sourceAmp])
489 # TODO: DM-21904
490 # Iterating over all other entries in ratioDict[targetChip] will yield
491 # inter-chip terms.
493 for targetAmp in combinedRatios:
494 for sourceAmp in combinedRatios[targetAmp]:
495 self.log.info("Read %d pixels for %s -> %s",
496 len(combinedRatios[targetAmp][sourceAmp]),
497 targetAmp, sourceAmp)
498 if len(combinedRatios[targetAmp][sourceAmp]) > 1:
499 self.debugRatios('reduce', combinedRatios, targetAmp, sourceAmp)
501 if self.config.fluxOrder == 0:
502 self.log.info("Fitting crosstalk coefficients.")
503 calib = self.measureCrosstalkCoefficients(combinedRatios,
504 self.config.rejIter, self.config.rejSigma)
505 else:
506 raise NotImplementedError("Non-linear crosstalk terms are not yet supported.")
508 self.log.info("Number of valid coefficients: %d", np.sum(calib.coeffValid))
510 if self.config.doFiltering:
511 # This step will apply the calculated validity values to
512 # censor poorly measured coefficients.
513 self.log.info("Filtering measured crosstalk to remove invalid solutions.")
514 calib = self.filterCrosstalkCalib(calib)
516 # Populate the remainder of the calibration information.
517 calib.hasCrosstalk = True
518 calib.interChip = {}
519 calib._detectorName = calibChip
520 if camera:
521 for chip in camera:
522 if chip.getName() == calibChip:
523 calib._detectorSerial = chip.getSerial()
524 calib._instrument = instrument
525 calib.updateMetadata()
527 # Make an IsrProvenance().
528 provenance = IsrProvenance(calibType="CROSSTALK")
529 provenance._detectorName = calibChip
530 if inputDims:
531 provenance.fromDataIds(inputDims)
532 provenance._instrument = instrument
533 provenance.updateMetadata()
535 return pipeBase.Struct(
536 outputCrosstalk=calib,
537 outputProvenance=provenance,
538 )
540 def measureCrosstalkCoefficients(self, ratios, rejIter, rejSigma):
541 """Measure crosstalk coefficients from the ratios.
543 Given a list of ratios for each target/source amp combination,
544 we measure a sigma clipped mean and error.
546 The coefficient errors returned are the standard deviation of
547 the final set of clipped input ratios.
549 Parameters
550 ----------
551 ratios : `dict` of `dict` of `numpy.ndarray`
552 Catalog of arrays of ratios.
553 rejIter : `int`
554 Number of rejection iterations.
555 rejSigma : `float`
556 Rejection threshold (sigma).
558 Returns
559 -------
560 calib : `lsst.ip.isr.CrosstalkCalib`
561 The output crosstalk calibration.
563 Notes
564 -----
565 The lsstDebug.Info() method can be rewritten for __name__ =
566 `lsst.ip.isr.measureCrosstalk`, and supports the parameters:
568 debug.display['measure'] : `bool`
569 Display the CDF of the combined ratio measurements for
570 a pair of source/target amplifiers from the final set of
571 clipped input ratios.
572 """
573 calib = CrosstalkCalib(nAmp=len(ratios))
575 # Calibration stores coefficients as a numpy ndarray.
576 ordering = list(ratios.keys())
577 for ii, jj in itertools.product(range(calib.nAmp), range(calib.nAmp)):
578 if ii == jj:
579 values = [0.0]
580 else:
581 values = np.array(ratios[ordering[ii]][ordering[jj]])
582 values = values[np.abs(values) < 1.0] # Discard unreasonable values
584 calib.coeffNum[ii][jj] = len(values)
586 if len(values) == 0:
587 self.log.warn("No values for matrix element %d,%d" % (ii, jj))
588 calib.coeffs[ii][jj] = np.nan
589 calib.coeffErr[ii][jj] = np.nan
590 calib.coeffValid[ii][jj] = False
591 else:
592 if ii != jj:
593 for rej in range(rejIter):
594 lo, med, hi = np.percentile(values, [25.0, 50.0, 75.0])
595 sigma = 0.741*(hi - lo)
596 good = np.abs(values - med) < rejSigma*sigma
597 if good.sum() == len(good):
598 break
599 values = values[good]
601 calib.coeffs[ii][jj] = np.mean(values)
602 if calib.coeffNum[ii][jj] == 1:
603 calib.coeffErr[ii][jj] = np.nan
604 else:
605 correctionFactor = self.sigmaClipCorrection(rejSigma)
606 calib.coeffErr[ii][jj] = np.std(values) * correctionFactor
607 calib.coeffValid[ii][jj] = (np.abs(calib.coeffs[ii][jj]) >
608 calib.coeffErr[ii][jj] / np.sqrt(calib.coeffNum[ii][jj]))
610 if calib.coeffNum[ii][jj] > 1:
611 self.debugRatios('measure', ratios, ordering[ii], ordering[jj],
612 calib.coeffs[ii][jj], calib.coeffValid[ii][jj])
614 return calib
616 @staticmethod
617 def sigmaClipCorrection(nSigClip):
618 """Correct measured sigma to account for clipping.
620 If we clip our input data and then measure sigma, then the
621 measured sigma is smaller than the true value because real
622 points beyond the clip threshold have been removed. This is a
623 small (1.5% at nSigClip=3) effect when nSigClip >~ 3, but the
624 default parameters for measure crosstalk use nSigClip=2.0.
625 This causes the measured sigma to be about 15% smaller than
626 real. This formula corrects the issue, for the symmetric case
627 (upper clip threshold equal to lower clip threshold).
629 Parameters
630 ----------
631 nSigClip : `float`
632 Number of sigma the measurement was clipped by.
634 Returns
635 -------
636 scaleFactor : `float`
637 Scale factor to increase the measured sigma by.
639 """
640 varFactor = 1.0 + (2 * nSigClip * norm.pdf(nSigClip)) / (norm.cdf(nSigClip) - norm.cdf(-nSigClip))
641 return 1.0 / np.sqrt(varFactor)
643 @staticmethod
644 def filterCrosstalkCalib(inCalib):
645 """Apply valid constraints to the measured values.
647 Any measured coefficient that is determined to be invalid is
648 set to zero, and has the error set to nan. The validation is
649 determined by checking that the measured coefficient is larger
650 than the calculated standard error of the mean.
652 Parameters
653 ----------
654 inCalib : `lsst.ip.isr.CrosstalkCalib`
655 Input calibration to filter.
657 Returns
658 -------
659 outCalib : `lsst.ip.isr.CrosstalkCalib`
660 Filtered calibration.
661 """
662 outCalib = CrosstalkCalib()
663 outCalib.numAmps = inCalib.numAmps
665 outCalib.coeffs = inCalib.coeffs
666 outCalib.coeffs[~inCalib.coeffValid] = 0.0
668 outCalib.coeffErr = inCalib.coeffErr
669 outCalib.coeffErr[~inCalib.coeffValid] = np.nan
671 outCalib.coeffNum = inCalib.coeffNum
672 outCalib.coeffValid = inCalib.coeffValid
674 return outCalib
676 def debugRatios(self, stepname, ratios, i, j, coeff=0.0, valid=False):
677 """Utility function to examine the final CT ratio set.
679 Parameters
680 ----------
681 stepname : `str`
682 State of processing to view.
683 ratios : `dict` of `dict` of `np.ndarray`
684 Array of measured CT ratios, indexed by source/victim
685 amplifier.
686 i : `str`
687 Index of the source amplifier.
688 j : `str`
689 Index of the target amplifier.
690 coeff : `float`, optional
691 Coefficient calculated to plot along with the simple mean.
692 valid : `bool`, optional
693 Validity to be added to the plot title.
694 """
695 frame = getDebugFrame(self._display, stepname)
696 if frame:
697 if i == j or ratios is None or len(ratios) < 1:
698 pass
700 ratioList = ratios[i][j]
701 if ratioList is None or len(ratioList) < 1:
702 pass
704 mean = np.mean(ratioList)
705 std = np.std(ratioList)
706 import matplotlib.pyplot as plt
707 figure = plt.figure(1)
708 figure.clear()
709 plt.hist(x=ratioList, bins=len(ratioList),
710 cumulative=True, color='b', density=True, histtype='step')
711 plt.xlabel("Measured pixel ratio")
712 plt.ylabel(f"CDF: n={len(ratioList)}")
713 plt.xlim(np.percentile(ratioList, [1.0, 99]))
714 plt.axvline(x=mean, color="k")
715 plt.axvline(x=coeff, color='g')
716 plt.axvline(x=(std / np.sqrt(len(ratioList))), color='r')
717 plt.axvline(x=-(std / np.sqrt(len(ratioList))), color='r')
718 plt.title(f"(Source {i} -> Target {j}) mean: {mean:.2g} coeff: {coeff:.2g} valid: {valid}")
719 figure.show()
721 prompt = "Press Enter to continue: "
722 while True:
723 ans = input(prompt).lower()
724 if ans in ("", "c",):
725 break
726 elif ans in ("pdb", "p",):
727 import pdb
728 pdb.set_trace()
729 plt.close()
732class MeasureCrosstalkConfig(Config):
733 extract = ConfigurableField(
734 target=CrosstalkExtractTask,
735 doc="Task to measure pixel ratios.",
736 )
737 solver = ConfigurableField(
738 target=CrosstalkSolveTask,
739 doc="Task to convert ratio lists to crosstalk coefficients.",
740 )
743class MeasureCrosstalkTask(pipeBase.CmdLineTask):
744 """Measure intra-detector crosstalk.
746 See also
747 --------
748 lsst.ip.isr.crosstalk.CrosstalkCalib
749 lsst.cp.pipe.measureCrosstalk.CrosstalkExtractTask
750 lsst.cp.pipe.measureCrosstalk.CrosstalkSolveTask
752 Notes
753 -----
754 The crosstalk this method measures assumes that when a bright
755 pixel is found in one detector amplifier, all other detector
756 amplifiers may see a signal change in the same pixel location
757 (relative to the readout amplifier) as these other pixels are read
758 out at the same time.
760 After processing each input exposure through a limited set of ISR
761 stages, bright unmasked pixels above the threshold are identified.
762 The potential CT signal is found by taking the ratio of the
763 appropriate background-subtracted pixel value on the other
764 amplifiers to the input value on the source amplifier. If the
765 source amplifier has a large number of bright pixels as well, the
766 background level may be elevated, leading to poor ratio
767 measurements.
769 The set of ratios found between each pair of amplifiers across all
770 input exposures is then gathered to produce the final CT
771 coefficients. The sigma-clipped mean and sigma are returned from
772 these sets of ratios, with the coefficient to supply to the ISR
773 CrosstalkTask() being the multiplicative inverse of these values.
775 This Task simply calls the pipetask versions of the measure
776 crosstalk code.
777 """
778 ConfigClass = MeasureCrosstalkConfig
779 _DefaultName = "measureCrosstalk"
781 # Let's use this instead of messing with parseAndRun.
782 RunnerClass = DataRefListRunner
784 def __init__(self, **kwargs):
785 super().__init__(**kwargs)
786 self.makeSubtask("extract")
787 self.makeSubtask("solver")
789 def runDataRef(self, dataRefList):
790 """Run extract task on each of inputs in the dataRef list, then pass
791 that to the solver task.
793 Parameters
794 ----------
795 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
796 Data references for exposures for detectors to process.
798 Returns
799 -------
800 results : `lsst.pipe.base.Struct`
801 The results struct containing:
803 ``outputCrosstalk`` : `lsst.ip.isr.CrosstalkCalib`
804 Final crosstalk calibration.
805 ``outputProvenance`` : `lsst.ip.isr.IsrProvenance`
806 Provenance data for the new calibration.
808 Raises
809 ------
810 RuntimeError
811 Raised if multiple target detectors are supplied.
812 """
813 dataRef = dataRefList[0]
814 camera = dataRef.get("camera")
816 ratios = []
817 activeChip = None
818 for dataRef in dataRefList:
819 exposure = dataRef.get("postISRCCD")
820 if activeChip:
821 if exposure.getDetector().getName() != activeChip:
822 raise RuntimeError("Too many input detectors supplied!")
823 else:
824 activeChip = exposure.getDetector().getName()
826 self.extract.debugView("extract", exposure)
827 result = self.extract.run(exposure)
828 ratios.append(result.outputRatios)
830 finalResults = self.solver.run(ratios, camera=camera)
831 dataRef.put(finalResults.outputCrosstalk, "crosstalk")
833 return finalResults