Coverage for python/lsst/cp/pipe/ptc.py : 10%

Hot-keys on this page
r m x p toggle line displays
j k next/prev highlighted chunk
0 (zero) top of page
1 (one) first highlighted chunk
1# This file is part of cp_pipe.
2#
3# Developed for the LSST Data Management System.
4# This product includes software developed by the LSST Project
5# (https://www.lsst.org).
6# See the COPYRIGHT file at the top-level directory of this distribution
7# for details of code ownership.
8#
9# This program is free software: you can redistribute it and/or modify
10# it under the terms of the GNU General Public License as published by
11# the Free Software Foundation, either version 3 of the License, or
12# (at your option) any later version.
13#
14# This program is distributed in the hope that it will be useful,
15# but WITHOUT ANY WARRANTY; without even the implied warranty of
16# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17# GNU General Public License for more details.
18#
19# You should have received a copy of the GNU General Public License
20# along with this program. If not, see <https://www.gnu.org/licenses/>.
21#
22import numpy as np
23import matplotlib.pyplot as plt
24from collections import Counter
26import lsst.afw.math as afwMath
27import lsst.pex.config as pexConfig
28import lsst.pipe.base as pipeBase
29from .utils import (fitLeastSq, fitBootstrap, funcPolynomial, funcAstier)
30from scipy.optimize import least_squares
32import datetime
34from .astierCovPtcUtils import (fftSize, CovFft, computeCovDirect, fitData)
35from .linearity import LinearitySolveTask
37from lsst.pipe.tasks.getRepositoryData import DataRefListRunner
39__all__ = ['MeasurePhotonTransferCurveTask',
40 'MeasurePhotonTransferCurveTaskConfig',
41 'PhotonTransferCurveDataset']
44class MeasurePhotonTransferCurveTaskConfig(pexConfig.Config):
45 """Config class for photon transfer curve measurement task"""
46 ccdKey = pexConfig.Field(
47 dtype=str,
48 doc="The key by which to pull a detector from a dataId, e.g. 'ccd' or 'detector'.",
49 default='ccd',
50 )
51 ptcFitType = pexConfig.ChoiceField(
52 dtype=str,
53 doc="Fit PTC to approximation in Astier+19 (Equation 16) or to a polynomial.",
54 default="POLYNOMIAL",
55 allowed={
56 "POLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegree' to set 'n').",
57 "EXPAPPROXIMATION": "Approximation in Astier+19 (Eq. 16).",
58 "FULLCOVARIANCE": "Full covariances model in Astier+19 (Eq. 20)"
59 }
60 )
61 sigmaClipFullFitCovariancesAstier = pexConfig.Field(
62 dtype=float,
63 doc="sigma clip for full model fit for FULLCOVARIANCE ptcFitType ",
64 default=5.0,
65 )
66 maxIterFullFitCovariancesAstier = pexConfig.Field(
67 dtype=int,
68 doc="Maximum number of iterations in full model fit for FULLCOVARIANCE ptcFitType",
69 default=3,
70 )
71 maximumRangeCovariancesAstier = pexConfig.Field(
72 dtype=int,
73 doc="Maximum range of covariances as in Astier+19",
74 default=8,
75 )
76 covAstierRealSpace = pexConfig.Field(
77 dtype=bool,
78 doc="Calculate covariances in real space or via FFT? (see appendix A of Astier+19).",
79 default=False,
80 )
81 polynomialFitDegree = pexConfig.Field(
82 dtype=int,
83 doc="Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.",
84 default=3,
85 )
86 linearity = pexConfig.ConfigurableField(
87 target=LinearitySolveTask,
88 doc="Task to solve the linearity."
89 )
91 doCreateLinearizer = pexConfig.Field(
92 dtype=bool,
93 doc="Calculate non-linearity and persist linearizer?",
94 default=False,
95 )
97 binSize = pexConfig.Field(
98 dtype=int,
99 doc="Bin the image by this factor in both dimensions.",
100 default=1,
101 )
102 minMeanSignal = pexConfig.Field(
103 dtype=float,
104 doc="Minimum value (inclusive) of mean signal (in DN) above which to consider.",
105 default=0,
106 )
107 maxMeanSignal = pexConfig.Field(
108 dtype=float,
109 doc="Maximum value (inclusive) of mean signal (in DN) below which to consider.",
110 default=9e6,
111 )
112 initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField(
113 dtype=float,
114 doc="Initially exclude data points with a variance that are more than a factor of this from being"
115 " linear in the positive direction, from the PTC fit. Note that these points will also be"
116 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
117 " to allow an accurate determination of the sigmas for said iterative fit.",
118 default=0.12,
119 min=0.0,
120 max=1.0,
121 )
122 initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField(
123 dtype=float,
124 doc="Initially exclude data points with a variance that are more than a factor of this from being"
125 " linear in the negative direction, from the PTC fit. Note that these points will also be"
126 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
127 " to allow an accurate determination of the sigmas for said iterative fit.",
128 default=0.25,
129 min=0.0,
130 max=1.0,
131 )
132 sigmaCutPtcOutliers = pexConfig.Field(
133 dtype=float,
134 doc="Sigma cut for outlier rejection in PTC.",
135 default=5.0,
136 )
137 maskNameList = pexConfig.ListField(
138 dtype=str,
139 doc="Mask list to exclude from statistics calculations.",
140 default=['SUSPECT', 'BAD', 'NO_DATA'],
141 )
142 nSigmaClipPtc = pexConfig.Field(
143 dtype=float,
144 doc="Sigma cut for afwMath.StatisticsControl()",
145 default=5.5,
146 )
147 nIterSigmaClipPtc = pexConfig.Field(
148 dtype=int,
149 doc="Number of sigma-clipping iterations for afwMath.StatisticsControl()",
150 default=1,
151 )
152 maxIterationsPtcOutliers = pexConfig.Field(
153 dtype=int,
154 doc="Maximum number of iterations for outlier rejection in PTC.",
155 default=2,
156 )
157 doFitBootstrap = pexConfig.Field(
158 dtype=bool,
159 doc="Use bootstrap for the PTC fit parameters and errors?.",
160 default=False,
161 )
162 instrumentName = pexConfig.Field(
163 dtype=str,
164 doc="Instrument name.",
165 default='',
166 )
169class PhotonTransferCurveDataset:
170 """A simple class to hold the output data from the PTC task.
172 The dataset is made up of a dictionary for each item, keyed by the
173 amplifiers' names, which much be supplied at construction time.
175 New items cannot be added to the class to save accidentally saving to the
176 wrong property, and the class can be frozen if desired.
178 inputExpIdPairs records the exposures used to produce the data.
179 When fitPtc() or fitCovariancesAstier() is run, a mask is built up, which is by definition
180 always the same length as inputExpIdPairs, rawExpTimes, rawMeans
181 and rawVars, and is a list of bools, which are incrementally set to False
182 as points are discarded from the fits.
184 PTC fit parameters for polynomials are stored in a list in ascending order
185 of polynomial term, i.e. par[0]*x^0 + par[1]*x + par[2]*x^2 etc
186 with the length of the list corresponding to the order of the polynomial
187 plus one.
189 Parameters
190 ----------
191 ampNames : `list`
192 List with the names of the amplifiers of the detector at hand.
194 ptcFitType : `str`
195 Type of model fitted to the PTC: "POLYNOMIAL", "EXPAPPROXIMATION", or "FULLCOVARIANCE".
197 Returns
198 -------
199 `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
200 Output dataset from MeasurePhotonTransferCurveTask.
201 """
203 def __init__(self, ampNames, ptcFitType):
204 # add items to __dict__ directly because __setattr__ is overridden
206 # instance variables
207 self.__dict__["ptcFitType"] = ptcFitType
208 self.__dict__["ampNames"] = ampNames
209 self.__dict__["badAmps"] = []
211 # raw data variables
212 # expIdMask is the mask produced after outlier rejection. The mask produced by "FULLCOVARIANCE"
213 # may differ from the one produced in the other two PTC fit types.
214 self.__dict__["inputExpIdPairs"] = {ampName: [] for ampName in ampNames}
215 self.__dict__["expIdMask"] = {ampName: [] for ampName in ampNames}
216 self.__dict__["rawExpTimes"] = {ampName: [] for ampName in ampNames}
217 self.__dict__["rawMeans"] = {ampName: [] for ampName in ampNames}
218 self.__dict__["rawVars"] = {ampName: [] for ampName in ampNames}
220 # Gain and noise
221 self.__dict__["gain"] = {ampName: -1. for ampName in ampNames}
222 self.__dict__["gainErr"] = {ampName: -1. for ampName in ampNames}
223 self.__dict__["noise"] = {ampName: -1. for ampName in ampNames}
224 self.__dict__["noiseErr"] = {ampName: -1. for ampName in ampNames}
226 # if ptcFitTye in ["POLYNOMIAL", "EXPAPPROXIMATION"]
227 # fit information
228 self.__dict__["ptcFitPars"] = {ampName: [] for ampName in ampNames}
229 self.__dict__["ptcFitParsError"] = {ampName: [] for ampName in ampNames}
230 self.__dict__["ptcFitReducedChiSquared"] = {ampName: [] for ampName in ampNames}
232 # if ptcFitTye in ["FULLCOVARIANCE"]
233 # "covariancesTuple" is a numpy recarray with entries of the form
234 # ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName']
235 # "covariancesFits" has CovFit objects that fit the measured covariances to Eq. 20 of Astier+19.
236 # In "covariancesFitsWithNoB", "b"=0 in the model described by Eq. 20 of Astier+19.
237 self.__dict__["covariancesTuple"] = {ampName: [] for ampName in ampNames}
238 self.__dict__["covariancesFitsWithNoB"] = {ampName: [] for ampName in ampNames}
239 self.__dict__["covariancesFits"] = {ampName: [] for ampName in ampNames}
240 self.__dict__["aMatrix"] = {ampName: [] for ampName in ampNames}
241 self.__dict__["bMatrix"] = {ampName: [] for ampName in ampNames}
243 # "final" means that the "raw" vectors above had "expIdMask" applied.
244 self.__dict__["finalVars"] = {ampName: [] for ampName in ampNames}
245 self.__dict__["finalModelVars"] = {ampName: [] for ampName in ampNames}
246 self.__dict__["finalMeans"] = {ampName: [] for ampName in ampNames}
248 def __setattr__(self, attribute, value):
249 """Protect class attributes"""
250 if attribute not in self.__dict__:
251 raise AttributeError(f"{attribute} is not already a member of PhotonTransferCurveDataset, which"
252 " does not support setting of new attributes.")
253 else:
254 self.__dict__[attribute] = value
256 def getExpIdsUsed(self, ampName):
257 """Get the exposures used, i.e. not discarded, for a given amp.
259 If no mask has been created yet, all exposures are returned.
260 """
261 if len(self.expIdMask[ampName]) == 0:
262 return self.inputExpIdPairs[ampName]
264 # if the mask exists it had better be the same length as the expIdPairs
265 assert len(self.expIdMask[ampName]) == len(self.inputExpIdPairs[ampName])
267 pairs = self.inputExpIdPairs[ampName]
268 mask = self.expIdMask[ampName]
269 # cast to bool required because numpy
270 return [(exp1, exp2) for ((exp1, exp2), m) in zip(pairs, mask) if bool(m) is True]
272 def getGoodAmps(self):
273 return [amp for amp in self.ampNames if amp not in self.badAmps]
276class MeasurePhotonTransferCurveTask(pipeBase.CmdLineTask):
277 """A class to calculate, fit, and plot a PTC from a set of flat pairs.
279 The Photon Transfer Curve (var(signal) vs mean(signal)) is a standard tool
280 used in astronomical detectors characterization (e.g., Janesick 2001,
281 Janesick 2007). If ptcFitType is "EXPAPPROXIMATION" or "POLYNOMIAL", this task calculates the
282 PTC from a series of pairs of flat-field images; each pair taken at identical exposure
283 times. The difference image of each pair is formed to eliminate fixed pattern noise,
284 and then the variance of the difference image and the mean of the average image
285 are used to produce the PTC. An n-degree polynomial or the approximation in Equation
286 16 of Astier+19 ("The Shape of the Photon Transfer Curve of CCD sensors",
287 arXiv:1905.08677) can be fitted to the PTC curve. These models include
288 parameters such as the gain (e/DN) and readout noise.
290 Linearizers to correct for signal-chain non-linearity are also calculated.
291 The `Linearizer` class, in general, can support per-amp linearizers, but in this
292 task this is not supported.
294 If ptcFitType is "FULLCOVARIANCE", the covariances of the difference images are calculated via the
295 DFT methods described in Astier+19 and the variances for the PTC are given by the cov[0,0] elements
296 at each signal level. The full model in Equation 20 of Astier+19 is fit to the PTC to get the gain
297 and the noise.
299 Parameters
300 ----------
302 *args: `list`
303 Positional arguments passed to the Task constructor. None used at this
304 time.
305 **kwargs: `dict`
306 Keyword arguments passed on to the Task constructor. None used at this
307 time.
309 """
311 RunnerClass = DataRefListRunner
312 ConfigClass = MeasurePhotonTransferCurveTaskConfig
313 _DefaultName = "measurePhotonTransferCurve"
315 def __init__(self, *args, **kwargs):
316 pipeBase.CmdLineTask.__init__(self, *args, **kwargs)
317 self.makeSubtask("linearity")
318 plt.interactive(False) # stop windows popping up when plotting. When headless, use 'agg' backend too
319 self.config.validate()
320 self.config.freeze()
322 @pipeBase.timeMethod
323 def runDataRef(self, dataRefList):
324 """Run the Photon Transfer Curve (PTC) measurement task.
326 For a dataRef (which is each detector here),
327 and given a list of exposure pairs (postISR) at different exposure times,
328 measure the PTC.
330 Parameters
331 ----------
332 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
333 Data references for exposures for detectors to process.
334 """
335 if len(dataRefList) < 2:
336 raise RuntimeError("Insufficient inputs to combine.")
338 # setup necessary objects
339 dataRef = dataRefList[0]
341 detNum = dataRef.dataId[self.config.ccdKey]
342 camera = dataRef.get('camera')
343 detector = camera[dataRef.dataId[self.config.ccdKey]]
345 amps = detector.getAmplifiers()
346 ampNames = [amp.getName() for amp in amps]
347 datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType)
349 # Get the pairs of flat indexed by expTime
350 expPairs = self.makePairs(dataRefList)
351 expIds = []
352 for (exp1, exp2) in expPairs.values():
353 id1 = exp1.getInfo().getVisitInfo().getExposureId()
354 id2 = exp2.getInfo().getVisitInfo().getExposureId()
355 expIds.append((id1, id2))
356 self.log.info(f"Measuring PTC using {expIds} exposures for detector {detector.getId()}")
357 tupleRecords = []
358 allTags = []
359 for expTime, (exp1, exp2) in expPairs.items():
360 expId1 = exp1.getInfo().getVisitInfo().getExposureId()
361 expId2 = exp2.getInfo().getVisitInfo().getExposureId()
362 tupleRows = []
363 nAmpsNan = 0
364 for ampNumber, amp in enumerate(detector):
365 ampName = amp.getName()
366 # covAstier: (i, j, var (cov[0,0]), cov, npix)
367 doRealSpace = self.config.covAstierRealSpace
368 muDiff, varDiff, covAstier = self.measureMeanVarCov(exp1, exp2, region=amp.getBBox(),
369 covAstierRealSpace=doRealSpace)
370 if np.isnan(muDiff) or np.isnan(varDiff) or (covAstier is None):
371 msg = (f"NaN mean or var, or None cov in amp {ampNumber} in exposure pair {expId1},"
372 f" {expId2} of detector {detNum}.")
373 self.log.warn(msg)
374 nAmpsNan += 1
375 continue
376 tags = ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName']
377 if (muDiff <= self.config.minMeanSignal) or (muDiff >= self.config.maxMeanSignal):
378 continue
379 datasetPtc.rawExpTimes[ampName].append(expTime)
380 datasetPtc.rawMeans[ampName].append(muDiff)
381 datasetPtc.rawVars[ampName].append(varDiff)
382 datasetPtc.inputExpIdPairs[ampName].append((expId1, expId2))
384 tupleRows += [(muDiff, ) + covRow + (ampNumber, expTime, ampName) for covRow in covAstier]
385 if nAmpsNan == len(ampNames):
386 msg = f"NaN mean in all amps of exposure pair {expId1}, {expId2} of detector {detNum}."
387 self.log.warn(msg)
388 continue
389 allTags += tags
390 tupleRecords += tupleRows
391 covariancesWithTags = np.core.records.fromrecords(tupleRecords, names=allTags)
393 if self.config.ptcFitType in ["FULLCOVARIANCE", ]:
394 # Calculate covariances and fit them, including the PTC, to Astier+19 full model (Eq. 20)
395 datasetPtc = self.fitCovariancesAstier(datasetPtc, covariancesWithTags)
396 elif self.config.ptcFitType in ["EXPAPPROXIMATION", "POLYNOMIAL"]:
397 # Fit the PTC to a polynomial or to Astier+19 exponential approximation (Eq. 16)
398 # Fill up PhotonTransferCurveDataset object.
399 datasetPtc = self.fitPtc(datasetPtc, self.config.ptcFitType)
401 # Fit a poynomial to calculate non-linearity and persist linearizer.
402 if self.config.doCreateLinearizer:
403 # Fit (non)linearity of signal vs time curve.
404 # Fill up PhotonTransferCurveDataset object.
405 # Fill up array for LUT linearizer (tableArray).
406 # Produce coefficients for Polynomial ans Squared linearizers.
407 # Build linearizer objects.
408 dimensions = {'camera': camera.getName(), 'detector': detector.getId()}
409 linearityResults = self.linearity.run(datasetPtc, camera, dimensions)
410 linearizer = linearityResults.outputLinearizer
412 butler = dataRef.getButler()
413 self.log.info("Writing linearizer:")
415 detName = detector.getName()
416 now = datetime.datetime.utcnow()
417 calibDate = now.strftime("%Y-%m-%d")
419 butler.put(linearizer, datasetType='Linearizer', dataId={'detector': detNum,
420 'detectorName': detName, 'calibDate': calibDate})
421 self.log.info(f"Writing PTC data.")
422 dataRef.put(datasetPtc, datasetType="photonTransferCurveDataset")
424 return pipeBase.Struct(exitStatus=0)
426 def makePairs(self, dataRefList):
427 """Produce a list of flat pairs indexed by exposure time.
429 Parameters
430 ----------
431 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
432 Data references for exposures for detectors to process.
434 Return
435 ------
436 flatPairs : `dict` [`float`, `lsst.afw.image.exposure.exposure.ExposureF`]
437 Dictionary that groups flat-field exposures that have the same exposure time (seconds).
439 Notes
440 -----
441 We use the difference of one pair of flat-field images taken at the same exposure time when
442 calculating the PTC to reduce Fixed Pattern Noise. If there are > 2 flat-field images with the
443 same exposure time, the first two are kept and the rest discarded.
444 """
446 # Organize exposures by observation date.
447 expDict = {}
448 for dataRef in dataRefList:
449 try:
450 tempFlat = dataRef.get("postISRCCD")
451 except RuntimeError:
452 self.log.warn(f"postISR exposure could not be retrieved. Ignoring flat.")
453 continue
454 expDate = tempFlat.getInfo().getVisitInfo().getDate().get()
455 expDict.setdefault(expDate, tempFlat)
456 sortedExps = {k: expDict[k] for k in sorted(expDict)}
458 flatPairs = {}
459 for exp in sortedExps:
460 tempFlat = sortedExps[exp]
461 expTime = tempFlat.getInfo().getVisitInfo().getExposureTime()
462 listAtExpTime = flatPairs.setdefault(expTime, [])
463 if len(listAtExpTime) < 2:
464 listAtExpTime.append(tempFlat)
465 if len(listAtExpTime) > 2:
466 self.log.warn("More than 2 exposures found at expTime {expTime}. Dropping exposures "
467 f"{listAtExpTime[2:]}.")
469 for (key, value) in flatPairs.items():
470 if len(value) < 2:
471 flatPairs.pop(key)
472 self.log.warn("Only one exposure found at expTime {key}. Dropping exposure {value}.")
473 return flatPairs
475 def fitCovariancesAstier(self, dataset, covariancesWithTagsArray):
476 """Fit measured flat covariances to full model in Astier+19.
478 Parameters
479 ----------
480 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
481 The dataset containing information such as the means, variances and exposure times.
483 covariancesWithTagsArray : `numpy.recarray`
484 Tuple with at least (mu, cov, var, i, j, npix), where:
485 mu : 0.5*(m1 + m2), where:
486 mu1: mean value of flat1
487 mu2: mean value of flat2
488 cov: covariance value at lag(i, j)
489 var: variance(covariance value at lag(0, 0))
490 i: lag dimension
491 j: lag dimension
492 npix: number of pixels used for covariance calculation.
494 Returns
495 -------
496 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
497 This is the same dataset as the input paramter, however, it has been modified
498 to include information such as the fit vectors and the fit parameters. See
499 the class `PhotonTransferCurveDatase`.
500 """
502 covFits, covFitsNoB = fitData(covariancesWithTagsArray, maxMu=self.config.maxMeanSignal,
503 r=self.config.maximumRangeCovariancesAstier,
504 nSigmaFullFit=self.config.sigmaClipFullFitCovariancesAstier,
505 maxIterFullFit=self.config.maxIterFullFitCovariancesAstier)
507 dataset.covariancesTuple = covariancesWithTagsArray
508 dataset.covariancesFits = covFits
509 dataset.covariancesFitsWithNoB = covFitsNoB
510 dataset = self.getOutputPtcDataCovAstier(dataset, covFits)
512 return dataset
514 def getOutputPtcDataCovAstier(self, dataset, covFits):
515 """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects.
517 Parameters
518 ----------
519 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
520 The dataset containing information such as the means, variances and exposure times.
522 covFits: `dict`
523 Dictionary of CovFit objects, with amp names as keys.
525 Returns
526 -------
527 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
528 This is the same dataset as the input paramter, however, it has been modified
529 to include extra information such as the mask 1D array, gains, reoudout noise, measured signal,
530 measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier.
531 See the class `PhotonTransferCurveDatase`.
532 """
534 for i, amp in enumerate(covFits):
535 fit = covFits[amp]
536 (meanVecFinal, varVecFinal, varVecModel,
537 wc, varMask) = fit.getFitData(0, 0, divideByMu=False, returnMasked=True)
538 gain = fit.getGain()
539 dataset.expIdMask[amp] = varMask
540 dataset.gain[amp] = gain
541 dataset.gainErr[amp] = fit.getGainErr()
542 dataset.noise[amp] = np.sqrt(np.fabs(fit.getRon()))
543 dataset.noiseErr[amp] = fit.getRonErr()
544 dataset.finalVars[amp].append(varVecFinal/(gain**2))
545 dataset.finalModelVars[amp].append(varVecModel/(gain**2))
546 dataset.finalMeans[amp].append(meanVecFinal/gain)
547 dataset.aMatrix[amp].append(fit.getA())
548 dataset.bMatrix[amp].append(fit.getB())
550 return dataset
552 def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False):
553 """Calculate the mean of each of two exposures and the variance and covariance of their difference.
555 The variance is calculated via afwMath, and the covariance via the methods in Astier+19 (appendix A).
556 In theory, var = covariance[0,0]. This should be validated, and in the future, we may decide to just
557 keep one (covariance).
559 Parameters
560 ----------
561 exposure1 : `lsst.afw.image.exposure.exposure.ExposureF`
562 First exposure of flat field pair.
564 exposure2 : `lsst.afw.image.exposure.exposure.ExposureF`
565 Second exposure of flat field pair.
567 region : `lsst.geom.Box2I`, optional
568 Region of each exposure where to perform the calculations (e.g, an amplifier).
570 covAstierRealSpace : `bool`, optional
571 Should the covariannces in Astier+19 be calculated in real space or via FFT?
572 See Appendix A of Astier+19.
574 Returns
575 -------
576 mu : `float` or `NaN`
577 0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in
578 both exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
580 varDiff : `float` or `NaN`
581 Half of the clipped variance of the difference of the regions inthe two input
582 exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
584 covDiffAstier : `list` or `NaN`
585 List with tuples of the form (dx, dy, var, cov, npix), where:
586 dx : `int`
587 Lag in x
588 dy : `int`
589 Lag in y
590 var : `float`
591 Variance at (dx, dy).
592 cov : `float`
593 Covariance at (dx, dy).
594 nPix : `int`
595 Number of pixel pairs used to evaluate var and cov.
596 If either mu1 or m2 are NaN's, the returned value is NaN.
597 """
599 if region is not None:
600 im1Area = exposure1.maskedImage[region]
601 im2Area = exposure2.maskedImage[region]
602 else:
603 im1Area = exposure1.maskedImage
604 im2Area = exposure2.maskedImage
606 if self.config.binSize > 1:
607 im1Area = afwMath.binImage(im1Area, self.config.binSize)
608 im2Area = afwMath.binImage(im2Area, self.config.binSize)
610 im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList)
611 im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
612 self.config.nIterSigmaClipPtc,
613 im1MaskVal)
614 im1StatsCtrl.setNanSafe(True)
615 im1StatsCtrl.setAndMask(im1MaskVal)
617 im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList)
618 im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
619 self.config.nIterSigmaClipPtc,
620 im2MaskVal)
621 im2StatsCtrl.setNanSafe(True)
622 im2StatsCtrl.setAndMask(im2MaskVal)
624 # Clipped mean of images; then average of mean.
625 mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue()
626 mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue()
627 if np.isnan(mu1) or np.isnan(mu2):
628 return np.nan, np.nan, None
629 mu = 0.5*(mu1 + mu2)
631 # Take difference of pairs
632 # symmetric formula: diff = (mu2*im1-mu1*im2)/(0.5*(mu1+mu2))
633 temp = im2Area.clone()
634 temp *= mu1
635 diffIm = im1Area.clone()
636 diffIm *= mu2
637 diffIm -= temp
638 diffIm /= mu
640 diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList)
641 diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
642 self.config.nIterSigmaClipPtc,
643 diffImMaskVal)
644 diffImStatsCtrl.setNanSafe(True)
645 diffImStatsCtrl.setAndMask(diffImMaskVal)
647 varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue())
649 # Get the mask and identify good pixels as '1', and the rest as '0'.
650 w1 = np.where(im1Area.getMask().getArray() == 0, 1, 0)
651 w2 = np.where(im2Area.getMask().getArray() == 0, 1, 0)
653 w12 = w1*w2
654 wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0)
655 w = w12*wDiff
657 maxRangeCov = self.config.maximumRangeCovariancesAstier
658 if covAstierRealSpace:
659 covDiffAstier = computeCovDirect(diffIm.getImage().getArray(), w, maxRangeCov)
660 else:
661 shapeDiff = diffIm.getImage().getArray().shape
662 fftShape = (fftSize(shapeDiff[0] + maxRangeCov), fftSize(shapeDiff[1]+maxRangeCov))
663 c = CovFft(diffIm.getImage().getArray(), w, fftShape, maxRangeCov)
664 covDiffAstier = c.reportCovFft(maxRangeCov)
666 return mu, varDiff, covDiffAstier
668 def computeCovDirect(self, diffImage, weightImage, maxRange):
669 """Compute covariances of diffImage in real space.
671 For lags larger than ~25, it is slower than the FFT way.
672 Taken from https://github.com/PierreAstier/bfptc/
674 Parameters
675 ----------
676 diffImage : `numpy.array`
677 Image to compute the covariance of.
679 weightImage : `numpy.array`
680 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
682 maxRange : `int`
683 Last index of the covariance to be computed.
685 Returns
686 -------
687 outList : `list`
688 List with tuples of the form (dx, dy, var, cov, npix), where:
689 dx : `int`
690 Lag in x
691 dy : `int`
692 Lag in y
693 var : `float`
694 Variance at (dx, dy).
695 cov : `float`
696 Covariance at (dx, dy).
697 nPix : `int`
698 Number of pixel pairs used to evaluate var and cov.
699 """
700 outList = []
701 var = 0
702 # (dy,dx) = (0,0) has to be first
703 for dy in range(maxRange + 1):
704 for dx in range(0, maxRange + 1):
705 if (dx*dy > 0):
706 cov1, nPix1 = self.covDirectValue(diffImage, weightImage, dx, dy)
707 cov2, nPix2 = self.covDirectValue(diffImage, weightImage, dx, -dy)
708 cov = 0.5*(cov1 + cov2)
709 nPix = nPix1 + nPix2
710 else:
711 cov, nPix = self.covDirectValue(diffImage, weightImage, dx, dy)
712 if (dx == 0 and dy == 0):
713 var = cov
714 outList.append((dx, dy, var, cov, nPix))
716 return outList
718 def covDirectValue(self, diffImage, weightImage, dx, dy):
719 """Compute covariances of diffImage in real space at lag (dx, dy).
721 Taken from https://github.com/PierreAstier/bfptc/ (c.f., appendix of Astier+19).
723 Parameters
724 ----------
725 diffImage : `numpy.array`
726 Image to compute the covariance of.
728 weightImage : `numpy.array`
729 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
731 dx : `int`
732 Lag in x.
734 dy : `int`
735 Lag in y.
737 Returns
738 -------
739 cov : `float`
740 Covariance at (dx, dy)
742 nPix : `int`
743 Number of pixel pairs used to evaluate var and cov.
744 """
745 (nCols, nRows) = diffImage.shape
746 # switching both signs does not change anything:
747 # it just swaps im1 and im2 below
748 if (dx < 0):
749 (dx, dy) = (-dx, -dy)
750 # now, we have dx >0. We have to distinguish two cases
751 # depending on the sign of dy
752 if dy >= 0:
753 im1 = diffImage[dy:, dx:]
754 w1 = weightImage[dy:, dx:]
755 im2 = diffImage[:nCols - dy, :nRows - dx]
756 w2 = weightImage[:nCols - dy, :nRows - dx]
757 else:
758 im1 = diffImage[:nCols + dy, dx:]
759 w1 = weightImage[:nCols + dy, dx:]
760 im2 = diffImage[-dy:, :nRows - dx]
761 w2 = weightImage[-dy:, :nRows - dx]
762 # use the same mask for all 3 calculations
763 wAll = w1*w2
764 # do not use mean() because weightImage=0 pixels would then count
765 nPix = wAll.sum()
766 im1TimesW = im1*wAll
767 s1 = im1TimesW.sum()/nPix
768 s2 = (im2*wAll).sum()/nPix
769 p = (im1TimesW*im2).sum()/nPix
770 cov = p - s1*s2
772 return cov, nPix
774 @staticmethod
775 def _initialParsForPolynomial(order):
776 assert(order >= 2)
777 pars = np.zeros(order, dtype=np.float)
778 pars[0] = 10
779 pars[1] = 1
780 pars[2:] = 0.0001
781 return pars
783 @staticmethod
784 def _boundsForPolynomial(initialPars):
785 lowers = [np.NINF for p in initialPars]
786 uppers = [np.inf for p in initialPars]
787 lowers[1] = 0 # no negative gains
788 return (lowers, uppers)
790 @staticmethod
791 def _boundsForAstier(initialPars):
792 lowers = [np.NINF for p in initialPars]
793 uppers = [np.inf for p in initialPars]
794 return (lowers, uppers)
796 @staticmethod
797 def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative):
798 """Return a boolean array to mask bad points.
800 A linear function has a constant ratio, so find the median
801 value of the ratios, and exclude the points that deviate
802 from that by more than a factor of maxDeviationPositive/negative.
803 Asymmetric deviations are supported as we expect the PTC to turn
804 down as the flux increases, but sometimes it anomalously turns
805 upwards just before turning over, which ruins the fits, so it
806 is wise to be stricter about restricting positive outliers than
807 negative ones.
809 Too high and points that are so bad that fit will fail will be included
810 Too low and the non-linear points will be excluded, biasing the NL fit."""
811 ratios = [b/a for (a, b) in zip(means, variances)]
812 medianRatio = np.median(ratios)
813 ratioDeviations = [(r/medianRatio)-1 for r in ratios]
815 # so that it doesn't matter if the deviation is expressed as positive or negative
816 maxDeviationPositive = abs(maxDeviationPositive)
817 maxDeviationNegative = -1. * abs(maxDeviationNegative)
819 goodPoints = np.array([True if (r < maxDeviationPositive and r > maxDeviationNegative)
820 else False for r in ratioDeviations])
821 return goodPoints
823 def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9):
824 """"""
825 nBad = Counter(array)[0]
826 if nBad == 0:
827 return array
829 if warn:
830 msg = f"Found {nBad} zeros in array at elements {[x for x in np.where(array==0)[0]]}"
831 self.log.warn(msg)
833 array[array == 0] = substituteValue
834 return array
836 def fitPtc(self, dataset, ptcFitType):
837 """Fit the photon transfer curve to a polynimial or to Astier+19 approximation.
839 Fit the photon transfer curve with either a polynomial of the order
840 specified in the task config, or using the Astier approximation.
842 Sigma clipping is performed iteratively for the fit, as well as an
843 initial clipping of data points that are more than
844 config.initialNonLinearityExclusionThreshold away from lying on a
845 straight line. This other step is necessary because the photon transfer
846 curve turns over catastrophically at very high flux (because saturation
847 drops the variance to ~0) and these far outliers cause the initial fit
848 to fail, meaning the sigma cannot be calculated to perform the
849 sigma-clipping.
851 Parameters
852 ----------
853 dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
854 The dataset containing the means, variances and exposure times
856 ptcFitType : `str`
857 Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or
858 'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC
860 Returns
861 -------
862 dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
863 This is the same dataset as the input paramter, however, it has been modified
864 to include information such as the fit vectors and the fit parameters. See
865 the class `PhotonTransferCurveDatase`.
866 """
868 def errFunc(p, x, y):
869 return ptcFunc(p, x) - y
871 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers
872 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers
874 for i, ampName in enumerate(dataset.ampNames):
875 timeVecOriginal = np.array(dataset.rawExpTimes[ampName])
876 meanVecOriginal = np.array(dataset.rawMeans[ampName])
877 varVecOriginal = np.array(dataset.rawVars[ampName])
878 varVecOriginal = self._makeZeroSafe(varVecOriginal)
880 mask = ((meanVecOriginal >= self.config.minMeanSignal) &
881 (meanVecOriginal <= self.config.maxMeanSignal))
883 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal,
884 self.config.initialNonLinearityExclusionThresholdPositive,
885 self.config.initialNonLinearityExclusionThresholdNegative)
886 if not (mask.any() and goodPoints.any()):
887 msg = (f"\nSERIOUS: All points in either mask: {mask} or goodPoints: {goodPoints} are bad."
888 f"Setting {ampName} to BAD.")
889 self.log.warn(msg)
890 # The first and second parameters of initial fit are discarded (bias and gain)
891 # for the final NL coefficients
892 dataset.badAmps.append(ampName)
893 dataset.gain[ampName] = np.nan
894 dataset.gainErr[ampName] = np.nan
895 dataset.noise[ampName] = np.nan
896 dataset.noiseErr[ampName] = np.nan
897 dataset.ptcFitPars[ampName] = np.nan
898 dataset.ptcFitParsError[ampName] = np.nan
899 dataset.ptcFitReducedChiSquared[ampName] = np.nan
900 continue
902 mask = mask & goodPoints
904 if ptcFitType == 'EXPAPPROXIMATION':
905 ptcFunc = funcAstier
906 parsIniPtc = [-1e-9, 1.0, 10.] # a00, gain, noise
907 bounds = self._boundsForAstier(parsIniPtc)
908 if ptcFitType == 'POLYNOMIAL':
909 ptcFunc = funcPolynomial
910 parsIniPtc = self._initialParsForPolynomial(self.config.polynomialFitDegree + 1)
911 bounds = self._boundsForPolynomial(parsIniPtc)
913 # Before bootstrap fit, do an iterative fit to get rid of outliers
914 count = 1
915 while count <= maxIterationsPtcOutliers:
916 # Note that application of the mask actually shrinks the array
917 # to size rather than setting elements to zero (as we want) so
918 # always update mask itself and re-apply to the original data
919 meanTempVec = meanVecOriginal[mask]
920 varTempVec = varVecOriginal[mask]
921 res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec))
922 pars = res.x
924 # change this to the original from the temp because the masks are ANDed
925 # meaning once a point is masked it's always masked, and the masks must
926 # always be the same length for broadcasting
927 sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal)
928 newMask = np.array([True if np.abs(r) < sigmaCutPtcOutliers else False for r in sigResids])
929 mask = mask & newMask
930 if not (mask.any() and newMask.any()):
931 msg = (f"\nSERIOUS: All points in either mask: {mask} or newMask: {newMask} are bad. "
932 f"Setting {ampName} to BAD.")
933 self.log.warn(msg)
934 # The first and second parameters of initial fit are discarded (bias and gain)
935 # for the final NL coefficients
936 dataset.badAmps.append(ampName)
937 dataset.gain[ampName] = np.nan
938 dataset.gainErr[ampName] = np.nan
939 dataset.noise[ampName] = np.nan
940 dataset.noiseErr[ampName] = np.nan
941 dataset.ptcFitPars[ampName] = np.nan
942 dataset.ptcFitParsError[ampName] = np.nan
943 dataset.ptcFitReducedChiSquared[ampName] = np.nan
944 break
945 nDroppedTotal = Counter(mask)[False]
946 self.log.debug(f"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}")
947 count += 1
948 # objects should never shrink
949 assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal))
951 if not (mask.any() and newMask.any()):
952 continue
953 dataset.expIdMask[ampName] = mask # store the final mask
954 parsIniPtc = pars
955 meanVecFinal = meanVecOriginal[mask]
956 varVecFinal = varVecOriginal[mask]
958 if Counter(mask)[False] > 0:
959 self.log.info((f"Number of points discarded in PTC of amplifier {ampName}:" +
960 f" {Counter(mask)[False]} out of {len(meanVecOriginal)}"))
962 if (len(meanVecFinal) < len(parsIniPtc)):
963 msg = (f"\nSERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of"
964 f"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.")
965 self.log.warn(msg)
966 # The first and second parameters of initial fit are discarded (bias and gain)
967 # for the final NL coefficients
968 dataset.badAmps.append(ampName)
969 dataset.gain[ampName] = np.nan
970 dataset.gainErr[ampName] = np.nan
971 dataset.noise[ampName] = np.nan
972 dataset.noiseErr[ampName] = np.nan
973 dataset.ptcFitPars[ampName] = np.nan
974 dataset.ptcFitParsError[ampName] = np.nan
975 dataset.ptcFitReducedChiSquared[ampName] = np.nan
976 continue
978 # Fit the PTC
979 if self.config.doFitBootstrap:
980 parsFit, parsFitErr, reducedChiSqPtc = fitBootstrap(parsIniPtc, meanVecFinal,
981 varVecFinal, ptcFunc,
982 weightsY=1./np.sqrt(varVecFinal))
983 else:
984 parsFit, parsFitErr, reducedChiSqPtc = fitLeastSq(parsIniPtc, meanVecFinal,
985 varVecFinal, ptcFunc,
986 weightsY=1./np.sqrt(varVecFinal))
987 dataset.ptcFitPars[ampName] = parsFit
988 dataset.ptcFitParsError[ampName] = parsFitErr
989 dataset.ptcFitReducedChiSquared[ampName] = reducedChiSqPtc
991 if ptcFitType == 'EXPAPPROXIMATION':
992 ptcGain = parsFit[1]
993 ptcGainErr = parsFitErr[1]
994 ptcNoise = np.sqrt(np.fabs(parsFit[2]))
995 ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2]))
996 if ptcFitType == 'POLYNOMIAL':
997 ptcGain = 1./parsFit[1]
998 ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1])
999 ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain
1000 ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain
1001 dataset.gain[ampName] = ptcGain
1002 dataset.gainErr[ampName] = ptcGainErr
1003 dataset.noise[ampName] = ptcNoise
1004 dataset.noiseErr[ampName] = ptcNoiseErr
1005 if not len(dataset.ptcFitType) == 0:
1006 dataset.ptcFitType = ptcFitType
1008 return dataset