23 import matplotlib.pyplot
as plt
24 from collections
import Counter
29 from .utils
import (fitLeastSq, fitBootstrap, funcPolynomial, funcAstier)
30 from scipy.optimize
import least_squares
34 from .astierCovPtcUtils
import (fftSize, CovFft, computeCovDirect, fitData)
35 from .linearity
import LinearitySolveTask
36 from .photodiode
import getBOTphotodiodeData
38 from lsst.pipe.tasks.getRepositoryData
import DataRefListRunner
40 __all__ = [
'MeasurePhotonTransferCurveTask',
41 'MeasurePhotonTransferCurveTaskConfig',
42 'PhotonTransferCurveDataset']
46 """Config class for photon transfer curve measurement task"""
47 ccdKey = pexConfig.Field(
49 doc=
"The key by which to pull a detector from a dataId, e.g. 'ccd' or 'detector'.",
52 ptcFitType = pexConfig.ChoiceField(
54 doc=
"Fit PTC to approximation in Astier+19 (Equation 16) or to a polynomial.",
57 "POLYNOMIAL":
"n-degree polynomial (use 'polynomialFitDegree' to set 'n').",
58 "EXPAPPROXIMATION":
"Approximation in Astier+19 (Eq. 16).",
59 "FULLCOVARIANCE":
"Full covariances model in Astier+19 (Eq. 20)"
62 sigmaClipFullFitCovariancesAstier = pexConfig.Field(
64 doc=
"sigma clip for full model fit for FULLCOVARIANCE ptcFitType ",
67 maxIterFullFitCovariancesAstier = pexConfig.Field(
69 doc=
"Maximum number of iterations in full model fit for FULLCOVARIANCE ptcFitType",
72 maximumRangeCovariancesAstier = pexConfig.Field(
74 doc=
"Maximum range of covariances as in Astier+19",
77 covAstierRealSpace = pexConfig.Field(
79 doc=
"Calculate covariances in real space or via FFT? (see appendix A of Astier+19).",
82 polynomialFitDegree = pexConfig.Field(
84 doc=
"Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.",
87 linearity = pexConfig.ConfigurableField(
88 target=LinearitySolveTask,
89 doc=
"Task to solve the linearity."
92 doCreateLinearizer = pexConfig.Field(
94 doc=
"Calculate non-linearity and persist linearizer?",
98 binSize = pexConfig.Field(
100 doc=
"Bin the image by this factor in both dimensions.",
103 minMeanSignal = pexConfig.Field(
105 doc=
"Minimum value (inclusive) of mean signal (in DN) above which to consider.",
108 maxMeanSignal = pexConfig.Field(
110 doc=
"Maximum value (inclusive) of mean signal (in DN) below which to consider.",
113 initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField(
115 doc=
"Initially exclude data points with a variance that are more than a factor of this from being"
116 " linear in the positive direction, from the PTC fit. Note that these points will also be"
117 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
118 " to allow an accurate determination of the sigmas for said iterative fit.",
123 initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField(
125 doc=
"Initially exclude data points with a variance that are more than a factor of this from being"
126 " linear in the negative direction, from the PTC fit. Note that these points will also be"
127 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
128 " to allow an accurate determination of the sigmas for said iterative fit.",
133 sigmaCutPtcOutliers = pexConfig.Field(
135 doc=
"Sigma cut for outlier rejection in PTC.",
138 maskNameList = pexConfig.ListField(
140 doc=
"Mask list to exclude from statistics calculations.",
141 default=[
'SUSPECT',
'BAD',
'NO_DATA'],
143 nSigmaClipPtc = pexConfig.Field(
145 doc=
"Sigma cut for afwMath.StatisticsControl()",
148 nIterSigmaClipPtc = pexConfig.Field(
150 doc=
"Number of sigma-clipping iterations for afwMath.StatisticsControl()",
153 maxIterationsPtcOutliers = pexConfig.Field(
155 doc=
"Maximum number of iterations for outlier rejection in PTC.",
158 doFitBootstrap = pexConfig.Field(
160 doc=
"Use bootstrap for the PTC fit parameters and errors?.",
163 doPhotodiode = pexConfig.Field(
165 doc=
"Apply a correction based on the photodiode readings if available?",
168 photodiodeDataPath = pexConfig.Field(
170 doc=
"Gen2 only: path to locate the data photodiode data files.",
173 instrumentName = pexConfig.Field(
175 doc=
"Instrument name.",
181 """A simple class to hold the output data from the PTC task.
183 The dataset is made up of a dictionary for each item, keyed by the
184 amplifiers' names, which much be supplied at construction time.
186 New items cannot be added to the class to save accidentally saving to the
187 wrong property, and the class can be frozen if desired.
189 inputExpIdPairs records the exposures used to produce the data.
190 When fitPtc() or fitCovariancesAstier() is run, a mask is built up, which is by definition
191 always the same length as inputExpIdPairs, rawExpTimes, rawMeans
192 and rawVars, and is a list of bools, which are incrementally set to False
193 as points are discarded from the fits.
195 PTC fit parameters for polynomials are stored in a list in ascending order
196 of polynomial term, i.e. par[0]*x^0 + par[1]*x + par[2]*x^2 etc
197 with the length of the list corresponding to the order of the polynomial
203 List with the names of the amplifiers of the detector at hand.
206 Type of model fitted to the PTC: "POLYNOMIAL", "EXPAPPROXIMATION", or "FULLCOVARIANCE".
210 `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
211 Output dataset from MeasurePhotonTransferCurveTask.
218 self.__dict__[
"ptcFitType"] = ptcFitType
219 self.__dict__[
"ampNames"] = ampNames
220 self.__dict__[
"badAmps"] = []
225 self.__dict__[
"inputExpIdPairs"] = {ampName: []
for ampName
in ampNames}
226 self.__dict__[
"expIdMask"] = {ampName: []
for ampName
in ampNames}
227 self.__dict__[
"rawExpTimes"] = {ampName: []
for ampName
in ampNames}
228 self.__dict__[
"rawMeans"] = {ampName: []
for ampName
in ampNames}
229 self.__dict__[
"rawVars"] = {ampName: []
for ampName
in ampNames}
230 self.__dict__[
"photoCharge"] = {ampName: []
for ampName
in ampNames}
233 self.__dict__[
"gain"] = {ampName: -1.
for ampName
in ampNames}
234 self.__dict__[
"gainErr"] = {ampName: -1.
for ampName
in ampNames}
235 self.__dict__[
"noise"] = {ampName: -1.
for ampName
in ampNames}
236 self.__dict__[
"noiseErr"] = {ampName: -1.
for ampName
in ampNames}
240 self.__dict__[
"ptcFitPars"] = {ampName: []
for ampName
in ampNames}
241 self.__dict__[
"ptcFitParsError"] = {ampName: []
for ampName
in ampNames}
242 self.__dict__[
"ptcFitReducedChiSquared"] = {ampName: []
for ampName
in ampNames}
249 self.__dict__[
"covariancesTuple"] = {ampName: []
for ampName
in ampNames}
250 self.__dict__[
"covariancesFitsWithNoB"] = {ampName: []
for ampName
in ampNames}
251 self.__dict__[
"covariancesFits"] = {ampName: []
for ampName
in ampNames}
252 self.__dict__[
"aMatrix"] = {ampName: []
for ampName
in ampNames}
253 self.__dict__[
"bMatrix"] = {ampName: []
for ampName
in ampNames}
256 self.__dict__[
"finalVars"] = {ampName: []
for ampName
in ampNames}
257 self.__dict__[
"finalModelVars"] = {ampName: []
for ampName
in ampNames}
258 self.__dict__[
"finalMeans"] = {ampName: []
for ampName
in ampNames}
261 """Protect class attributes"""
262 if attribute
not in self.__dict__:
263 raise AttributeError(f
"{attribute} is not already a member of PhotonTransferCurveDataset, which"
264 " does not support setting of new attributes.")
266 self.__dict__[attribute] = value
269 """Get the exposures used, i.e. not discarded, for a given amp.
271 If no mask has been created yet, all exposures are returned.
273 if len(self.expIdMask[ampName]) == 0:
274 return self.inputExpIdPairs[ampName]
277 assert len(self.expIdMask[ampName]) == len(self.inputExpIdPairs[ampName])
279 pairs = self.inputExpIdPairs[ampName]
280 mask = self.expIdMask[ampName]
282 return [(exp1, exp2)
for ((exp1, exp2), m)
in zip(pairs, mask)
if bool(m)
is True]
285 return [amp
for amp
in self.ampNames
if amp
not in self.badAmps]
289 """A class to calculate, fit, and plot a PTC from a set of flat pairs.
291 The Photon Transfer Curve (var(signal) vs mean(signal)) is a standard tool
292 used in astronomical detectors characterization (e.g., Janesick 2001,
293 Janesick 2007). If ptcFitType is "EXPAPPROXIMATION" or "POLYNOMIAL", this task calculates the
294 PTC from a series of pairs of flat-field images; each pair taken at identical exposure
295 times. The difference image of each pair is formed to eliminate fixed pattern noise,
296 and then the variance of the difference image and the mean of the average image
297 are used to produce the PTC. An n-degree polynomial or the approximation in Equation
298 16 of Astier+19 ("The Shape of the Photon Transfer Curve of CCD sensors",
299 arXiv:1905.08677) can be fitted to the PTC curve. These models include
300 parameters such as the gain (e/DN) and readout noise.
302 Linearizers to correct for signal-chain non-linearity are also calculated.
303 The `Linearizer` class, in general, can support per-amp linearizers, but in this
304 task this is not supported.
306 If ptcFitType is "FULLCOVARIANCE", the covariances of the difference images are calculated via the
307 DFT methods described in Astier+19 and the variances for the PTC are given by the cov[0,0] elements
308 at each signal level. The full model in Equation 20 of Astier+19 is fit to the PTC to get the gain
315 Positional arguments passed to the Task constructor. None used at this
318 Keyword arguments passed on to the Task constructor. None used at this
323 RunnerClass = DataRefListRunner
324 ConfigClass = MeasurePhotonTransferCurveTaskConfig
325 _DefaultName =
"measurePhotonTransferCurve"
328 pipeBase.CmdLineTask.__init__(self, *args, **kwargs)
329 self.makeSubtask(
"linearity")
330 plt.interactive(
False)
331 self.config.validate()
336 """Run the Photon Transfer Curve (PTC) measurement task.
338 For a dataRef (which is each detector here),
339 and given a list of exposure pairs (postISR) at different exposure times,
344 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
345 Data references for exposures for detectors to process.
347 if len(dataRefList) < 2:
348 raise RuntimeError(
"Insufficient inputs to combine.")
351 dataRef = dataRefList[0]
353 detNum = dataRef.dataId[self.config.ccdKey]
354 camera = dataRef.get(
'camera')
355 detector = camera[dataRef.dataId[self.config.ccdKey]]
357 amps = detector.getAmplifiers()
358 ampNames = [amp.getName()
for amp
in amps]
364 for (exp1, exp2)
in expPairs.values():
365 id1 = exp1.getInfo().getVisitInfo().getExposureId()
366 id2 = exp2.getInfo().getVisitInfo().getExposureId()
367 expIds.append((id1, id2))
368 self.log.info(f
"Measuring PTC using {expIds} exposures for detector {detector.getId()}")
373 if self.config.doPhotodiode:
374 for (expId1, expId2)
in expIds:
376 for i, expId
in enumerate([expId1, expId2]):
381 dataRef.dataId[
'expId'] = expId//1000
382 if self.config.photodiodeDataPath:
387 charges[i] = photodiodeData.getCharge()
391 self.log.warn(f
"No photodiode data found for {expId}")
393 for ampName
in ampNames:
394 datasetPtc.photoCharge[ampName].append((charges[0], charges[1]))
398 for expTime, (exp1, exp2)
in expPairs.items():
399 expId1 = exp1.getInfo().getVisitInfo().getExposureId()
400 expId2 = exp2.getInfo().getVisitInfo().getExposureId()
403 for ampNumber, amp
in enumerate(detector):
404 ampName = amp.getName()
406 doRealSpace = self.config.covAstierRealSpace
407 muDiff, varDiff, covAstier = self.
measureMeanVarCov(exp1, exp2, region=amp.getBBox(),
408 covAstierRealSpace=doRealSpace)
409 if np.isnan(muDiff)
or np.isnan(varDiff)
or (covAstier
is None):
410 msg = (f
"NaN mean or var, or None cov in amp {ampNumber} in exposure pair {expId1},"
411 f
" {expId2} of detector {detNum}.")
415 tags = [
'mu',
'i',
'j',
'var',
'cov',
'npix',
'ext',
'expTime',
'ampName']
416 if (muDiff <= self.config.minMeanSignal)
or (muDiff >= self.config.maxMeanSignal):
418 datasetPtc.rawExpTimes[ampName].append(expTime)
419 datasetPtc.rawMeans[ampName].append(muDiff)
420 datasetPtc.rawVars[ampName].append(varDiff)
421 datasetPtc.inputExpIdPairs[ampName].append((expId1, expId2))
423 tupleRows += [(muDiff, ) + covRow + (ampNumber, expTime, ampName)
for covRow
in covAstier]
424 if nAmpsNan == len(ampNames):
425 msg = f
"NaN mean in all amps of exposure pair {expId1}, {expId2} of detector {detNum}."
429 tupleRecords += tupleRows
430 covariancesWithTags = np.core.records.fromrecords(tupleRecords, names=allTags)
432 if self.config.ptcFitType
in [
"FULLCOVARIANCE", ]:
435 elif self.config.ptcFitType
in [
"EXPAPPROXIMATION",
"POLYNOMIAL"]:
438 datasetPtc = self.
fitPtc(datasetPtc, self.config.ptcFitType)
441 if self.config.doCreateLinearizer:
447 dimensions = {
'camera': camera.getName(),
'detector': detector.getId()}
448 linearityResults = self.linearity.run(datasetPtc, camera, dimensions)
449 linearizer = linearityResults.outputLinearizer
451 butler = dataRef.getButler()
452 self.log.info(
"Writing linearizer:")
454 detName = detector.getName()
455 now = datetime.datetime.utcnow()
456 calibDate = now.strftime(
"%Y-%m-%d")
458 butler.put(linearizer, datasetType=
'Linearizer', dataId={
'detector': detNum,
459 'detectorName': detName,
'calibDate': calibDate})
460 self.log.info(
"Writing PTC data.")
461 dataRef.put(datasetPtc, datasetType=
"photonTransferCurveDataset")
463 return pipeBase.Struct(exitStatus=0)
466 """Produce a list of flat pairs indexed by exposure time.
470 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
471 Data references for exposures for detectors to process.
475 flatPairs : `dict` [`float`, `lsst.afw.image.exposure.exposure.ExposureF`]
476 Dictionary that groups flat-field exposures that have the same exposure time (seconds).
480 We use the difference of one pair of flat-field images taken at the same exposure time when
481 calculating the PTC to reduce Fixed Pattern Noise. If there are > 2 flat-field images with the
482 same exposure time, the first two are kept and the rest discarded.
487 for dataRef
in dataRefList:
489 tempFlat = dataRef.get(
"postISRCCD")
491 self.log.warn(
"postISR exposure could not be retrieved. Ignoring flat.")
493 expDate = tempFlat.getInfo().getVisitInfo().getDate().get()
494 expDict.setdefault(expDate, tempFlat)
495 sortedExps = {k: expDict[k]
for k
in sorted(expDict)}
498 for exp
in sortedExps:
499 tempFlat = sortedExps[exp]
500 expTime = tempFlat.getInfo().getVisitInfo().getExposureTime()
501 listAtExpTime = flatPairs.setdefault(expTime, [])
502 if len(listAtExpTime) < 2:
503 listAtExpTime.append(tempFlat)
504 if len(listAtExpTime) > 2:
505 self.log.warn(f
"More than 2 exposures found at expTime {expTime}. Dropping exposures "
506 f
"{listAtExpTime[2:]}.")
508 for (key, value)
in flatPairs.items():
511 self.log.warn(f
"Only one exposure found at expTime {key}. Dropping exposure {value}.")
515 """Fit measured flat covariances to full model in Astier+19.
519 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
520 The dataset containing information such as the means, variances and exposure times.
522 covariancesWithTagsArray : `numpy.recarray`
523 Tuple with at least (mu, cov, var, i, j, npix), where:
524 mu : 0.5*(m1 + m2), where:
525 mu1: mean value of flat1
526 mu2: mean value of flat2
527 cov: covariance value at lag(i, j)
528 var: variance(covariance value at lag(0, 0))
531 npix: number of pixels used for covariance calculation.
535 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
536 This is the same dataset as the input paramter, however, it has been modified
537 to include information such as the fit vectors and the fit parameters. See
538 the class `PhotonTransferCurveDatase`.
541 covFits, covFitsNoB =
fitData(covariancesWithTagsArray, maxMu=self.config.maxMeanSignal,
542 r=self.config.maximumRangeCovariancesAstier,
543 nSigmaFullFit=self.config.sigmaClipFullFitCovariancesAstier,
544 maxIterFullFit=self.config.maxIterFullFitCovariancesAstier)
546 dataset.covariancesTuple = covariancesWithTagsArray
547 dataset.covariancesFits = covFits
548 dataset.covariancesFitsWithNoB = covFitsNoB
554 """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects.
558 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
559 The dataset containing information such as the means, variances and exposure times.
562 Dictionary of CovFit objects, with amp names as keys.
566 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
567 This is the same dataset as the input paramter, however, it has been modified
568 to include extra information such as the mask 1D array, gains, reoudout noise, measured signal,
569 measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier.
570 See the class `PhotonTransferCurveDatase`.
573 for i, amp
in enumerate(covFits):
575 (meanVecFinal, varVecFinal, varVecModel,
576 wc, varMask) = fit.getFitData(0, 0, divideByMu=
False, returnMasked=
True)
578 dataset.expIdMask[amp] = varMask
579 dataset.gain[amp] = gain
580 dataset.gainErr[amp] = fit.getGainErr()
581 dataset.noise[amp] = np.sqrt(np.fabs(fit.getRon()))
582 dataset.noiseErr[amp] = fit.getRonErr()
583 dataset.finalVars[amp].append(varVecFinal/(gain**2))
584 dataset.finalModelVars[amp].append(varVecModel/(gain**2))
585 dataset.finalMeans[amp].append(meanVecFinal/gain)
586 dataset.aMatrix[amp].append(fit.getA())
587 dataset.bMatrix[amp].append(fit.getB())
592 """Calculate the mean of each of two exposures and the variance and covariance of their difference.
594 The variance is calculated via afwMath, and the covariance via the methods in Astier+19 (appendix A).
595 In theory, var = covariance[0,0]. This should be validated, and in the future, we may decide to just
596 keep one (covariance).
600 exposure1 : `lsst.afw.image.exposure.exposure.ExposureF`
601 First exposure of flat field pair.
603 exposure2 : `lsst.afw.image.exposure.exposure.ExposureF`
604 Second exposure of flat field pair.
606 region : `lsst.geom.Box2I`, optional
607 Region of each exposure where to perform the calculations (e.g, an amplifier).
609 covAstierRealSpace : `bool`, optional
610 Should the covariannces in Astier+19 be calculated in real space or via FFT?
611 See Appendix A of Astier+19.
615 mu : `float` or `NaN`
616 0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in
617 both exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
619 varDiff : `float` or `NaN`
620 Half of the clipped variance of the difference of the regions inthe two input
621 exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
623 covDiffAstier : `list` or `NaN`
624 List with tuples of the form (dx, dy, var, cov, npix), where:
630 Variance at (dx, dy).
632 Covariance at (dx, dy).
634 Number of pixel pairs used to evaluate var and cov.
635 If either mu1 or m2 are NaN's, the returned value is NaN.
638 if region
is not None:
639 im1Area = exposure1.maskedImage[region]
640 im2Area = exposure2.maskedImage[region]
642 im1Area = exposure1.maskedImage
643 im2Area = exposure2.maskedImage
645 if self.config.binSize > 1:
646 im1Area = afwMath.binImage(im1Area, self.config.binSize)
647 im2Area = afwMath.binImage(im2Area, self.config.binSize)
649 im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList)
650 im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
651 self.config.nIterSigmaClipPtc,
653 im1StatsCtrl.setNanSafe(
True)
654 im1StatsCtrl.setAndMask(im1MaskVal)
656 im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList)
657 im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
658 self.config.nIterSigmaClipPtc,
660 im2StatsCtrl.setNanSafe(
True)
661 im2StatsCtrl.setAndMask(im2MaskVal)
664 mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue()
665 mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue()
666 if np.isnan(mu1)
or np.isnan(mu2):
667 return np.nan, np.nan,
None
672 temp = im2Area.clone()
674 diffIm = im1Area.clone()
679 diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList)
680 diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
681 self.config.nIterSigmaClipPtc,
683 diffImStatsCtrl.setNanSafe(
True)
684 diffImStatsCtrl.setAndMask(diffImMaskVal)
686 varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue())
689 w1 = np.where(im1Area.getMask().getArray() == 0, 1, 0)
690 w2 = np.where(im2Area.getMask().getArray() == 0, 1, 0)
693 wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0)
696 maxRangeCov = self.config.maximumRangeCovariancesAstier
697 if covAstierRealSpace:
698 covDiffAstier =
computeCovDirect(diffIm.getImage().getArray(), w, maxRangeCov)
700 shapeDiff = diffIm.getImage().getArray().shape
701 fftShape = (
fftSize(shapeDiff[0] + maxRangeCov),
fftSize(shapeDiff[1]+maxRangeCov))
702 c =
CovFft(diffIm.getImage().getArray(), w, fftShape, maxRangeCov)
703 covDiffAstier = c.reportCovFft(maxRangeCov)
705 return mu, varDiff, covDiffAstier
708 """Compute covariances of diffImage in real space.
710 For lags larger than ~25, it is slower than the FFT way.
711 Taken from https://github.com/PierreAstier/bfptc/
715 diffImage : `numpy.array`
716 Image to compute the covariance of.
718 weightImage : `numpy.array`
719 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
722 Last index of the covariance to be computed.
727 List with tuples of the form (dx, dy, var, cov, npix), where:
733 Variance at (dx, dy).
735 Covariance at (dx, dy).
737 Number of pixel pairs used to evaluate var and cov.
742 for dy
in range(maxRange + 1):
743 for dx
in range(0, maxRange + 1):
746 cov2, nPix2 = self.
covDirectValue(diffImage, weightImage, dx, -dy)
747 cov = 0.5*(cov1 + cov2)
751 if (dx == 0
and dy == 0):
753 outList.append((dx, dy, var, cov, nPix))
758 """Compute covariances of diffImage in real space at lag (dx, dy).
760 Taken from https://github.com/PierreAstier/bfptc/ (c.f., appendix of Astier+19).
764 diffImage : `numpy.array`
765 Image to compute the covariance of.
767 weightImage : `numpy.array`
768 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
779 Covariance at (dx, dy)
782 Number of pixel pairs used to evaluate var and cov.
784 (nCols, nRows) = diffImage.shape
788 (dx, dy) = (-dx, -dy)
792 im1 = diffImage[dy:, dx:]
793 w1 = weightImage[dy:, dx:]
794 im2 = diffImage[:nCols - dy, :nRows - dx]
795 w2 = weightImage[:nCols - dy, :nRows - dx]
797 im1 = diffImage[:nCols + dy, dx:]
798 w1 = weightImage[:nCols + dy, dx:]
799 im2 = diffImage[-dy:, :nRows - dx]
800 w2 = weightImage[-dy:, :nRows - dx]
806 s1 = im1TimesW.sum()/nPix
807 s2 = (im2*wAll).sum()/nPix
808 p = (im1TimesW*im2).sum()/nPix
814 def _initialParsForPolynomial(order):
816 pars = np.zeros(order, dtype=np.float)
823 def _boundsForPolynomial(initialPars):
824 lowers = [np.NINF
for p
in initialPars]
825 uppers = [np.inf
for p
in initialPars]
827 return (lowers, uppers)
830 def _boundsForAstier(initialPars):
831 lowers = [np.NINF
for p
in initialPars]
832 uppers = [np.inf
for p
in initialPars]
833 return (lowers, uppers)
836 def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative):
837 """Return a boolean array to mask bad points.
839 A linear function has a constant ratio, so find the median
840 value of the ratios, and exclude the points that deviate
841 from that by more than a factor of maxDeviationPositive/negative.
842 Asymmetric deviations are supported as we expect the PTC to turn
843 down as the flux increases, but sometimes it anomalously turns
844 upwards just before turning over, which ruins the fits, so it
845 is wise to be stricter about restricting positive outliers than
848 Too high and points that are so bad that fit will fail will be included
849 Too low and the non-linear points will be excluded, biasing the NL fit."""
850 ratios = [b/a
for (a, b)
in zip(means, variances)]
851 medianRatio = np.median(ratios)
852 ratioDeviations = [(r/medianRatio)-1
for r
in ratios]
855 maxDeviationPositive = abs(maxDeviationPositive)
856 maxDeviationNegative = -1. * abs(maxDeviationNegative)
858 goodPoints = np.array([
True if (r < maxDeviationPositive
and r > maxDeviationNegative)
859 else False for r
in ratioDeviations])
862 def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9):
864 nBad = Counter(array)[0]
869 msg = f
"Found {nBad} zeros in array at elements {[x for x in np.where(array==0)[0]]}"
872 array[array == 0] = substituteValue
876 """Fit the photon transfer curve to a polynimial or to Astier+19 approximation.
878 Fit the photon transfer curve with either a polynomial of the order
879 specified in the task config, or using the Astier approximation.
881 Sigma clipping is performed iteratively for the fit, as well as an
882 initial clipping of data points that are more than
883 config.initialNonLinearityExclusionThreshold away from lying on a
884 straight line. This other step is necessary because the photon transfer
885 curve turns over catastrophically at very high flux (because saturation
886 drops the variance to ~0) and these far outliers cause the initial fit
887 to fail, meaning the sigma cannot be calculated to perform the
892 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
893 The dataset containing the means, variances and exposure times
896 Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or
897 'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC
901 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
902 This is the same dataset as the input paramter, however, it has been modified
903 to include information such as the fit vectors and the fit parameters. See
904 the class `PhotonTransferCurveDatase`.
907 def errFunc(p, x, y):
908 return ptcFunc(p, x) - y
910 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers
911 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers
913 for i, ampName
in enumerate(dataset.ampNames):
914 timeVecOriginal = np.array(dataset.rawExpTimes[ampName])
915 meanVecOriginal = np.array(dataset.rawMeans[ampName])
916 varVecOriginal = np.array(dataset.rawVars[ampName])
919 mask = ((meanVecOriginal >= self.config.minMeanSignal) &
920 (meanVecOriginal <= self.config.maxMeanSignal))
923 self.config.initialNonLinearityExclusionThresholdPositive,
924 self.config.initialNonLinearityExclusionThresholdNegative)
925 if not (mask.any()
and goodPoints.any()):
926 msg = (f
"\nSERIOUS: All points in either mask: {mask} or goodPoints: {goodPoints} are bad."
927 f
"Setting {ampName} to BAD.")
931 dataset.badAmps.append(ampName)
932 dataset.gain[ampName] = np.nan
933 dataset.gainErr[ampName] = np.nan
934 dataset.noise[ampName] = np.nan
935 dataset.noiseErr[ampName] = np.nan
936 dataset.ptcFitPars[ampName] = np.nan
937 dataset.ptcFitParsError[ampName] = np.nan
938 dataset.ptcFitReducedChiSquared[ampName] = np.nan
941 mask = mask & goodPoints
943 if ptcFitType ==
'EXPAPPROXIMATION':
945 parsIniPtc = [-1e-9, 1.0, 10.]
947 if ptcFitType ==
'POLYNOMIAL':
948 ptcFunc = funcPolynomial
954 while count <= maxIterationsPtcOutliers:
958 meanTempVec = meanVecOriginal[mask]
959 varTempVec = varVecOriginal[mask]
960 res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec))
966 sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal)
967 newMask = np.array([
True if np.abs(r) < sigmaCutPtcOutliers
else False for r
in sigResids])
968 mask = mask & newMask
969 if not (mask.any()
and newMask.any()):
970 msg = (f
"\nSERIOUS: All points in either mask: {mask} or newMask: {newMask} are bad. "
971 f
"Setting {ampName} to BAD.")
975 dataset.badAmps.append(ampName)
976 dataset.gain[ampName] = np.nan
977 dataset.gainErr[ampName] = np.nan
978 dataset.noise[ampName] = np.nan
979 dataset.noiseErr[ampName] = np.nan
980 dataset.ptcFitPars[ampName] = np.nan
981 dataset.ptcFitParsError[ampName] = np.nan
982 dataset.ptcFitReducedChiSquared[ampName] = np.nan
984 nDroppedTotal = Counter(mask)[
False]
985 self.log.debug(f
"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}")
988 assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal))
990 if not (mask.any()
and newMask.any()):
992 dataset.expIdMask[ampName] = mask
994 meanVecFinal = meanVecOriginal[mask]
995 varVecFinal = varVecOriginal[mask]
997 if Counter(mask)[
False] > 0:
998 self.log.info((f
"Number of points discarded in PTC of amplifier {ampName}:" +
999 f
" {Counter(mask)[False]} out of {len(meanVecOriginal)}"))
1001 if (len(meanVecFinal) < len(parsIniPtc)):
1002 msg = (f
"\nSERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of"
1003 f
"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.")
1007 dataset.badAmps.append(ampName)
1008 dataset.gain[ampName] = np.nan
1009 dataset.gainErr[ampName] = np.nan
1010 dataset.noise[ampName] = np.nan
1011 dataset.noiseErr[ampName] = np.nan
1012 dataset.ptcFitPars[ampName] = np.nan
1013 dataset.ptcFitParsError[ampName] = np.nan
1014 dataset.ptcFitReducedChiSquared[ampName] = np.nan
1018 if self.config.doFitBootstrap:
1019 parsFit, parsFitErr, reducedChiSqPtc =
fitBootstrap(parsIniPtc, meanVecFinal,
1020 varVecFinal, ptcFunc,
1021 weightsY=1./np.sqrt(varVecFinal))
1023 parsFit, parsFitErr, reducedChiSqPtc =
fitLeastSq(parsIniPtc, meanVecFinal,
1024 varVecFinal, ptcFunc,
1025 weightsY=1./np.sqrt(varVecFinal))
1026 dataset.ptcFitPars[ampName] = parsFit
1027 dataset.ptcFitParsError[ampName] = parsFitErr
1028 dataset.ptcFitReducedChiSquared[ampName] = reducedChiSqPtc
1030 if ptcFitType ==
'EXPAPPROXIMATION':
1031 ptcGain = parsFit[1]
1032 ptcGainErr = parsFitErr[1]
1033 ptcNoise = np.sqrt(np.fabs(parsFit[2]))
1034 ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2]))
1035 if ptcFitType ==
'POLYNOMIAL':
1036 ptcGain = 1./parsFit[1]
1037 ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1])
1038 ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain
1039 ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain
1040 dataset.gain[ampName] = ptcGain
1041 dataset.gainErr[ampName] = ptcGainErr
1042 dataset.noise[ampName] = ptcNoise
1043 dataset.noiseErr[ampName] = ptcNoiseErr
1044 if not len(dataset.ptcFitType) == 0:
1045 dataset.ptcFitType = ptcFitType