228 def run(self, inputExp, untrimmedInputExposure=None, ptc=None, serialOverscanResults=None,
229 parallelOverscanResults=None, doLegacyCtiStatistics=False, **kwargs):
230 """Task to run arbitrary statistics.
232 The statistics should be measured by individual methods, and
233 add to the dictionary in the return struct.
237 inputExp : `lsst.afw.image.Exposure`
238 The exposure to measure.
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:
247 Value or fit subtracted from the amplifier image data
248 (scalar or `lsst.afw.image.Image`).
250 Value or fit subtracted from the overscan image data
251 (scalar or `lsst.afw.image.Image`).
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
257 parallelOverscanResults : `list` [`lsst.pipe.base.Struct`], optional
258 List of parallel overscan results. Expected fields are:
261 Value or fit subtracted from the amplifier image data
262 (scalar or `lsst.afw.image.Image`).
264 Value or fit subtracted from the overscan image data
265 (scalar or `lsst.afw.image.Image`).
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
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).
277 Keyword arguments. Calibrations being passed in should
282 resultStruct : `lsst.pipe.base.Struct`
283 Contains the measured statistics as a dict stored in a
284 field named ``results``.
289 Raised if the amplifier gains could not be found.
292 if self.config.doCtiStatistics
and not doLegacyCtiStatistics:
293 if untrimmedInputExposure
is None:
294 raise RuntimeError(
"Must pass untrimmed exposure if doCtiStatistics==True.")
297 detector = inputExp.getDetector()
300 elif detector
is not None:
301 gains = {amp.getName(): amp.getGain()
for amp
in detector.getAmplifiers()}
303 raise RuntimeError(
"No source of gains provided.")
306 if self.config.doCtiStatistics:
307 if doLegacyCtiStatistics:
308 self.log.warning(
"Calculating the legacy CTI statistics is not recommended.")
311 ctiResults = self.
measureCti(inputExp, untrimmedInputExposure, gains)
313 bandingResults =
None
314 if self.config.doBandingStatistics:
315 bandingResults = self.
measureBanding(inputExp, serialOverscanResults)
317 projectionResults =
None
318 if self.config.doProjectionStatistics:
321 divisaderoResults =
None
322 if self.config.doDivisaderoStatistics:
325 calibDistributionResults =
None
326 if self.config.doCopyCalibDistributionStatistics:
329 biasShiftResults =
None
330 if self.config.doBiasShiftStatistics:
333 ampCorrelationResults =
None
334 if self.config.doAmplifierCorrelationStatistics:
337 mjd = inputExp.getMetadata().get(
"MJD",
None)
339 return pipeBase.Struct(
340 results={
"CTI": ctiResults,
341 "BANDING": bandingResults,
342 "PROJECTION": projectionResults,
343 "CALIBDIST": calibDistributionResults,
344 "BIASSHIFT": biasShiftResults,
345 "AMPCORR": ampCorrelationResults,
347 'DIVISADERO': divisaderoResults,
352 """Task to measure CTI statistics.
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.
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.
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.
380 detector = inputExp.getDetector()
381 image = inputExp.image
382 untrimmedImage = untrimmedInputExp.image
386 imageUnits = inputExp.getMetadata().get(
"LSST ISR UNITS")
387 untrimmedImageUnits = untrimmedInputExp.getMetadata().get(
"LSST ISR UNITS")
390 applyGainToImage =
True if imageUnits ==
"adu" else False
391 applyGainToUntrimmedImage =
True if untrimmedImageUnits ==
"adu" else False
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():
398 readoutCorner = amp.getReadoutCorner()
400 ampStats[
"INPUT_GAIN"] = float(gains[amp.getName()])
401 turnoff = inputExp.metadata[f
"LSST ISR PTCTURNOFF {amp.getName()}"]
402 ampStats[
"INPUT_PTCTURNOFF"] = float(turnoff)
405 dataRegion = image[amp.getBBox()]
406 serialOverscanImage = untrimmedImage[amp.getRawSerialOverscanBBox()]
407 parallelOverscanImage = untrimmedImage[amp.getRawSerialOverscanBBox()]
410 ampStats[
"IMAGE_MEAN"] = afwMath.makeStatistics(dataRegion, self.
statType,
414 colA = afwMath.makeStatistics(
415 dataRegion.array[:, 0],
419 colZ = afwMath.makeStatistics(
420 dataRegion.array[:, -1],
426 rowA = afwMath.makeStatistics(
427 dataRegion.array[0, :],
431 rowZ = afwMath.makeStatistics(
432 dataRegion.array[-1, :],
439 if readoutCorner
in (ReadoutCorner.LR, ReadoutCorner.UR):
440 ampStats[
"FIRST_COLUMN_MEAN"] = colZ
441 ampStats[
"LAST_COLUMN_MEAN"] = colA
443 ampStats[
"FIRST_COLUMN_MEAN"] = colA
444 ampStats[
"LAST_COLUMN_MEAN"] = colZ
448 if readoutCorner
in (ReadoutCorner.UR, ReadoutCorner.UL):
449 ampStats[
"FIRST_ROW_MEAN"] = rowZ
450 ampStats[
"LAST_ROW_MEAN"] = rowA
452 ampStats[
"FIRST_ROW_MEAN"] = rowA
453 ampStats[
"LAST_ROW_MEAN"] = rowZ
457 nSerialOverscanCols = amp.getRawSerialOverscanBBox().getWidth()
458 nParallelOverscanRows = amp.getRawParallelOverscanBBox().getHeight()
463 for idx
in range(0, nSerialOverscanCols):
467 serialOverscanColMean = afwMath.makeStatistics(
468 serialOverscanImage.array[:, idx],
475 columnValues.append(serialOverscanColMean)
479 if readoutCorner
in (ReadoutCorner.LR, ReadoutCorner.UR):
480 ampStats[
"SERIAL_OVERSCAN_COLUMNS"] = list(reversed(columns))
481 ampStats[
"SERIAL_OVERSCAN_VALUES"] = list(reversed(columnValues))
483 ampStats[
"SERIAL_OVERSCAN_COLUMNS"] = columns
484 ampStats[
"SERIAL_OVERSCAN_VALUES"] = columnValues
489 for idx
in range(0, nParallelOverscanRows):
490 parallelOverscanRowMean = afwMath.makeStatistics(
491 parallelOverscanImage.array[idx, :],
498 rowValues.append(parallelOverscanRowMean)
502 if readoutCorner
in (ReadoutCorner.UR, ReadoutCorner.UL):
503 ampStats[
"PARALLEL_OVERSCAN_ROWS"] = list(reversed(rows))
504 ampStats[
"PARALLEL_OVERSCAN_VALUES"] = list(reversed(rowValues))
506 ampStats[
"PARALLEL_OVERSCAN_ROWS"] = rows
507 ampStats[
"PARALLEL_OVERSCAN_VALUES"] = rowValues
509 outputStats[amp.getName()] = ampStats
514 """Task to measure CTI statistics.
518 inputExp : `lsst.afw.image.Exposure`
520 serialOverscans : `list` [`lsst.pipe.base.Struct`]
521 List of serial overscan results (expects base units of adu).
525 Value or fit subtracted from the amplifier image data
526 (scalar or `lsst.afw.image.Image`).
528 Value or fit subtracted from the overscan image data
529 (scalar or `lsst.afw.image.Image`).
531 Image of the overscan region with the overscan
532 correction applied (`lsst.afw.image.Image`). This
533 quantity is used to estimate the amplifier read noise
535 gains : `dict` [`str` `float`]
536 Dictionary of per-amplifier gains, indexed by amplifier name.
540 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
541 Dictionary of measurements, keyed by amplifier name and
542 statistics segment. Everything in units based on electron.
544 self.log.warning(
"Using the legacy version of CTI statistics may give wrong CTI")
548 detector = inputExp.getDetector()
549 image = inputExp.image
553 imageUnits = inputExp.getMetadata().get(
"LSST ISR UNITS")
555 if imageUnits ==
"adu":
562 assert len(serialOverscans) == len(detector.getAmplifiers())
564 with gainContext(inputExp, image, applyGain, gains, isTrimmed=isTrimmed):
565 for ampIter, amp
in enumerate(detector.getAmplifiers()):
567 readoutCorner = amp.getReadoutCorner()
569 ampStats[
"INPUT_GAIN"] = float(gains[amp.getName()])
570 turnoff = inputExp.metadata[f
"LSST ISR PTCTURNOFF {amp.getName()}"]
571 ampStats[
"INPUT_PTCTURNOFF"] = float(turnoff)
574 dataRegion = image[amp.getBBox()
if isTrimmed
else amp.getRawDataBBox()]
577 ampStats[
"IMAGE_MEAN"] = afwMath.makeStatistics(dataRegion, self.
statType,
581 pixelA = afwMath.makeStatistics(dataRegion.array[:, 0],
584 pixelZ = afwMath.makeStatistics(dataRegion.array[:, -1],
590 if readoutCorner
in (ReadoutCorner.LR, ReadoutCorner.UR):
591 ampStats[
"FIRST_MEAN"] = pixelZ
592 ampStats[
"LAST_MEAN"] = pixelA
594 ampStats[
"FIRST_MEAN"] = pixelA
595 ampStats[
"LAST_MEAN"] = pixelZ
598 if serialOverscans[ampIter]
is None:
601 self.log.warning(
"No overscan information available for ISR statistics for amp %s.",
603 nCols = amp.getRawSerialOverscanBBox().getWidth()
604 ampStats[
"OVERSCAN_COLUMNS"] = np.full((nCols, ), np.nan)
605 ampStats[
"OVERSCAN_VALUES"] = np.full((nCols, ), np.nan)
607 overscanImage = serialOverscans[ampIter].overscanImage
611 for column
in range(0, overscanImage.getWidth()):
618 nRows = amp.getRawSerialOverscanBBox().getHeight()
619 osMean = afwMath.makeStatistics(overscanImage.image.array[:nRows, column],
621 columns.append(column)
624 values.append(osMean * gains[amp.getName()])
628 if readoutCorner
in (ReadoutCorner.LR, ReadoutCorner.UR):
629 ampStats[
"OVERSCAN_COLUMNS"] = list(reversed(columns))
630 ampStats[
"OVERSCAN_VALUES"] = list(reversed(values))
632 ampStats[
"OVERSCAN_COLUMNS"] = columns
633 ampStats[
"OVERSCAN_VALUES"] = values
635 outputStats[amp.getName()] = ampStats
660 """Task to measure banding statistics.
664 inputExp : `lsst.afw.image.Exposure`
666 overscans : `list` [`lsst.pipe.base.Struct`]
667 List of overscan results. Expected fields are:
670 Value or fit subtracted from the amplifier image data
671 (scalar or `lsst.afw.image.Image`).
673 Value or fit subtracted from the overscan image data
674 (scalar or `lsst.afw.image.Image`).
676 Image of the overscan region with the overscan
677 correction applied (`lsst.afw.image.Image`). This
678 quantity is used to estimate the amplifier read noise
683 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
684 Dictionary of measurements, keyed by amplifier name and
689 detector = inputExp.getDetector()
690 kernel = self.
makeKernel(self.config.bandingKernelSize)
692 outputStats[
"AMP_BANDING"] = []
693 for amp, overscanData
in zip(detector.getAmplifiers(), overscans):
694 overscanFit = np.array(overscanData.overscanFit)
695 overscanArray = overscanData.overscanImage.image.array
696 rawOverscan = np.mean(overscanArray + overscanFit, axis=1)
698 smoothedOverscan = np.convolve(rawOverscan, kernel, mode=
"valid")
700 low, high = np.quantile(smoothedOverscan, [self.config.bandingFractionLow,
701 self.config.bandingFractionHigh])
702 outputStats[
"AMP_BANDING"].append(float(high - low))
704 if self.config.bandingUseHalfDetector:
705 fullLength = len(outputStats[
"AMP_BANDING"])
706 outputStats[
"DET_BANDING"] = float(np.nanmedian(outputStats[
"AMP_BANDING"][0:fullLength//2]))
708 outputStats[
"DET_BANDING"] = float(np.nanmedian(outputStats[
"AMP_BANDING"]))
713 """Task to measure metrics from image slicing.
717 inputExp : `lsst.afw.image.Exposure`
719 overscans : `list` [`lsst.pipe.base.Struct`]
720 List of overscan results. Expected fields are:
723 Value or fit subtracted from the amplifier image data
724 (scalar or `lsst.afw.image.Image`).
726 Value or fit subtracted from the overscan image data
727 (scalar or `lsst.afw.image.Image`).
729 Image of the overscan region with the overscan
730 correction applied (`lsst.afw.image.Image`). This
731 quantity is used to estimate the amplifier read noise
736 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
737 Dictionary of measurements, keyed by amplifier name and
742 detector = inputExp.getDetector()
743 kernel = self.
makeKernel(self.config.projectionKernelSize)
745 outputStats[
"AMP_VPROJECTION"] = {}
746 outputStats[
"AMP_HPROJECTION"] = {}
747 convolveMode =
"valid"
748 if self.config.doProjectionFft:
749 outputStats[
"AMP_VFFT_REAL"] = {}
750 outputStats[
"AMP_VFFT_IMAG"] = {}
751 outputStats[
"AMP_HFFT_REAL"] = {}
752 outputStats[
"AMP_HFFT_IMAG"] = {}
753 convolveMode =
"same"
755 for amp
in detector.getAmplifiers():
756 ampArray = inputExp.image[amp.getBBox()].array
758 horizontalProjection = np.mean(ampArray, axis=0)
759 verticalProjection = np.mean(ampArray, axis=1)
761 horizontalProjection = np.convolve(horizontalProjection, kernel, mode=convolveMode)
762 verticalProjection = np.convolve(verticalProjection, kernel, mode=convolveMode)
764 outputStats[
"AMP_HPROJECTION"][amp.getName()] = horizontalProjection.tolist()
765 outputStats[
"AMP_VPROJECTION"][amp.getName()] = verticalProjection.tolist()
767 if self.config.doProjectionFft:
768 horizontalWindow = np.ones_like(horizontalProjection)
769 verticalWindow = np.ones_like(verticalProjection)
770 if self.config.projectionFftWindow ==
"NONE":
772 elif self.config.projectionFftWindow ==
"HAMMING":
773 horizontalWindow = hamming(len(horizontalProjection))
774 verticalWindow = hamming(len(verticalProjection))
775 elif self.config.projectionFftWindow ==
"HANN":
776 horizontalWindow = hann(len(horizontalProjection))
777 verticalWindow = hann(len(verticalProjection))
778 elif self.config.projectionFftWindow ==
"GAUSSIAN":
779 horizontalWindow = gaussian(len(horizontalProjection))
780 verticalWindow = gaussian(len(verticalProjection))
782 raise RuntimeError(f
"Invalid window function: {self.config.projectionFftWindow}")
784 horizontalFFT = np.fft.rfft(np.multiply(horizontalProjection, horizontalWindow))
785 verticalFFT = np.fft.rfft(np.multiply(verticalProjection, verticalWindow))
787 outputStats[
"AMP_HFFT_REAL"][amp.getName()] = np.real(horizontalFFT).tolist()
788 outputStats[
"AMP_HFFT_IMAG"][amp.getName()] = np.imag(horizontalFFT).tolist()
789 outputStats[
"AMP_VFFT_REAL"][amp.getName()] = np.real(verticalFFT).tolist()
790 outputStats[
"AMP_VFFT_IMAG"][amp.getName()] = np.imag(verticalFFT).tolist()
795 """Copy calibration statistics for this exposure.
799 inputExp : `lsst.afw.image.Exposure`
800 The exposure being processed.
802 Keyword arguments with calibrations.
806 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
807 Dictionary of measurements, keyed by amplifier name and
813 for amp
in inputExp.getDetector():
816 for calibType
in (
"bias",
"dark",
"flat"):
817 if kwargs.get(calibType,
None)
is not None:
818 metadata = kwargs[calibType].getMetadata()
819 for pct
in self.config.expectedDistributionLevels:
820 key = f
"LSST CALIB {calibType.upper()} {amp.getName()} DISTRIBUTION {pct}-PCT"
821 ampStats[key] = metadata.get(key, np.nan)
823 for calibType
in (
"defects"):
824 if kwargs.get(calibType,
None)
is not None:
825 metadata = kwargs[calibType].getMetadata()
826 for key
in (f
"LSST CALIB {calibType.upper()} {amp.getName()} N_HOT",
827 f
"LSST CALIB {calibType.upper()} {amp.getName()} N_COLD"):
828 ampStats[key] = metadata.get(key, np.nan)
829 outputStats[amp.getName()] = ampStats
832 for calibType
in (
"defects"):
833 if kwargs.get(calibType,
None)
is not None:
834 metadata = kwargs[calibType].getMetadata()
835 for key
in (f
"LSST CALIB {calibType.upper()} N_BAD_COLUMNS"):
836 outputStats[
"detector"][key] = metadata.get(key, np.nan)
841 """Measure number of bias shifts from overscan data.
845 inputExp : `lsst.afw.image.Exposure`
847 overscans : `list` [`lsst.pipe.base.Struct`]
848 List of overscan results. Expected fields are:
851 Value or fit subtracted from the amplifier image data
852 (scalar or `lsst.afw.image.Image`).
854 Value or fit subtracted from the overscan image data
855 (scalar or `lsst.afw.image.Image`).
857 Image of the overscan region with the overscan
858 correction applied (`lsst.afw.image.Image`). This
859 quantity is used to estimate the amplifier read noise
864 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
865 Dictionary of measurements, keyed by amplifier name and
870 Based on eo_pipe implementation:
871 https://github.com/lsst-camera-dh/eo_pipe/blob/main/python/lsst/eo/pipe/biasShiftsTask.py # noqa: E501 W505
875 detector = inputExp.getDetector()
876 for amp, overscans
in zip(detector, serialOverscanResults):
882 if overscans
is None:
883 ampStats[
"LOCAL_NOISE"] = 0.0
884 ampStats[
"BIAS_SHIFTS"] = []
885 outputStats[amp.getName()] = ampStats
889 rawOverscan = overscans.overscanImage.image.array + overscans.overscanFit
892 rawOverscan = np.mean(rawOverscan[:, self.config.biasShiftColumnSkip:], axis=1)
896 ampStats[
"LOCAL_NOISE"] = float(noise)
897 ampStats[
"BIAS_SHIFTS"] = shift_peaks
899 outputStats[amp.getName()] = ampStats
903 """Scan overscan data for shifts.
907 overscanData : `list` [`float`]
908 Overscan data to search for shifts.
913 Noise estimated from Butterworth filtered overscan data.
914 peaks : `list` [`float`, `float`, `int`, `int`]
915 Shift peak information, containing the convolved peak
916 value, the raw peak value, and the lower and upper bounds
917 of the region checked.
919 numerator, denominator = butter(self.config.biasShiftFilterOrder,
920 self.config.biasShiftCutoff,
921 btype=
"high", analog=
False)
922 noise = np.std(filtfilt(numerator, denominator, overscanData))
923 kernel = np.concatenate([np.arange(self.config.biasShiftWindow),
924 np.arange(-self.config.biasShiftWindow + 1, 0)])
925 kernel = kernel/np.sum(kernel[:self.config.biasShiftWindow])
927 convolved = np.convolve(overscanData, kernel, mode=
"valid")
928 convolved = np.pad(convolved, (self.config.biasShiftWindow - 1, self.config.biasShiftWindow))
930 shift_check = np.abs(convolved)/noise
931 shift_mask = shift_check > self.config.biasShiftThreshold
932 shift_mask[:self.config.biasShiftRowSkip] =
False
934 shift_regions = np.flatnonzero(np.diff(np.r_[np.int8(0),
935 shift_mask.view(np.int8),
936 np.int8(0)])).reshape(-1, 2)
938 for region
in shift_regions:
939 region_peak = np.argmax(shift_check[region[0]:region[1]]) + region[0]
942 [float(convolved[region_peak]), float(region_peak),
943 int(region[0]), int(region[1])])
944 return noise, shift_peaks
947 """Determine if a region is flat.
952 Row with possible peak.
954 Value at the possible peak.
955 overscanData : `list` [`float`]
956 Overscan data used to fit around the possible peak.
961 Indicates if the region is flat, and so the peak is valid.
963 prerange = np.arange(shiftRow - self.config.biasShiftWindow, shiftRow)
964 postrange = np.arange(shiftRow, shiftRow + self.config.biasShiftWindow)
966 preFit = linregress(prerange, overscanData[prerange])
967 postFit = linregress(postrange, overscanData[postrange])
970 preTrend = (2*preFit[0]*len(prerange) < shiftPeak)
971 postTrend = (2*postFit[0]*len(postrange) < shiftPeak)
973 preTrend = (2*preFit[0]*len(prerange) > shiftPeak)
974 postTrend = (2*postFit[0]*len(postrange) > shiftPeak)
976 return (preTrend
and postTrend)
979 """Measure correlations between amplifier segments.
983 inputExp : `lsst.afw.image.Exposure`
985 overscans : `list` [`lsst.pipe.base.Struct`]
986 List of overscan results. Expected fields are:
989 Value or fit subtracted from the amplifier image data
990 (scalar or `lsst.afw.image.Image`).
992 Value or fit subtracted from the overscan image data
993 (scalar or `lsst.afw.image.Image`).
995 Image of the overscan region with the overscan
996 correction applied (`lsst.afw.image.Image`). This
997 quantity is used to estimate the amplifier read noise
1002 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
1003 Dictionary of measurements, keyed by amplifier name and
1008 Based on eo_pipe implementation:
1009 https://github.com/lsst-camera-dh/eo_pipe/blob/main/python/lsst/eo/pipe/raft_level_correlations.py # noqa: E501 W505
1011 detector = inputExp.getDetector()
1014 for ampId, overscan
in enumerate(serialOverscanResults):
1018 if overscan
is None:
1019 for ampId2, overscan2
in enumerate(serialOverscanResults):
1024 'detector': detector.getId(),
1025 'detectorComp': detector.getId(),
1026 'ampName': detector[ampId].getName(),
1027 'ampComp': detector[ampId2].getName(),
1028 'imageCorr': float(imC),
1029 'overscanCorr': float(osC),
1033 rawOverscan = overscan.overscanImage.image.array + overscan.overscanFit
1034 rawOverscan = rawOverscan.ravel()
1036 ampImage = inputExp[detector[ampId].getBBox()]
1037 ampImage = ampImage.image.array.ravel()
1039 for ampId2, overscan2
in enumerate(serialOverscanResults):
1046 if overscan2
is None:
1050 rawOverscan2 = overscan2.overscanImage.image.array + overscan2.overscanFit
1051 rawOverscan2 = rawOverscan2.ravel()
1053 osC = np.corrcoef(rawOverscan, rawOverscan2)[0, 1]
1055 ampImage2 = inputExp[detector[ampId2].getBBox()]
1056 ampImage2 = ampImage2.image.array.ravel()
1058 imC = np.corrcoef(ampImage, ampImage2)[0, 1]
1060 'detector': detector.getId(),
1061 'detectorComp': detector.getId(),
1062 'ampName': detector[ampId].getName(),
1063 'ampComp': detector[ampId2].getName(),
1064 'imageCorr': float(imC),
1065 'overscanCorr': float(osC),
1071 """Measure Max Divisadero Tearing effect per amp.
1075 inputExp : `lsst.afw.image.Exposure`
1076 Exposure to measure. Usually a flat.
1078 The flat will be selected from here.
1082 outputStats : `dict` [`str`, [`dict` [`str`, `float`]]]
1083 Dictionary of measurements, keyed by amplifier name and
1085 Measurements include
1087 - DIVISADERO_PROFILE: Robust mean of rows between
1088 divisaderoProjection<Maximum|Minumum> on readout edge of ccd
1089 normalized by a linear fit to the same rows.
1090 - DIVISADERO_MAX_PAIR: Tuple of maximum of the absolute values of
1091 the DIVISADERO_PROFILE, for number of pixels (specified by
1092 divisaderoNumImpactPixels on left and right side of amp.
1093 - DIVISADERO_MAX: Maximum of the absolute values of the
1094 the DIVISADERO_PROFILE, for the divisaderoNumImpactPixels on
1095 boundaries of neighboring amps (including the pixels in those
1096 neighborboring amps).
1100 for amp
in inputExp.getDetector():
1102 ampArray = inputExp.image[amp.getBBox()].array
1104 if amp.getReadoutCorner().name
in (
'UL',
'UR'):
1105 minRow = amp.getBBox().getHeight() - self.config.divisaderoProjectionMaximum
1106 maxRow = amp.getBBox().getHeight() - self.config.divisaderoProjectionMinimum
1108 minRow = self.config.divisaderoProjectionMinimum
1109 maxRow = self.config.divisaderoProjectionMaximum
1111 segment = slice(minRow, maxRow)
1112 projection, _, _ = astropy.stats.sigma_clipped_stats(ampArray[segment, :], axis=0)
1115 projection /= np.median(projection)
1116 columns = np.arange(len(projection))
1118 segment = slice(self.config.divisaderoEdgePixels, -self.config.divisaderoEdgePixels)
1119 model = np.polyfit(columns[segment], projection[segment], 1)
1120 modelProjection = model[0] * columns + model[1]
1121 divisaderoProfile = projection / modelProjection
1124 leftMax = np.nanmax(np.abs(divisaderoProfile[0:self.config.divisaderoNumImpactPixels] - 1.0))
1125 rightMax = np.nanmax(np.abs(divisaderoProfile[-self.config.divisaderoNumImpactPixels:] - 1.0))
1127 ampStats[
'DIVISADERO_PROFILE'] = np.array(divisaderoProfile).tolist()
1128 ampStats[
'DIVISADERO_MAX_PAIR'] = [leftMax, rightMax]
1129 outputStats[amp.getName()] = ampStats
1131 detector = inputExp.getDetector()
1132 xCenters = [amp.getBBox().getCenterX()
for amp
in detector]
1133 yCenters = [amp.getBBox().getCenterY()
for amp
in detector]
1134 xIndices = np.ceil(xCenters / np.min(xCenters) / 2).astype(int) - 1
1135 yIndices = np.ceil(yCenters / np.min(yCenters) / 2).astype(int) - 1
1136 ampIds = np.zeros((len(set(yIndices)), len(set(xIndices))), dtype=int)
1137 for ampId, xIndex, yIndex
in zip(np.arange(len(detector)), xIndices, yIndices):
1138 ampIds[yIndex, xIndex] = ampId
1142 for i, amp
in enumerate(detector):
1143 y, x = np.where(ampIds == i)
1144 end = ampIds.shape[1] - 1
1147 thisAmpsPair = outputStats[amp.getName()][
'DIVISADERO_MAX_PAIR']
1151 myMax = thisAmpsPair[1]
1153 neighborMax = outputStats[detector[ampIds[yInd, 1]].getName()][
'DIVISADERO_MAX_PAIR'][0]
1156 myMax = thisAmpsPair[0]
1158 neighborMax = outputStats[detector[ampIds[yInd, end - 1]].getName()][
'DIVISADERO_MAX_PAIR'][1]
1161 myMax = max(thisAmpsPair)
1162 leftName = detector[ampIds[yInd, max(xInd - 1, 0)]].getName()
1163 rightName = detector[ampIds[yInd, min(xInd + 1, ampIds.shape[1] - 1)]].getName()
1166 neighborMax = max(outputStats[leftName][
'DIVISADERO_MAX_PAIR'][1],
1167 outputStats[rightName][
'DIVISADERO_MAX_PAIR'][0])
1169 divisaderoMax = max([myMax, neighborMax])
1170 outputStats[amp.getName()][
'DIVISADERO_MAX'] = divisaderoMax