LSST Applications g0fba68d861+5616995c1c,g1ebb85f214+2420ccdea7,g1fd858c14a+44c57a1f81,g21d47ad084+8e51fce9ac,g262e1987ae+1a7d68eb3b,g2cef7863aa+3bd8df3d95,g35bb328faa+fcb1d3bbc8,g36ff55ed5b+2420ccdea7,g47891489e3+5c6313fe9a,g53246c7159+fcb1d3bbc8,g646c943bdb+dbb9921566,g67b6fd64d1+5c6313fe9a,g6bd32b75b5+2420ccdea7,g74acd417e5+37fc0c974d,g786e29fd12+cf7ec2a62a,g86c591e316+6e13bcb9e9,g87389fa792+1e0a283bba,g89139ef638+5c6313fe9a,g90f42f885a+fce05a46d3,g9125e01d80+fcb1d3bbc8,g93e38de9ac+5345a64125,g95a1e89356+47d08a1cc6,g97be763408+bba861c665,ga9e4eb89a6+85210110a1,gb0b61e0e8e+1f27f70249,gb58c049af0+f03b321e39,gb89ab40317+5c6313fe9a,gc4e39d7843+4e09c98c3d,gd16ba4ae74+5402bcf54a,gd8ff7fe66e+2420ccdea7,gd9a9a58781+fcb1d3bbc8,gdab6d2f7ff+37fc0c974d,gde280f09ee+604b327636,ge278dab8ac+50e2446c94,ge410e46f29+5c6313fe9a,gef3c2e6661+6b480e0fb7,gf67bdafdda+5c6313fe9a,gffca2db377+fcb1d3bbc8,v29.2.0.rc1
LSST Data Management Base Package
Loading...
Searching...
No Matches
isrStatistics.py
Go to the documentation of this file.
1# This file is part of ip_isr.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21
22__all__ = ["IsrStatisticsTaskConfig", "IsrStatisticsTask"]
23
24import numpy as np
25import astropy.stats
26from scipy.signal.windows import hamming, hann, gaussian
27from scipy.signal import butter, filtfilt
28from scipy.stats import linregress
29
30import lsst.afw.math as afwMath
31import lsst.afw.image as afwImage
32import lsst.pipe.base as pipeBase
33import lsst.pex.config as pexConfig
34
35from lsst.afw.cameraGeom import ReadoutCorner
36from .isrFunctions import gainContext, isTrimmedExposure
37
38
39class IsrStatisticsTaskConfig(pexConfig.Config):
40 """Image statistics options.
41 """
42 doCtiStatistics = pexConfig.Field(
43 dtype=bool,
44 doc="Measure CTI statistics from image and overscans?",
45 default=False,
46 )
47 doApplyGainsForCtiStatistics = pexConfig.Field(
48 dtype=bool,
49 doc="Apply gain to the overscan region when measuring CTI statistics?",
50 default=True,
51 # TODO: DM-46721
52 deprecated="This field is no longer used. Will be removed after v28.",
53 )
54
55 doBandingStatistics = pexConfig.Field(
56 dtype=bool,
57 doc="Measure image banding metric?",
58 default=False,
59 )
60 bandingKernelSize = pexConfig.Field(
61 dtype=int,
62 doc="Width of box for boxcar smoothing for banding metric.",
63 default=3,
64 check=lambda x: x == 0 or x % 2 != 0,
65 )
66 bandingFractionLow = pexConfig.Field(
67 dtype=float,
68 doc="Fraction of values to exclude from low samples.",
69 default=0.1,
70 check=lambda x: x >= 0.0 and x <= 1.0
71 )
72 bandingFractionHigh = pexConfig.Field(
73 dtype=float,
74 doc="Fraction of values to exclude from high samples.",
75 default=0.9,
76 check=lambda x: x >= 0.0 and x <= 1.0,
77 )
78 bandingUseHalfDetector = pexConfig.Field(
79 dtype=float,
80 doc="Use only the first half set of amplifiers.",
81 default=True,
82 )
83
84 doProjectionStatistics = pexConfig.Field(
85 dtype=bool,
86 doc="Measure projection metric?",
87 default=False,
88 )
89 projectionKernelSize = pexConfig.Field(
90 dtype=int,
91 doc="Width of box for boxcar smoothing of projections.",
92 default=0,
93 check=lambda x: x == 0 or x % 2 != 0,
94 )
95 doProjectionFft = pexConfig.Field(
96 dtype=bool,
97 doc="Generate FFTs from the image projections?",
98 default=False,
99 )
100 projectionFftWindow = pexConfig.ChoiceField(
101 dtype=str,
102 doc="Type of windowing to use prior to calculating FFT.",
103 default="HAMMING",
104 allowed={
105 "HAMMING": "Hamming window.",
106 "HANN": "Hann window.",
107 "GAUSSIAN": "Gaussian window.",
108 "NONE": "No window."
109 }
110 )
111
112 doDivisaderoStatistics = pexConfig.Field(
113 dtype=bool,
114 doc="Measure divisadero tearing statistics?",
115 default=False,
116 )
117 divisaderoEdgePixels = pexConfig.Field(
118 dtype=int,
119 doc="Number of edge pixels excluded from divisadero linear fit.",
120 default=25,
121 )
122 divisaderoNumImpactPixels = pexConfig.Field(
123 dtype=int,
124 doc="Number of edge pixels to examine for divisadero tearing.",
125 default=2,
126 )
127 divisaderoProjectionMinimum = pexConfig.Field(
128 dtype=int,
129 doc="Minimum row to consider when taking robust mean of columns.",
130 default=10,
131 )
132 divisaderoProjectionMaximum = pexConfig.Field(
133 dtype=int,
134 doc="Maximum row to consider when taking robust mean of columns",
135 default=210,
136 )
137
138 doCopyCalibDistributionStatistics = pexConfig.Field(
139 dtype=bool,
140 doc="Copy calibration distribution statistics to output?",
141 default=False,
142 )
143 expectedDistributionLevels = pexConfig.ListField(
144 dtype=float,
145 doc="Percentile levels expected in the calibration header.",
146 default=[0, 5, 16, 50, 84, 95, 100],
147 )
148
149 doBiasShiftStatistics = pexConfig.Field(
150 dtype=bool,
151 doc="Measure number of image shifts in overscan?",
152 default=False,
153 )
154 biasShiftFilterOrder = pexConfig.Field(
155 dtype=int,
156 doc="Filter order for Butterworth highpass filter.",
157 default=5,
158 )
159 biasShiftCutoff = pexConfig.Field(
160 dtype=float,
161 doc="Cutoff frequency for highpass filter.",
162 default=1.0/15.0,
163 )
164 biasShiftWindow = pexConfig.Field(
165 dtype=int,
166 doc="Filter window size in pixels for highpass filter.",
167 default=30,
168 )
169 biasShiftThreshold = pexConfig.Field(
170 dtype=float,
171 doc="S/N threshold for bias shift detection.",
172 default=3.0,
173 )
174 biasShiftRowSkip = pexConfig.Field(
175 dtype=int,
176 doc="Number of rows to skip for the bias shift detection.",
177 default=30,
178 )
179 biasShiftColumnSkip = pexConfig.Field(
180 dtype=int,
181 doc="Number of columns to skip when averaging the overscan region.",
182 default=3,
183 )
184
185 doAmplifierCorrelationStatistics = pexConfig.Field(
186 dtype=bool,
187 doc="Measure amplifier correlations?",
188 default=False,
189 )
190
191 stat = pexConfig.Field(
192 dtype=str,
193 default="MEANCLIP",
194 doc="Statistic name to use to measure regions.",
195 )
196 nSigmaClip = pexConfig.Field(
197 dtype=float,
198 default=3.0,
199 doc="Clipping threshold for background",
200 )
201 nIter = pexConfig.Field(
202 dtype=int,
203 default=3,
204 doc="Clipping iterations for background",
205 )
206 badMask = pexConfig.ListField(
207 dtype=str,
208 default=["BAD", "INTRP", "SAT"],
209 doc="Mask planes to ignore when identifying source pixels."
210 )
211
212
213class IsrStatisticsTask(pipeBase.Task):
214 """Task to measure arbitrary statistics on ISR processed exposures.
215
216 The goal is to wrap a number of optional measurements that are
217 useful for calibration production and detector stability.
218 """
219 ConfigClass = IsrStatisticsTaskConfig
220 _DefaultName = "isrStatistics"
221
222 def __init__(self, statControl=None, **kwargs):
223 super().__init__(**kwargs)
224 self.statControl = afwMath.StatisticsControl(self.config.nSigmaClip, self.config.nIter,
225 afwImage.Mask.getPlaneBitMask(self.config.badMask))
227
228 def run(self, inputExp, untrimmedInputExposure=None, ptc=None, serialOverscanResults=None,
229 parallelOverscanResults=None, doLegacyCtiStatistics=False, **kwargs):
230 """Task to run arbitrary statistics.
231
232 The statistics should be measured by individual methods, and
233 add to the dictionary in the return struct.
234
235 Parameters
236 ----------
237 inputExp : `lsst.afw.image.Exposure`
238 The exposure to measure.
239 untrimmedInputExp :
240 The exposure to measure overscan statistics from.
241 ptc : `lsst.ip.isr.PtcDataset`, optional
242 A PTC object containing gains to use.
243 serialOverscanResults : `list` [`lsst.pipe.base.Struct`], optional
244 List of serial overscan results. Expected fields are:
245
246 ``imageFit``
247 Value or fit subtracted from the amplifier image data
248 (scalar or `lsst.afw.image.Image`).
249 ``overscanFit``
250 Value or fit subtracted from the overscan image data
251 (scalar or `lsst.afw.image.Image`).
252 ``overscanImage``
253 Image of the overscan region with the overscan
254 correction applied (`lsst.afw.image.Image`). This
255 quantity is used to estimate the amplifier read noise
256 empirically.
257 parallelOverscanResults : `list` [`lsst.pipe.base.Struct`], optional
258 List of parallel overscan results. Expected fields are:
259
260 ``imageFit``
261 Value or fit subtracted from the amplifier image data
262 (scalar or `lsst.afw.image.Image`).
263 ``overscanFit``
264 Value or fit subtracted from the overscan image data
265 (scalar or `lsst.afw.image.Image`).
266 ``overscanImage``
267 Image of the overscan region with the overscan
268 correction applied (`lsst.afw.image.Image`). This
269 quantity is used to estimate the amplifier read noise
270 empirically.
271 doLegacyCtiStatistics : `bool`, optional
272 Use the older version of measureCti (not recommended).
273 This should be True if and only if this task is called
274 from IsrTask. TODO: Deprecate legacy CTI + CTI correction
275 from IsrTask (DM-48757).
276 **kwargs :
277 Keyword arguments. Calibrations being passed in should
278 have an entry here.
279
280 Returns
281 -------
282 resultStruct : `lsst.pipe.base.Struct`
283 Contains the measured statistics as a dict stored in a
284 field named ``results``.
285
286 Raises
287 ------
288 RuntimeError
289 Raised if the amplifier gains could not be found.
290 """
291 # Verify inputs
292 if self.config.doCtiStatistics and not doLegacyCtiStatistics:
293 if untrimmedInputExposure is None:
294 raise RuntimeError("Must pass untrimmed exposure if doCtiStatistics==True.")
295
296 # Find gains.
297 detector = inputExp.getDetector()
298 if ptc is not None:
299 gains = ptc.gain
300 elif detector is not None:
301 gains = {amp.getName(): amp.getGain() for amp in detector.getAmplifiers()}
302 else:
303 raise RuntimeError("No source of gains provided.")
304
305 ctiResults = None
306 if self.config.doCtiStatistics:
307 if doLegacyCtiStatistics:
308 self.log.warning("Calculating the legacy CTI statistics is not recommended.")
309 ctiResults = self.measureCtiLegacy(inputExp, serialOverscanResults, gains)
310 else:
311 ctiResults = self.measureCti(inputExp, untrimmedInputExposure, gains)
312
313 bandingResults = None
314 if self.config.doBandingStatistics:
315 bandingResults = self.measureBanding(inputExp, serialOverscanResults)
316
317 projectionResults = None
318 if self.config.doProjectionStatistics:
319 projectionResults = self.measureProjectionStatistics(inputExp, serialOverscanResults)
320
321 divisaderoResults = None
322 if self.config.doDivisaderoStatistics:
323 divisaderoResults = self.measureDivisaderoStatistics(inputExp, **kwargs)
324
325 calibDistributionResults = None
326 if self.config.doCopyCalibDistributionStatistics:
327 calibDistributionResults = self.copyCalibDistributionStatistics(inputExp, **kwargs)
328
329 biasShiftResults = None
330 if self.config.doBiasShiftStatistics:
331 biasShiftResults = self.measureBiasShifts(inputExp, serialOverscanResults)
332
333 ampCorrelationResults = None
334 if self.config.doAmplifierCorrelationStatistics:
335 ampCorrelationResults = self.measureAmpCorrelations(inputExp, serialOverscanResults)
336
337 mjd = inputExp.getMetadata().get("MJD", None)
338
339 return pipeBase.Struct(
340 results={"CTI": ctiResults,
341 "BANDING": bandingResults,
342 "PROJECTION": projectionResults,
343 "CALIBDIST": calibDistributionResults,
344 "BIASSHIFT": biasShiftResults,
345 "AMPCORR": ampCorrelationResults,
346 "MJD": mjd,
347 'DIVISADERO': divisaderoResults,
348 },
349 )
350
351 def measureCti(self, inputExp, untrimmedInputExp, gains):
352 """Task to measure CTI statistics.
353
354 Parameters
355 ----------
356 inputExp : `lsst.afw.image.Exposure`
357 The exposure to measure.
358 untrimmedInputExp : `lsst.afw.image.Exposure`
359 Exposure to measure overscan from.
360 gains : `dict` [`str` `float`]
361 Dictionary of per-amplifier gains, indexed by amplifier name.
362
363 Returns
364 -------
365 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
366 Dictionary of measurements, keyed by amplifier name and
367 statistics segment. Everything in units based on electron.
368
369 Notes
370 -------
371 The input exposure is needed because it contains the last imaging
372 pixel, with defects applied. And the untrimmed input exposure is
373 needed because it contains the overscan regions. It needs to be
374 this way because the defect masking code requires that the image
375 be trimmed, but we need the image with defects masked to measure
376 the CTI from the last imaging pixel.
377 """
378 outputStats = {}
379
380 detector = inputExp.getDetector()
381 image = inputExp.image
382 untrimmedImage = untrimmedInputExp.image
383
384 # It only makes sense to measure CTI in electron units.
385 # Make it so.
386 imageUnits = inputExp.getMetadata().get("LSST ISR UNITS")
387 untrimmedImageUnits = untrimmedInputExp.getMetadata().get("LSST ISR UNITS")
388
389 # Ensure we always have the same units.
390 applyGainToImage = True if imageUnits == "adu" else False
391 applyGainToUntrimmedImage = True if untrimmedImageUnits == "adu" else False
392
393 with gainContext(inputExp, image, applyGainToImage, gains, isTrimmed=True):
394 with gainContext(untrimmedInputExp, untrimmedImage, applyGainToUntrimmedImage,
395 gains, isTrimmed=False):
396 for amp in detector.getAmplifiers():
397 ampStats = {}
398 readoutCorner = amp.getReadoutCorner()
399
400 # Full data region.
401 dataRegion = image[amp.getBBox()]
402 serialOverscanImage = untrimmedImage[amp.getRawSerialOverscanBBox()]
403 parallelOverscanImage = untrimmedImage[amp.getRawSerialOverscanBBox()]
404
405 # Get the mean of the image
406 ampStats["IMAGE_MEAN"] = afwMath.makeStatistics(dataRegion, self.statType,
407 self.statControl).getValue()
408
409 # First and last image columns.
411 dataRegion.array[:, 0],
412 self.statType,
413 self.statControl,
414 ).getValue()
416 dataRegion.array[:, -1],
417 self.statType,
418 self.statControl,
419 ).getValue()
420
421 # First and last image rows.
423 dataRegion.array[0, :],
424 self.statType,
425 self.statControl,
426 ).getValue()
428 dataRegion.array[-1, :],
429 self.statType,
430 self.statControl,
431 ).getValue()
432
433 # We want these relative to the readout corner. If that's
434 # on the right side, we need to swap them.
435 if readoutCorner in (ReadoutCorner.LR, ReadoutCorner.UR):
436 ampStats["FIRST_COLUMN_MEAN"] = colZ
437 ampStats["LAST_COLUMN_MEAN"] = colA
438 else:
439 ampStats["FIRST_COLUMN_MEAN"] = colA
440 ampStats["LAST_COLUMN_MEAN"] = colZ
441
442 # We want these relative to the readout corner. If that's
443 # on the top, we need to swap them.
444 if readoutCorner in (ReadoutCorner.UR, ReadoutCorner.UL):
445 ampStats["FIRST_ROW_MEAN"] = rowZ
446 ampStats["LAST_ROW_MEAN"] = rowA
447 else:
448 ampStats["FIRST_ROW_MEAN"] = rowA
449 ampStats["LAST_ROW_MEAN"] = rowZ
450
451 # Measure the columns of the serial overscan and the rows
452 # of the parallel overscan.
453 nSerialOverscanCols = amp.getRawSerialOverscanBBox().getWidth()
454 nParallelOverscanRows = amp.getRawParallelOverscanBBox().getHeight()
455
456 # Calculate serial overscan statistics
457 columns = []
458 columnValues = []
459 for idx in range(0, nSerialOverscanCols):
460 # If overscan.doParallelOverscan=True, the
461 # overscanImage will contain both the serial
462 # and parallel overscan regions.
463 serialOverscanColMean = afwMath.makeStatistics(
464 serialOverscanImage.array[:, idx],
465 self.statType,
466 self.statControl,
467 ).getValue()
468 columns.append(idx)
469 # The overscan input is always in adu, but it only
470 # makes sense to measure CTI in electron units.
471 columnValues.append(serialOverscanColMean)
472
473 # We want these relative to the readout corner. If that's
474 # on the right side, we need to swap them.
475 if readoutCorner in (ReadoutCorner.LR, ReadoutCorner.UR):
476 ampStats["SERIAL_OVERSCAN_COLUMNS"] = list(reversed(columns))
477 ampStats["SERIAL_OVERSCAN_VALUES"] = list(reversed(columnValues))
478 else:
479 ampStats["SERIAL_OVERSCAN_COLUMNS"] = columns
480 ampStats["SERIAL_OVERSCAN_VALUES"] = columnValues
481
482 # Calculate parallel overscan statistics
483 rows = []
484 rowValues = []
485 for idx in range(0, nParallelOverscanRows):
486 parallelOverscanRowMean = afwMath.makeStatistics(
487 parallelOverscanImage.array[idx, :],
488 self.statType,
489 self.statControl,
490 ).getValue()
491 rows.append(idx)
492 # The overscan input is always in adu, but it only
493 # makes sense to measure CTI in electron units.
494 rowValues.append(parallelOverscanRowMean)
495
496 # We want these relative to the readout corner. If that's
497 # on the top, we need to swap them.
498 if readoutCorner in (ReadoutCorner.UR, ReadoutCorner.UL):
499 ampStats["PARALLEL_OVERSCAN_ROWS"] = list(reversed(rows))
500 ampStats["PARALLEL_OVERSCAN_VALUES"] = list(reversed(rowValues))
501 else:
502 ampStats["PARALLEL_OVERSCAN_ROWS"] = rows
503 ampStats["PARALLEL_OVERSCAN_VALUES"] = rowValues
504
505 outputStats[amp.getName()] = ampStats
506
507 return outputStats
508
509 def measureCtiLegacy(self, inputExp, serialOverscans, gains):
510 """Task to measure CTI statistics.
511
512 Parameters
513 ----------
514 inputExp : `lsst.afw.image.Exposure`
515 Exposure to measure.
516 serialOverscans : `list` [`lsst.pipe.base.Struct`]
517 List of serial overscan results (expects base units of adu).
518 Expected fields are:
519
520 ``imageFit``
521 Value or fit subtracted from the amplifier image data
522 (scalar or `lsst.afw.image.Image`).
523 ``overscanFit``
524 Value or fit subtracted from the overscan image data
525 (scalar or `lsst.afw.image.Image`).
526 ``overscanImage``
527 Image of the overscan region with the overscan
528 correction applied (`lsst.afw.image.Image`). This
529 quantity is used to estimate the amplifier read noise
530 empirically.
531 gains : `dict` [`str` `float`]
532 Dictionary of per-amplifier gains, indexed by amplifier name.
533
534 Returns
535 -------
536 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
537 Dictionary of measurements, keyed by amplifier name and
538 statistics segment. Everything in units based on electron.
539 """
540 self.log.warning("Using the legacy version of CTI statistics may give wrong CTI")
541
542 outputStats = {}
543
544 detector = inputExp.getDetector()
545 image = inputExp.image
546
547 # It only makes sense to measure CTI in electron units.
548 # Make it so.
549 imageUnits = inputExp.getMetadata().get("LSST ISR UNITS")
550 applyGain = False
551 if imageUnits == "adu":
552 applyGain = True
553
554 # Check if the image is trimmed.
555 isTrimmed = isTrimmedExposure(inputExp)
556
557 # Ensure we have the same number of overscans as amplifiers.
558 assert len(serialOverscans) == len(detector.getAmplifiers())
559
560 with gainContext(inputExp, image, applyGain, gains, isTrimmed=isTrimmed):
561 for ampIter, amp in enumerate(detector.getAmplifiers()):
562 ampStats = {}
563 readoutCorner = amp.getReadoutCorner()
564
565 # Full data region.
566 dataRegion = image[amp.getBBox() if isTrimmed else amp.getRawDataBBox()]
567
568 # Get the mean of the image
569 ampStats["IMAGE_MEAN"] = afwMath.makeStatistics(dataRegion, self.statType,
570 self.statControl).getValue()
571
572 # First and last image columns.
573 pixelA = afwMath.makeStatistics(dataRegion.array[:, 0],
574 self.statType,
575 self.statControl).getValue()
576 pixelZ = afwMath.makeStatistics(dataRegion.array[:, -1],
577 self.statType,
578 self.statControl).getValue()
579
580 # We want these relative to the readout corner. If that's
581 # on the right side, we need to swap them.
582 if readoutCorner in (ReadoutCorner.LR, ReadoutCorner.UR):
583 ampStats["FIRST_MEAN"] = pixelZ
584 ampStats["LAST_MEAN"] = pixelA
585 else:
586 ampStats["FIRST_MEAN"] = pixelA
587 ampStats["LAST_MEAN"] = pixelZ
588
589 # Measure the columns of the overscan.
590 if serialOverscans[ampIter] is None:
591 # The amplifier is likely entirely bad, and needs to
592 # be skipped.
593 self.log.warning("No overscan information available for ISR statistics for amp %s.",
594 amp.getName())
595 nCols = amp.getRawSerialOverscanBBox().getWidth()
596 ampStats["OVERSCAN_COLUMNS"] = np.full((nCols, ), np.nan)
597 ampStats["OVERSCAN_VALUES"] = np.full((nCols, ), np.nan)
598 else:
599 overscanImage = serialOverscans[ampIter].overscanImage
600
601 columns = []
602 values = []
603 for column in range(0, overscanImage.getWidth()):
604 # If overscan.doParallelOverscan=True, the
605 # overscanImage will contain both the serial
606 # and parallel overscan regions.
607 # Only the serial CTI correction has been
608 # implemented, so we must select only the
609 # serial overscan rows for a given column.
610 nRows = amp.getRawSerialOverscanBBox().getHeight()
611 osMean = afwMath.makeStatistics(overscanImage.image.array[:nRows, column],
612 self.statType, self.statControl).getValue()
613 columns.append(column)
614 # The overscan input is always in adu, but it only
615 # makes sense to measure CTI in electron units.
616 values.append(osMean * gains[amp.getName()])
617
618 # We want these relative to the readout corner. If that's
619 # on the right side, we need to swap them.
620 if readoutCorner in (ReadoutCorner.LR, ReadoutCorner.UR):
621 ampStats["OVERSCAN_COLUMNS"] = list(reversed(columns))
622 ampStats["OVERSCAN_VALUES"] = list(reversed(values))
623 else:
624 ampStats["OVERSCAN_COLUMNS"] = columns
625 ampStats["OVERSCAN_VALUES"] = values
626
627 outputStats[amp.getName()] = ampStats
628
629 return outputStats
630
631 @staticmethod
632 def makeKernel(kernelSize):
633 """Make a boxcar smoothing kernel.
634
635 Parameters
636 ----------
637 kernelSize : `int`
638 Size of the kernel in pixels.
639
640 Returns
641 -------
642 kernel : `np.array`
643 Kernel for boxcar smoothing.
644 """
645 if kernelSize > 0:
646 kernel = np.full(kernelSize, 1.0 / kernelSize)
647 else:
648 kernel = np.array([1.0])
649 return kernel
650
651 def measureBanding(self, inputExp, overscans):
652 """Task to measure banding statistics.
653
654 Parameters
655 ----------
656 inputExp : `lsst.afw.image.Exposure`
657 Exposure to measure.
658 overscans : `list` [`lsst.pipe.base.Struct`]
659 List of overscan results. Expected fields are:
660
661 ``imageFit``
662 Value or fit subtracted from the amplifier image data
663 (scalar or `lsst.afw.image.Image`).
664 ``overscanFit``
665 Value or fit subtracted from the overscan image data
666 (scalar or `lsst.afw.image.Image`).
667 ``overscanImage``
668 Image of the overscan region with the overscan
669 correction applied (`lsst.afw.image.Image`). This
670 quantity is used to estimate the amplifier read noise
671 empirically.
672
673 Returns
674 -------
675 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
676 Dictionary of measurements, keyed by amplifier name and
677 statistics segment.
678 """
679 outputStats = {}
680
681 detector = inputExp.getDetector()
682 kernel = self.makeKernel(self.config.bandingKernelSize)
683
684 outputStats["AMP_BANDING"] = []
685 for amp, overscanData in zip(detector.getAmplifiers(), overscans):
686 overscanFit = np.array(overscanData.overscanFit)
687 overscanArray = overscanData.overscanImage.image.array
688 rawOverscan = np.mean(overscanArray + overscanFit, axis=1)
689
690 smoothedOverscan = np.convolve(rawOverscan, kernel, mode="valid")
691
692 low, high = np.quantile(smoothedOverscan, [self.config.bandingFractionLow,
693 self.config.bandingFractionHigh])
694 outputStats["AMP_BANDING"].append(float(high - low))
695
696 if self.config.bandingUseHalfDetector:
697 fullLength = len(outputStats["AMP_BANDING"])
698 outputStats["DET_BANDING"] = float(np.nanmedian(outputStats["AMP_BANDING"][0:fullLength//2]))
699 else:
700 outputStats["DET_BANDING"] = float(np.nanmedian(outputStats["AMP_BANDING"]))
701
702 return outputStats
703
704 def measureProjectionStatistics(self, inputExp, overscans):
705 """Task to measure metrics from image slicing.
706
707 Parameters
708 ----------
709 inputExp : `lsst.afw.image.Exposure`
710 Exposure to measure.
711 overscans : `list` [`lsst.pipe.base.Struct`]
712 List of overscan results. Expected fields are:
713
714 ``imageFit``
715 Value or fit subtracted from the amplifier image data
716 (scalar or `lsst.afw.image.Image`).
717 ``overscanFit``
718 Value or fit subtracted from the overscan image data
719 (scalar or `lsst.afw.image.Image`).
720 ``overscanImage``
721 Image of the overscan region with the overscan
722 correction applied (`lsst.afw.image.Image`). This
723 quantity is used to estimate the amplifier read noise
724 empirically.
725
726 Returns
727 -------
728 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
729 Dictionary of measurements, keyed by amplifier name and
730 statistics segment.
731 """
732 outputStats = {}
733
734 detector = inputExp.getDetector()
735 kernel = self.makeKernel(self.config.projectionKernelSize)
736
737 outputStats["AMP_VPROJECTION"] = {}
738 outputStats["AMP_HPROJECTION"] = {}
739 convolveMode = "valid"
740 if self.config.doProjectionFft:
741 outputStats["AMP_VFFT_REAL"] = {}
742 outputStats["AMP_VFFT_IMAG"] = {}
743 outputStats["AMP_HFFT_REAL"] = {}
744 outputStats["AMP_HFFT_IMAG"] = {}
745 convolveMode = "same"
746
747 for amp in detector.getAmplifiers():
748 ampArray = inputExp.image[amp.getBBox()].array
749
750 horizontalProjection = np.mean(ampArray, axis=0)
751 verticalProjection = np.mean(ampArray, axis=1)
752
753 horizontalProjection = np.convolve(horizontalProjection, kernel, mode=convolveMode)
754 verticalProjection = np.convolve(verticalProjection, kernel, mode=convolveMode)
755
756 outputStats["AMP_HPROJECTION"][amp.getName()] = horizontalProjection.tolist()
757 outputStats["AMP_VPROJECTION"][amp.getName()] = verticalProjection.tolist()
758
759 if self.config.doProjectionFft:
760 horizontalWindow = np.ones_like(horizontalProjection)
761 verticalWindow = np.ones_like(verticalProjection)
762 if self.config.projectionFftWindow == "NONE":
763 pass
764 elif self.config.projectionFftWindow == "HAMMING":
765 horizontalWindow = hamming(len(horizontalProjection))
766 verticalWindow = hamming(len(verticalProjection))
767 elif self.config.projectionFftWindow == "HANN":
768 horizontalWindow = hann(len(horizontalProjection))
769 verticalWindow = hann(len(verticalProjection))
770 elif self.config.projectionFftWindow == "GAUSSIAN":
771 horizontalWindow = gaussian(len(horizontalProjection))
772 verticalWindow = gaussian(len(verticalProjection))
773 else:
774 raise RuntimeError(f"Invalid window function: {self.config.projectionFftWindow}")
775
776 horizontalFFT = np.fft.rfft(np.multiply(horizontalProjection, horizontalWindow))
777 verticalFFT = np.fft.rfft(np.multiply(verticalProjection, verticalWindow))
778
779 outputStats["AMP_HFFT_REAL"][amp.getName()] = np.real(horizontalFFT).tolist()
780 outputStats["AMP_HFFT_IMAG"][amp.getName()] = np.imag(horizontalFFT).tolist()
781 outputStats["AMP_VFFT_REAL"][amp.getName()] = np.real(verticalFFT).tolist()
782 outputStats["AMP_VFFT_IMAG"][amp.getName()] = np.imag(verticalFFT).tolist()
783
784 return outputStats
785
786 def copyCalibDistributionStatistics(self, inputExp, **kwargs):
787 """Copy calibration statistics for this exposure.
788
789 Parameters
790 ----------
791 inputExp : `lsst.afw.image.Exposure`
792 The exposure being processed.
793 **kwargs :
794 Keyword arguments with calibrations.
795
796 Returns
797 -------
798 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
799 Dictionary of measurements, keyed by amplifier name and
800 statistics segment.
801 """
802 outputStats = {}
803
804 # Amp level elements
805 for amp in inputExp.getDetector():
806 ampStats = {}
807
808 for calibType in ("bias", "dark", "flat"):
809 if kwargs.get(calibType, None) is not None:
810 metadata = kwargs[calibType].getMetadata()
811 for pct in self.config.expectedDistributionLevels:
812 key = f"LSST CALIB {calibType.upper()} {amp.getName()} DISTRIBUTION {pct}-PCT"
813 ampStats[key] = metadata.get(key, np.nan)
814
815 for calibType in ("defects"):
816 if kwargs.get(calibType, None) is not None:
817 metadata = kwargs[calibType].getMetadata()
818 for key in (f"LSST CALIB {calibType.upper()} {amp.getName()} N_HOT",
819 f"LSST CALIB {calibType.upper()} {amp.getName()} N_COLD"):
820 ampStats[key] = metadata.get(key, np.nan)
821 outputStats[amp.getName()] = ampStats
822
823 # Detector level elements
824 for calibType in ("defects"):
825 if kwargs.get(calibType, None) is not None:
826 metadata = kwargs[calibType].getMetadata()
827 for key in (f"LSST CALIB {calibType.upper()} N_BAD_COLUMNS"):
828 outputStats["detector"][key] = metadata.get(key, np.nan)
829
830 return outputStats
831
832 def measureBiasShifts(self, inputExp, serialOverscanResults):
833 """Measure number of bias shifts from overscan data.
834
835 Parameters
836 ----------
837 inputExp : `lsst.afw.image.Exposure`
838 Exposure to measure.
839 overscans : `list` [`lsst.pipe.base.Struct`]
840 List of overscan results. Expected fields are:
841
842 ``imageFit``
843 Value or fit subtracted from the amplifier image data
844 (scalar or `lsst.afw.image.Image`).
845 ``overscanFit``
846 Value or fit subtracted from the overscan image data
847 (scalar or `lsst.afw.image.Image`).
848 ``overscanImage``
849 Image of the overscan region with the overscan
850 correction applied (`lsst.afw.image.Image`). This
851 quantity is used to estimate the amplifier read noise
852 empirically.
853
854 Returns
855 -------
856 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
857 Dictionary of measurements, keyed by amplifier name and
858 statistics segment.
859
860 Notes
861 -----
862 Based on eop_pipe implementation:
863 https://github.com/lsst-camera-dh/eo_pipe/blob/main/python/lsst/eo/pipe/biasShiftsTask.py # noqa: E501 W505
864 """
865 outputStats = {}
866
867 detector = inputExp.getDetector()
868 for amp, overscans in zip(detector, serialOverscanResults):
869 ampStats = {}
870 # Add fit back to data
871 rawOverscan = overscans.overscanImage.image.array + overscans.overscanFit
872
873 # Collapse array, skipping first three columns
874 rawOverscan = np.mean(rawOverscan[:, self.config.biasShiftColumnSkip:], axis=1)
875
876 # Scan for shifts
877 noise, shift_peaks = self._scan_for_shifts(rawOverscan)
878 ampStats["LOCAL_NOISE"] = float(noise)
879 ampStats["BIAS_SHIFTS"] = shift_peaks
880
881 outputStats[amp.getName()] = ampStats
882 return outputStats
883
884 def _scan_for_shifts(self, overscanData):
885 """Scan overscan data for shifts.
886
887 Parameters
888 ----------
889 overscanData : `list` [`float`]
890 Overscan data to search for shifts.
891
892 Returns
893 -------
894 noise : `float`
895 Noise estimated from Butterworth filtered overscan data.
896 peaks : `list` [`float`, `float`, `int`, `int`]
897 Shift peak information, containing the convolved peak
898 value, the raw peak value, and the lower and upper bounds
899 of the region checked.
900 """
901 numerator, denominator = butter(self.config.biasShiftFilterOrder,
902 self.config.biasShiftCutoff,
903 btype="high", analog=False)
904 noise = np.std(filtfilt(numerator, denominator, overscanData))
905 kernel = np.concatenate([np.arange(self.config.biasShiftWindow),
906 np.arange(-self.config.biasShiftWindow + 1, 0)])
907 kernel = kernel/np.sum(kernel[:self.config.biasShiftWindow])
908
909 convolved = np.convolve(overscanData, kernel, mode="valid")
910 convolved = np.pad(convolved, (self.config.biasShiftWindow - 1, self.config.biasShiftWindow))
911
912 shift_check = np.abs(convolved)/noise
913 shift_mask = shift_check > self.config.biasShiftThreshold
914 shift_mask[:self.config.biasShiftRowSkip] = False
915
916 shift_regions = np.flatnonzero(np.diff(np.r_[np.int8(0),
917 shift_mask.view(np.int8),
918 np.int8(0)])).reshape(-1, 2)
919 shift_peaks = []
920 for region in shift_regions:
921 region_peak = np.argmax(shift_check[region[0]:region[1]]) + region[0]
922 if self._satisfies_flatness(region_peak, convolved[region_peak], overscanData):
923 shift_peaks.append(
924 [float(convolved[region_peak]), float(region_peak),
925 int(region[0]), int(region[1])])
926 return noise, shift_peaks
927
928 def _satisfies_flatness(self, shiftRow, shiftPeak, overscanData):
929 """Determine if a region is flat.
930
931 Parameters
932 ----------
933 shiftRow : `int`
934 Row with possible peak.
935 shiftPeak : `float`
936 Value at the possible peak.
937 overscanData : `list` [`float`]
938 Overscan data used to fit around the possible peak.
939
940 Returns
941 -------
942 isFlat : `bool`
943 Indicates if the region is flat, and so the peak is valid.
944 """
945 prerange = np.arange(shiftRow - self.config.biasShiftWindow, shiftRow)
946 postrange = np.arange(shiftRow, shiftRow + self.config.biasShiftWindow)
947
948 preFit = linregress(prerange, overscanData[prerange])
949 postFit = linregress(postrange, overscanData[postrange])
950
951 if shiftPeak > 0:
952 preTrend = (2*preFit[0]*len(prerange) < shiftPeak)
953 postTrend = (2*postFit[0]*len(postrange) < shiftPeak)
954 else:
955 preTrend = (2*preFit[0]*len(prerange) > shiftPeak)
956 postTrend = (2*postFit[0]*len(postrange) > shiftPeak)
957
958 return (preTrend and postTrend)
959
960 def measureAmpCorrelations(self, inputExp, serialOverscanResults):
961 """Measure correlations between amplifier segments.
962
963 Parameters
964 ----------
965 inputExp : `lsst.afw.image.Exposure`
966 Exposure to measure.
967 overscans : `list` [`lsst.pipe.base.Struct`]
968 List of overscan results. Expected fields are:
969
970 ``imageFit``
971 Value or fit subtracted from the amplifier image data
972 (scalar or `lsst.afw.image.Image`).
973 ``overscanFit``
974 Value or fit subtracted from the overscan image data
975 (scalar or `lsst.afw.image.Image`).
976 ``overscanImage``
977 Image of the overscan region with the overscan
978 correction applied (`lsst.afw.image.Image`). This
979 quantity is used to estimate the amplifier read noise
980 empirically.
981
982 Returns
983 -------
984 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
985 Dictionary of measurements, keyed by amplifier name and
986 statistics segment.
987
988 Notes
989 -----
990 Based on eo_pipe implementation:
991 https://github.com/lsst-camera-dh/eo_pipe/blob/main/python/lsst/eo/pipe/raft_level_correlations.py # noqa: E501 W505
992 """
993 detector = inputExp.getDetector()
994 rows = []
995
996 for ampId, overscan in enumerate(serialOverscanResults):
997 rawOverscan = overscan.overscanImage.image.array + overscan.overscanFit
998 rawOverscan = rawOverscan.ravel()
999
1000 ampImage = inputExp[detector[ampId].getBBox()]
1001 ampImage = ampImage.image.array.ravel()
1002
1003 for ampId2, overscan2 in enumerate(serialOverscanResults):
1004 osC = 1.0
1005 imC = 1.0
1006 if ampId2 != ampId:
1007 rawOverscan2 = overscan2.overscanImage.image.array + overscan2.overscanFit
1008 rawOverscan2 = rawOverscan2.ravel()
1009
1010 osC = np.corrcoef(rawOverscan, rawOverscan2)[0, 1]
1011
1012 ampImage2 = inputExp[detector[ampId2].getBBox()]
1013 ampImage2 = ampImage2.image.array.ravel()
1014
1015 imC = np.corrcoef(ampImage, ampImage2)[0, 1]
1016 rows.append(
1017 {'detector': detector.getId(),
1018 'detectorComp': detector.getId(),
1019 'ampName': detector[ampId].getName(),
1020 'ampComp': detector[ampId2].getName(),
1021 'imageCorr': float(imC),
1022 'overscanCorr': float(osC),
1023 }
1024 )
1025 return rows
1026
1027 def measureDivisaderoStatistics(self, inputExp, **kwargs):
1028 """Measure Max Divisadero Tearing effect per amp.
1029
1030 Parameters
1031 ----------
1032 inputExp : `lsst.afw.image.Exposure`
1033 Exposure to measure. Usually a flat.
1034 **kwargs :
1035 The flat will be selected from here.
1036
1037 Returns
1038 -------
1039 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
1040 Dictionary of measurements, keyed by amplifier name and
1041 statistics segment.
1042 Measurements include
1043
1044 - DIVISADERO_PROFILE: Robust mean of rows between
1045 divisaderoProjection<Maximum|Minumum> on readout edge of ccd
1046 normalized by a linear fit to the same rows.
1047 - DIVISADERO_MAX_PAIR: Tuple of maximum of the absolute values of
1048 the DIVISADERO_PROFILE, for number of pixels (specified by
1049 divisaderoNumImpactPixels on left and right side of amp.
1050 - DIVISADERO_MAX: Maximum of the absolute values of the
1051 the DIVISADERO_PROFILE, for the divisaderoNumImpactPixels on
1052 boundaries of neighboring amps (including the pixels in those
1053 neighborboring amps).
1054 """
1055 outputStats = {}
1056
1057 for amp in inputExp.getDetector():
1058 # Copy unneeded if we do not ever modify the array by flipping
1059 ampArray = inputExp.image[amp.getBBox()].array
1060 # slice the outer top or bottom of the amp: the readout side
1061 if amp.getReadoutCorner().name in ('UL', 'UR'):
1062 minRow = amp.getBBox().getHeight() - self.config.divisaderoProjectionMaximum
1063 maxRow = amp.getBBox().getHeight() - self.config.divisaderoProjectionMinimum
1064 else:
1065 minRow = self.config.divisaderoProjectionMinimum
1066 maxRow = self.config.divisaderoProjectionMaximum
1067
1068 segment = slice(minRow, maxRow)
1069 projection, _, _ = astropy.stats.sigma_clipped_stats(ampArray[segment, :], axis=0)
1070
1071 ampStats = {}
1072 projection /= np.median(projection)
1073 columns = np.arange(len(projection))
1074
1075 segment = slice(self.config.divisaderoEdgePixels, -self.config.divisaderoEdgePixels)
1076 model = np.polyfit(columns[segment], projection[segment], 1)
1077 modelProjection = model[0] * columns + model[1]
1078 divisaderoProfile = projection / modelProjection
1079
1080 # look for max at the edges:
1081 leftMax = np.nanmax(np.abs(divisaderoProfile[0:self.config.divisaderoNumImpactPixels] - 1.0))
1082 rightMax = np.nanmax(np.abs(divisaderoProfile[-self.config.divisaderoNumImpactPixels:] - 1.0))
1083
1084 ampStats['DIVISADERO_PROFILE'] = np.array(divisaderoProfile).tolist()
1085 ampStats['DIVISADERO_MAX_PAIR'] = [leftMax, rightMax]
1086 outputStats[amp.getName()] = ampStats
1087
1088 detector = inputExp.getDetector()
1089 xCenters = [amp.getBBox().getCenterX() for amp in detector]
1090 yCenters = [amp.getBBox().getCenterY() for amp in detector]
1091 xIndices = np.ceil(xCenters / np.min(xCenters) / 2).astype(int) - 1
1092 yIndices = np.ceil(yCenters / np.min(yCenters) / 2).astype(int) - 1
1093 ampIds = np.zeros((len(set(yIndices)), len(set(xIndices))), dtype=int)
1094 for ampId, xIndex, yIndex in zip(np.arange(len(detector)), xIndices, yIndices):
1095 ampIds[yIndex, xIndex] = ampId
1096
1097 # Loop over amps again because the DIVISIDERO_MAX will be the max
1098 # of the profile on its boundary with its neighboring amps
1099 for i, amp in enumerate(detector):
1100 y, x = np.where(ampIds == i)
1101 end = ampIds.shape[1] - 1
1102 xInd = x[0]
1103 yInd = y[0]
1104 thisAmpsPair = outputStats[amp.getName()]['DIVISADERO_MAX_PAIR']
1105
1106 if x == 0:
1107 # leftmost amp: take the max of your right side and
1108 myMax = thisAmpsPair[1]
1109 # your neighbor's left side
1110 neighborMax = outputStats[detector[ampIds[yInd, 1]].getName()]['DIVISADERO_MAX_PAIR'][0]
1111 elif x == end:
1112 # rightmost amp: take the max of your left side and
1113 myMax = thisAmpsPair[0]
1114 # your neighbor's right side
1115 neighborMax = outputStats[detector[ampIds[yInd, end - 1]].getName()]['DIVISADERO_MAX_PAIR'][1]
1116 else:
1117 # Middle amp: take the max of both your own sides and the
1118 myMax = max(thisAmpsPair)
1119 leftName = detector[ampIds[yInd, max(xInd - 1, 0)]].getName()
1120 rightName = detector[ampIds[yInd, min(xInd + 1, ampIds.shape[1] - 1)]].getName()
1121 # right side of the neighbor to your left
1122 # and left side of your neighbor to your right
1123 neighborMax = max(outputStats[leftName]['DIVISADERO_MAX_PAIR'][1],
1124 outputStats[rightName]['DIVISADERO_MAX_PAIR'][0])
1125
1126 divisaderoMax = max([myMax, neighborMax])
1127 outputStats[amp.getName()]['DIVISADERO_MAX'] = divisaderoMax
1128
1129 return outputStats
Pass parameters to a Statistics object.
Definition Statistics.h:83
measureDivisaderoStatistics(self, inputExp, **kwargs)
measureBiasShifts(self, inputExp, serialOverscanResults)
measureCti(self, inputExp, untrimmedInputExp, gains)
measureAmpCorrelations(self, inputExp, serialOverscanResults)
copyCalibDistributionStatistics(self, inputExp, **kwargs)
measureProjectionStatistics(self, inputExp, overscans)
run(self, inputExp, untrimmedInputExposure=None, ptc=None, serialOverscanResults=None, parallelOverscanResults=None, doLegacyCtiStatistics=False, **kwargs)
_satisfies_flatness(self, shiftRow, shiftPeak, overscanData)
measureCtiLegacy(self, inputExp, serialOverscans, gains)
__init__(self, statControl=None, **kwargs)
Statistics makeStatistics(lsst::afw::image::Image< Pixel > const &img, lsst::afw::image::Mask< image::MaskPixel > const &msk, int const flags, StatisticsControl const &sctrl=StatisticsControl())
Handle a watered-down front-end to the constructor (no variance)
Definition Statistics.h:361
Property stringToStatisticsProperty(std::string const property)
Conversion function to switch a string to a Property (see Statistics.h)
gainContext(exp, image, apply, gains=None, invert=False, isTrimmed=True)