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

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
36from .photodiode import getBOTphotodiodeData
38from lsst.pipe.tasks.getRepositoryData import DataRefListRunner
39from lsst.ip.isr import PhotonTransferCurveDataset
41import copy
43__all__ = ['MeasurePhotonTransferCurveTask',
44 'MeasurePhotonTransferCurveTaskConfig']
47class MeasurePhotonTransferCurveTaskConfig(pexConfig.Config):
48 """Config class for photon transfer curve measurement task"""
49 ccdKey = pexConfig.Field(
50 dtype=str,
51 doc="The key by which to pull a detector from a dataId, e.g. 'ccd' or 'detector'.",
52 default='ccd',
53 )
54 ptcFitType = pexConfig.ChoiceField(
55 dtype=str,
56 doc="Fit PTC to Eq. 16, Eq. 20 in Astier+19, or to a polynomial.",
57 default="POLYNOMIAL",
58 allowed={
59 "POLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegree' to set 'n').",
60 "EXPAPPROXIMATION": "Approximation in Astier+19 (Eq. 16).",
61 "FULLCOVARIANCE": "Full covariances model in Astier+19 (Eq. 20)"
62 }
63 )
64 maximumRangeCovariancesAstier = pexConfig.Field(
65 dtype=int,
66 doc="Maximum range of covariances as in Astier+19",
67 default=8,
68 )
69 covAstierRealSpace = pexConfig.Field(
70 dtype=bool,
71 doc="Calculate covariances in real space or via FFT? (see appendix A of Astier+19).",
72 default=False,
73 )
74 polynomialFitDegree = pexConfig.Field(
75 dtype=int,
76 doc="Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.",
77 default=3,
78 )
79 linearity = pexConfig.ConfigurableField(
80 target=LinearitySolveTask,
81 doc="Task to solve the linearity."
82 )
84 doCreateLinearizer = pexConfig.Field(
85 dtype=bool,
86 doc="Calculate non-linearity and persist linearizer?",
87 default=False,
88 )
90 binSize = pexConfig.Field(
91 dtype=int,
92 doc="Bin the image by this factor in both dimensions.",
93 default=1,
94 )
95 minMeanSignal = pexConfig.DictField(
96 keytype=str,
97 itemtype=float,
98 doc="Minimum values (inclusive) of mean signal (in ADU) above which to consider, per amp."
99 " The same cut is applied to all amps if this dictionary is of the form"
100 " {'ALL_AMPS': value}",
101 default={'ALL_AMPS': 0.0},
102 )
103 maxMeanSignal = pexConfig.DictField(
104 keytype=str,
105 itemtype=float,
106 doc="Maximum values (inclusive) of mean signal (in ADU) below which to consider, per amp."
107 " The same cut is applied to all amps if this dictionary is of the form"
108 " {'ALL_AMPS': value}",
109 default={'ALL_AMPS': 1e6},
110 )
111 initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField(
112 dtype=float,
113 doc="Initially exclude data points with a variance that are more than a factor of this from being"
114 " linear in the positive direction, from the PTC fit. Note that these points will also be"
115 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
116 " to allow an accurate determination of the sigmas for said iterative fit.",
117 default=0.12,
118 min=0.0,
119 max=1.0,
120 )
121 initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField(
122 dtype=float,
123 doc="Initially exclude data points with a variance that are more than a factor of this from being"
124 " linear in the negative direction, from the PTC fit. Note that these points will also be"
125 " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
126 " to allow an accurate determination of the sigmas for said iterative fit.",
127 default=0.25,
128 min=0.0,
129 max=1.0,
130 )
131 sigmaCutPtcOutliers = pexConfig.Field(
132 dtype=float,
133 doc="Sigma cut for outlier rejection in PTC.",
134 default=5.0,
135 )
136 maskNameList = pexConfig.ListField(
137 dtype=str,
138 doc="Mask list to exclude from statistics calculations.",
139 default=['SUSPECT', 'BAD', 'NO_DATA'],
140 )
141 nSigmaClipPtc = pexConfig.Field(
142 dtype=float,
143 doc="Sigma cut for afwMath.StatisticsControl()",
144 default=5.5,
145 )
146 nIterSigmaClipPtc = pexConfig.Field(
147 dtype=int,
148 doc="Number of sigma-clipping iterations for afwMath.StatisticsControl()",
149 default=1,
150 )
151 maxIterationsPtcOutliers = pexConfig.Field(
152 dtype=int,
153 doc="Maximum number of iterations for outlier rejection in PTC.",
154 default=2,
155 )
156 doFitBootstrap = pexConfig.Field(
157 dtype=bool,
158 doc="Use bootstrap for the PTC fit parameters and errors?.",
159 default=False,
160 )
161 doPhotodiode = pexConfig.Field(
162 dtype=bool,
163 doc="Apply a correction based on the photodiode readings if available?",
164 default=True,
165 )
166 photodiodeDataPath = pexConfig.Field(
167 dtype=str,
168 doc="Gen2 only: path to locate the data photodiode data files.",
169 default=""
170 )
171 instrumentName = pexConfig.Field(
172 dtype=str,
173 doc="Instrument name.",
174 default='',
175 )
178class MeasurePhotonTransferCurveTask(pipeBase.CmdLineTask):
179 """A class to calculate, fit, and plot a PTC from a set of flat pairs.
181 The Photon Transfer Curve (var(signal) vs mean(signal)) is a standard tool
182 used in astronomical detectors characterization (e.g., Janesick 2001,
183 Janesick 2007). If ptcFitType is "EXPAPPROXIMATION" or "POLYNOMIAL", this task calculates the
184 PTC from a series of pairs of flat-field images; each pair taken at identical exposure
185 times. The difference image of each pair is formed to eliminate fixed pattern noise,
186 and then the variance of the difference image and the mean of the average image
187 are used to produce the PTC. An n-degree polynomial or the approximation in Equation
188 16 of Astier+19 ("The Shape of the Photon Transfer Curve of CCD sensors",
189 arXiv:1905.08677) can be fitted to the PTC curve. These models include
190 parameters such as the gain (e/DN) and readout noise.
192 Linearizers to correct for signal-chain non-linearity are also calculated.
193 The `Linearizer` class, in general, can support per-amp linearizers, but in this
194 task this is not supported.
196 If ptcFitType is "FULLCOVARIANCE", the covariances of the difference images are calculated via the
197 DFT methods described in Astier+19 and the variances for the PTC are given by the cov[0,0] elements
198 at each signal level. The full model in Equation 20 of Astier+19 is fit to the PTC to get the gain
199 and the noise.
201 Parameters
202 ----------
204 *args: `list`
205 Positional arguments passed to the Task constructor. None used at this
206 time.
207 **kwargs: `dict`
208 Keyword arguments passed on to the Task constructor. None used at this
209 time.
211 """
213 RunnerClass = DataRefListRunner
214 ConfigClass = MeasurePhotonTransferCurveTaskConfig
215 _DefaultName = "measurePhotonTransferCurve"
217 def __init__(self, *args, **kwargs):
218 pipeBase.CmdLineTask.__init__(self, *args, **kwargs)
219 self.makeSubtask("linearity")
220 plt.interactive(False) # stop windows popping up when plotting. When headless, use 'agg' backend too
221 self.config.validate()
222 self.config.freeze()
224 @pipeBase.timeMethod
225 def runDataRef(self, dataRefList):
226 """Run the Photon Transfer Curve (PTC) measurement task.
228 For a dataRef (which is each detector here),
229 and given a list of exposure pairs (postISR) at different exposure times,
230 measure the PTC.
232 Parameters
233 ----------
234 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
235 Data references for exposures for detectors to process.
236 """
237 if len(dataRefList) < 2:
238 raise RuntimeError("Insufficient inputs to combine.")
240 # setup necessary objects
241 dataRef = dataRefList[0]
243 detNum = dataRef.dataId[self.config.ccdKey]
244 camera = dataRef.get('camera')
245 detector = camera[dataRef.dataId[self.config.ccdKey]]
247 amps = detector.getAmplifiers()
248 ampNames = [amp.getName() for amp in amps]
249 datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType)
251 # Get the pairs of flat indexed by expTime
252 expPairs = self.makePairs(dataRefList)
253 expIds = []
254 for (exp1, exp2) in expPairs.values():
255 id1 = exp1.getInfo().getVisitInfo().getExposureId()
256 id2 = exp2.getInfo().getVisitInfo().getExposureId()
257 expIds.append((id1, id2))
258 self.log.info(f"Measuring PTC using {expIds} exposures for detector {detector.getId()}")
260 # get photodiode data early so that logic can be put in to only use the
261 # data if all files are found, as partial corrections are not possible
262 # or at least require significant logic to deal with
263 if self.config.doPhotodiode:
264 for (expId1, expId2) in expIds:
265 charges = [-1, -1] # necessary to have a not-found value to keep lists in step
266 for i, expId in enumerate([expId1, expId2]):
267 # //1000 is a Gen2 only hack, working around the fact an
268 # exposure's ID is not the same as the expId in the
269 # registry. Currently expId is concatenated with the
270 # zero-padded detector ID. This will all go away in Gen3.
271 dataRef.dataId['expId'] = expId//1000
272 if self.config.photodiodeDataPath:
273 photodiodeData = getBOTphotodiodeData(dataRef, self.config.photodiodeDataPath)
274 else:
275 photodiodeData = getBOTphotodiodeData(dataRef)
276 if photodiodeData: # default path stored in function def to keep task clean
277 charges[i] = photodiodeData.getCharge()
278 else:
279 # full expId (not //1000) here, as that encodes the
280 # the detector number as so is fully qualifying
281 self.log.warn(f"No photodiode data found for {expId}")
283 for ampName in ampNames:
284 datasetPtc.photoCharge[ampName].append((charges[0], charges[1]))
285 else:
286 # Can't be an empty list, as initialized, because astropy.Table won't allow it
287 # when saving as fits
288 for ampName in ampNames:
289 datasetPtc.photoCharge[ampName] = np.repeat(np.nan, len(expIds))
291 for ampName in ampNames:
292 datasetPtc.inputExpIdPairs[ampName] = expIds
294 maxMeanSignalDict = {ampName: 1e6 for ampName in ampNames}
295 minMeanSignalDict = {ampName: 0.0 for ampName in ampNames}
296 for ampName in ampNames:
297 if 'ALL_AMPS' in self.config.maxMeanSignal:
298 maxMeanSignalDict[ampName] = self.config.maxMeanSignal['ALL_AMPS']
299 elif ampName in self.config.maxMeanSignal:
300 maxMeanSignalDict[ampName] = self.config.maxMeanSignal[ampName]
302 if 'ALL_AMPS' in self.config.minMeanSignal:
303 minMeanSignalDict[ampName] = self.config.minMeanSignal['ALL_AMPS']
304 elif ampName in self.config.minMeanSignal:
305 minMeanSignalDict[ampName] = self.config.minMeanSignal[ampName]
307 tupleRecords = []
308 allTags = []
309 for expTime, (exp1, exp2) in expPairs.items():
310 expId1 = exp1.getInfo().getVisitInfo().getExposureId()
311 expId2 = exp2.getInfo().getVisitInfo().getExposureId()
312 tupleRows = []
313 nAmpsNan = 0
314 tags = ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName']
315 for ampNumber, amp in enumerate(detector):
316 ampName = amp.getName()
317 # covAstier: (i, j, var (cov[0,0]), cov, npix)
318 doRealSpace = self.config.covAstierRealSpace
319 muDiff, varDiff, covAstier = self.measureMeanVarCov(exp1, exp2, region=amp.getBBox(),
320 covAstierRealSpace=doRealSpace)
322 if np.isnan(muDiff) or np.isnan(varDiff) or (covAstier is None):
323 msg = (f"NaN mean or var, or None cov in amp {ampName} in exposure pair {expId1},"
324 f" {expId2} of detector {detNum}.")
325 self.log.warn(msg)
326 nAmpsNan += 1
327 continue
328 if (muDiff <= minMeanSignalDict[ampName]) or (muDiff >= maxMeanSignalDict[ampName]):
329 continue
331 datasetPtc.rawExpTimes[ampName].append(expTime)
332 datasetPtc.rawMeans[ampName].append(muDiff)
333 datasetPtc.rawVars[ampName].append(varDiff)
335 tupleRows += [(muDiff, ) + covRow + (ampNumber, expTime, ampName) for covRow in covAstier]
336 if nAmpsNan == len(ampNames):
337 msg = f"NaN mean in all amps of exposure pair {expId1}, {expId2} of detector {detNum}."
338 self.log.warn(msg)
339 continue
340 allTags += tags
341 tupleRecords += tupleRows
342 covariancesWithTags = np.core.records.fromrecords(tupleRecords, names=allTags)
344 for ampName in datasetPtc.ampNames:
345 # Sort raw vectors by rawMeans index
346 index = np.argsort(datasetPtc.rawMeans[ampName])
347 datasetPtc.rawExpTimes[ampName] = np.array(datasetPtc.rawExpTimes[ampName])[index]
348 datasetPtc.rawMeans[ampName] = np.array(datasetPtc.rawMeans[ampName])[index]
349 datasetPtc.rawVars[ampName] = np.array(datasetPtc.rawVars[ampName])[index]
351 if self.config.ptcFitType in ["FULLCOVARIANCE", ]:
352 # Calculate covariances and fit them, including the PTC, to Astier+19 full model (Eq. 20)
353 # First, fit get the flat pairs that are masked, according to the regular PTC (C_00 vs mu)
354 # The points at these fluxes will also be masked when calculating the other covariances, C_ij)
355 newDatasetPtc = copy.copy(datasetPtc)
356 newDatasetPtc = self.fitPtc(newDatasetPtc, 'EXPAPPROXIMATION')
357 for ampName in datasetPtc.ampNames:
358 datasetPtc.expIdMask[ampName] = newDatasetPtc.expIdMask[ampName]
360 datasetPtc.fitType = "FULLCOVARIANCE"
361 datasetPtc = self.fitCovariancesAstier(datasetPtc, covariancesWithTags)
362 elif self.config.ptcFitType in ["EXPAPPROXIMATION", "POLYNOMIAL"]:
363 # Fit the PTC to a polynomial or to Astier+19 exponential approximation (Eq. 16)
364 # Fill up PhotonTransferCurveDataset object.
365 datasetPtc = self.fitPtc(datasetPtc, self.config.ptcFitType)
367 detName = detector.getName()
368 now = datetime.datetime.utcnow()
369 calibDate = now.strftime("%Y-%m-%d")
370 butler = dataRef.getButler()
372 datasetPtc.updateMetadata(setDate=True, camera=camera, detector=detector)
374 # Fit a poynomial to calculate non-linearity and persist linearizer.
375 if self.config.doCreateLinearizer:
376 # Fit (non)linearity of signal vs time curve.
377 # Fill up PhotonTransferCurveDataset object.
378 # Fill up array for LUT linearizer (tableArray).
379 # Produce coefficients for Polynomial and Squared linearizers.
380 # Build linearizer objects.
381 dimensions = {'camera': camera.getName(), 'detector': detector.getId()}
382 linearityResults = self.linearity.run(datasetPtc, camera, dimensions)
383 linearizer = linearityResults.outputLinearizer
385 self.log.info("Writing linearizer:")
387 detName = detector.getName()
388 now = datetime.datetime.utcnow()
389 calibDate = now.strftime("%Y-%m-%d")
391 butler.put(linearizer, datasetType='linearizer',
392 dataId={'detector': detNum, 'detectorName': detName, 'calibDate': calibDate})
394 self.log.info(f"Writing PTC data.")
395 butler.put(datasetPtc, datasetType='photonTransferCurveDataset', dataId={'detector': detNum,
396 'detectorName': detName, 'calibDate': calibDate})
398 return pipeBase.Struct(exitStatus=0)
400 def makePairs(self, dataRefList):
401 """Produce a list of flat pairs indexed by exposure time.
403 Parameters
404 ----------
405 dataRefList : `list` [`lsst.daf.peristence.ButlerDataRef`]
406 Data references for exposures for detectors to process.
408 Return
409 ------
410 flatPairs : `dict` [`float`, `lsst.afw.image.exposure.exposure.ExposureF`]
411 Dictionary that groups flat-field exposures that have the same exposure time (seconds).
413 Notes
414 -----
415 We use the difference of one pair of flat-field images taken at the same exposure time when
416 calculating the PTC to reduce Fixed Pattern Noise. If there are > 2 flat-field images with the
417 same exposure time, the first two are kept and the rest discarded.
418 """
420 # Organize exposures by observation date.
421 expDict = {}
422 for dataRef in dataRefList:
423 try:
424 tempFlat = dataRef.get("postISRCCD")
425 except RuntimeError:
426 self.log.warn("postISR exposure could not be retrieved. Ignoring flat.")
427 continue
428 expDate = tempFlat.getInfo().getVisitInfo().getDate().get()
429 expDict.setdefault(expDate, tempFlat)
430 sortedExps = {k: expDict[k] for k in sorted(expDict)}
432 flatPairs = {}
433 for exp in sortedExps:
434 tempFlat = sortedExps[exp]
435 expTime = tempFlat.getInfo().getVisitInfo().getExposureTime()
436 listAtExpTime = flatPairs.setdefault(expTime, [])
437 if len(listAtExpTime) >= 2:
438 self.log.warn(f"Already found 2 exposures at expTime {expTime}. "
439 f"Ignoring exposure {tempFlat.getInfo().getVisitInfo().getExposureId()}")
440 else:
441 listAtExpTime.append(tempFlat)
443 keysToDrop = []
444 for (key, value) in flatPairs.items():
445 if len(value) < 2:
446 keysToDrop.append(key)
448 if len(keysToDrop):
449 for key in keysToDrop:
450 self.log.warn(f"Only one exposure found at expTime {key}. Dropping exposure "
451 f"{flatPairs[key][0].getInfo().getVisitInfo().getExposureId()}.")
452 flatPairs.pop(key)
453 sortedFlatPairs = {k: flatPairs[k] for k in sorted(flatPairs)}
454 return sortedFlatPairs
456 def fitCovariancesAstier(self, dataset, covariancesWithTagsArray):
457 """Fit measured flat covariances to full model in Astier+19.
459 Parameters
460 ----------
461 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
462 The dataset containing information such as the means, variances and exposure times.
464 covariancesWithTagsArray : `numpy.recarray`
465 Tuple with at least (mu, cov, var, i, j, npix), where:
466 mu : 0.5*(m1 + m2), where:
467 mu1: mean value of flat1
468 mu2: mean value of flat2
469 cov: covariance value at lag(i, j)
470 var: variance(covariance value at lag(0, 0))
471 i: lag dimension
472 j: lag dimension
473 npix: number of pixels used for covariance calculation.
475 Returns
476 -------
477 dataset: `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
478 This is the same dataset as the input paramter, however, it has been modified
479 to include information such as the fit vectors and the fit parameters. See
480 the class `PhotonTransferCurveDatase`.
481 """
483 covFits, covFitsNoB = fitData(covariancesWithTagsArray,
484 r=self.config.maximumRangeCovariancesAstier,
485 expIdMask=dataset.expIdMask)
486 dataset = self.getOutputPtcDataCovAstier(dataset, covFits, covFitsNoB)
487 return dataset
489 def getOutputPtcDataCovAstier(self, dataset, covFits, covFitsNoB):
490 """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects.
492 Parameters
493 ----------
494 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
495 The dataset containing information such as the means, variances and exposure times.
497 covFits: `dict`
498 Dictionary of CovFit objects, with amp names as keys.
500 covFitsNoB : `dict`
501 Dictionary of CovFit objects, with amp names as keys, and 'b=0' in Eq. 20 of Astier+19.
503 Returns
504 -------
505 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
506 This is the same dataset as the input paramter, however, it has been modified
507 to include extra information such as the mask 1D array, gains, reoudout noise, measured signal,
508 measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier.
509 See the class `PhotonTransferCurveDatase`.
510 """
511 assert(len(covFits) == len(covFitsNoB))
513 for i, amp in enumerate(dataset.ampNames):
514 lenInputTimes = len(dataset.rawExpTimes[amp])
515 # Not used when ptcFitType is 'FULLCOVARIANCE'
516 dataset.ptcFitPars[amp] = np.nan
517 dataset.ptcFitParsError[amp] = np.nan
518 dataset.ptcFitChiSq[amp] = np.nan
519 if (amp in covFits and (covFits[amp].covParams is not None) and
520 (covFitsNoB[amp].covParams is not None)):
521 fit = covFits[amp]
522 fitNoB = covFitsNoB[amp]
523 # Save full covariances, covariances models, and their weights
524 dataset.covariances[amp] = fit.cov
525 dataset.covariancesModel[amp] = fit.evalCovModel()
526 dataset.covariancesSqrtWeights[amp] = fit.sqrtW
527 dataset.aMatrix[amp] = fit.getA()
528 dataset.bMatrix[amp] = fit.getB()
529 dataset.covariancesNoB[amp] = fitNoB.cov
530 dataset.covariancesModelNoB[amp] = fitNoB.evalCovModel()
531 dataset.covariancesSqrtWeightsNoB[amp] = fitNoB.sqrtW
532 dataset.aMatrixNoB[amp] = fitNoB.getA()
534 (meanVecFinal, varVecFinal, varVecModel,
535 wc, varMask) = fit.getFitData(0, 0, divideByMu=False)
536 gain = fit.getGain()
538 dataset.gain[amp] = gain
539 dataset.gainErr[amp] = fit.getGainErr()
540 dataset.noise[amp] = np.sqrt(fit.getRon())
541 dataset.noiseErr[amp] = fit.getRonErr()
543 padLength = lenInputTimes - len(varVecFinal)
544 dataset.finalVars[amp] = np.pad(varVecFinal/(gain**2), (0, padLength), 'constant',
545 constant_values=np.nan)
546 dataset.finalModelVars[amp] = np.pad(varVecModel/(gain**2), (0, padLength), 'constant',
547 constant_values=np.nan)
548 dataset.finalMeans[amp] = np.pad(meanVecFinal/gain, (0, padLength), 'constant',
549 constant_values=np.nan)
550 else:
551 # Bad amp
552 # Entries need to have proper dimensions so read/write with astropy.Table works.
553 matrixSide = self.config.maximumRangeCovariancesAstier
554 nanMatrix = np.full((matrixSide, matrixSide), np.nan)
555 listNanMatrix = np.full((lenInputTimes, matrixSide, matrixSide), np.nan)
557 dataset.covariances[amp] = listNanMatrix
558 dataset.covariancesModel[amp] = listNanMatrix
559 dataset.covariancesSqrtWeights[amp] = listNanMatrix
560 dataset.aMatrix[amp] = nanMatrix
561 dataset.bMatrix[amp] = nanMatrix
562 dataset.covariancesNoB[amp] = listNanMatrix
563 dataset.covariancesModelNoB[amp] = listNanMatrix
564 dataset.covariancesSqrtWeightsNoB[amp] = listNanMatrix
565 dataset.aMatrixNoB[amp] = nanMatrix
567 dataset.expIdMask[amp] = np.repeat(np.nan, lenInputTimes)
568 dataset.gain[amp] = np.nan
569 dataset.gainErr[amp] = np.nan
570 dataset.noise[amp] = np.nan
571 dataset.noiseErr[amp] = np.nan
572 dataset.finalVars[amp] = np.repeat(np.nan, lenInputTimes)
573 dataset.finalModelVars[amp] = np.repeat(np.nan, lenInputTimes)
574 dataset.finalMeans[amp] = np.repeat(np.nan, lenInputTimes)
576 return dataset
578 def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False):
579 """Calculate the mean of each of two exposures and the variance and covariance of their difference.
581 The variance is calculated via afwMath, and the covariance via the methods in Astier+19 (appendix A).
582 In theory, var = covariance[0,0]. This should be validated, and in the future, we may decide to just
583 keep one (covariance).
585 Parameters
586 ----------
587 exposure1 : `lsst.afw.image.exposure.exposure.ExposureF`
588 First exposure of flat field pair.
590 exposure2 : `lsst.afw.image.exposure.exposure.ExposureF`
591 Second exposure of flat field pair.
593 region : `lsst.geom.Box2I`, optional
594 Region of each exposure where to perform the calculations (e.g, an amplifier).
596 covAstierRealSpace : `bool`, optional
597 Should the covariannces in Astier+19 be calculated in real space or via FFT?
598 See Appendix A of Astier+19.
600 Returns
601 -------
602 mu : `float` or `NaN`
603 0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in
604 both exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
606 varDiff : `float` or `NaN`
607 Half of the clipped variance of the difference of the regions inthe two input
608 exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
610 covDiffAstier : `list` or `None`
611 List with tuples of the form (dx, dy, var, cov, npix), where:
612 dx : `int`
613 Lag in x
614 dy : `int`
615 Lag in y
616 var : `float`
617 Variance at (dx, dy).
618 cov : `float`
619 Covariance at (dx, dy).
620 nPix : `int`
621 Number of pixel pairs used to evaluate var and cov.
622 If either mu1 or m2 are NaN's, the returned value is None.
623 """
625 if region is not None:
626 im1Area = exposure1.maskedImage[region]
627 im2Area = exposure2.maskedImage[region]
628 else:
629 im1Area = exposure1.maskedImage
630 im2Area = exposure2.maskedImage
632 if self.config.binSize > 1:
633 im1Area = afwMath.binImage(im1Area, self.config.binSize)
634 im2Area = afwMath.binImage(im2Area, self.config.binSize)
636 im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList)
637 im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
638 self.config.nIterSigmaClipPtc,
639 im1MaskVal)
640 im1StatsCtrl.setNanSafe(True)
641 im1StatsCtrl.setAndMask(im1MaskVal)
643 im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList)
644 im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
645 self.config.nIterSigmaClipPtc,
646 im2MaskVal)
647 im2StatsCtrl.setNanSafe(True)
648 im2StatsCtrl.setAndMask(im2MaskVal)
650 # Clipped mean of images; then average of mean.
651 mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue()
652 mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue()
653 if np.isnan(mu1) or np.isnan(mu2):
654 return np.nan, np.nan, None
655 mu = 0.5*(mu1 + mu2)
657 # Take difference of pairs
658 # symmetric formula: diff = (mu2*im1-mu1*im2)/(0.5*(mu1+mu2))
659 temp = im2Area.clone()
660 temp *= mu1
661 diffIm = im1Area.clone()
662 diffIm *= mu2
663 diffIm -= temp
664 diffIm /= mu
666 diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList)
667 diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
668 self.config.nIterSigmaClipPtc,
669 diffImMaskVal)
670 diffImStatsCtrl.setNanSafe(True)
671 diffImStatsCtrl.setAndMask(diffImMaskVal)
673 varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue())
675 # Get the mask and identify good pixels as '1', and the rest as '0'.
676 w1 = np.where(im1Area.getMask().getArray() == 0, 1, 0)
677 w2 = np.where(im2Area.getMask().getArray() == 0, 1, 0)
679 w12 = w1*w2
680 wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0)
681 w = w12*wDiff
683 if np.sum(w) == 0:
684 return np.nan, np.nan, None
686 maxRangeCov = self.config.maximumRangeCovariancesAstier
687 if covAstierRealSpace:
688 covDiffAstier = computeCovDirect(diffIm.getImage().getArray(), w, maxRangeCov)
689 else:
690 shapeDiff = diffIm.getImage().getArray().shape
691 fftShape = (fftSize(shapeDiff[0] + maxRangeCov), fftSize(shapeDiff[1]+maxRangeCov))
692 c = CovFft(diffIm.getImage().getArray(), w, fftShape, maxRangeCov)
693 covDiffAstier = c.reportCovFft(maxRangeCov)
695 return mu, varDiff, covDiffAstier
697 def computeCovDirect(self, diffImage, weightImage, maxRange):
698 """Compute covariances of diffImage in real space.
700 For lags larger than ~25, it is slower than the FFT way.
701 Taken from https://github.com/PierreAstier/bfptc/
703 Parameters
704 ----------
705 diffImage : `numpy.array`
706 Image to compute the covariance of.
708 weightImage : `numpy.array`
709 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
711 maxRange : `int`
712 Last index of the covariance to be computed.
714 Returns
715 -------
716 outList : `list`
717 List with tuples of the form (dx, dy, var, cov, npix), where:
718 dx : `int`
719 Lag in x
720 dy : `int`
721 Lag in y
722 var : `float`
723 Variance at (dx, dy).
724 cov : `float`
725 Covariance at (dx, dy).
726 nPix : `int`
727 Number of pixel pairs used to evaluate var and cov.
728 """
729 outList = []
730 var = 0
731 # (dy,dx) = (0,0) has to be first
732 for dy in range(maxRange + 1):
733 for dx in range(0, maxRange + 1):
734 if (dx*dy > 0):
735 cov1, nPix1 = self.covDirectValue(diffImage, weightImage, dx, dy)
736 cov2, nPix2 = self.covDirectValue(diffImage, weightImage, dx, -dy)
737 cov = 0.5*(cov1 + cov2)
738 nPix = nPix1 + nPix2
739 else:
740 cov, nPix = self.covDirectValue(diffImage, weightImage, dx, dy)
741 if (dx == 0 and dy == 0):
742 var = cov
743 outList.append((dx, dy, var, cov, nPix))
745 return outList
747 def covDirectValue(self, diffImage, weightImage, dx, dy):
748 """Compute covariances of diffImage in real space at lag (dx, dy).
750 Taken from https://github.com/PierreAstier/bfptc/ (c.f., appendix of Astier+19).
752 Parameters
753 ----------
754 diffImage : `numpy.array`
755 Image to compute the covariance of.
757 weightImage : `numpy.array`
758 Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
760 dx : `int`
761 Lag in x.
763 dy : `int`
764 Lag in y.
766 Returns
767 -------
768 cov : `float`
769 Covariance at (dx, dy)
771 nPix : `int`
772 Number of pixel pairs used to evaluate var and cov.
773 """
774 (nCols, nRows) = diffImage.shape
775 # switching both signs does not change anything:
776 # it just swaps im1 and im2 below
777 if (dx < 0):
778 (dx, dy) = (-dx, -dy)
779 # now, we have dx >0. We have to distinguish two cases
780 # depending on the sign of dy
781 if dy >= 0:
782 im1 = diffImage[dy:, dx:]
783 w1 = weightImage[dy:, dx:]
784 im2 = diffImage[:nCols - dy, :nRows - dx]
785 w2 = weightImage[:nCols - dy, :nRows - dx]
786 else:
787 im1 = diffImage[:nCols + dy, dx:]
788 w1 = weightImage[:nCols + dy, dx:]
789 im2 = diffImage[-dy:, :nRows - dx]
790 w2 = weightImage[-dy:, :nRows - dx]
791 # use the same mask for all 3 calculations
792 wAll = w1*w2
793 # do not use mean() because weightImage=0 pixels would then count
794 nPix = wAll.sum()
795 im1TimesW = im1*wAll
796 s1 = im1TimesW.sum()/nPix
797 s2 = (im2*wAll).sum()/nPix
798 p = (im1TimesW*im2).sum()/nPix
799 cov = p - s1*s2
801 return cov, nPix
803 @staticmethod
804 def _initialParsForPolynomial(order):
805 assert(order >= 2)
806 pars = np.zeros(order, dtype=np.float)
807 pars[0] = 10
808 pars[1] = 1
809 pars[2:] = 0.0001
810 return pars
812 @staticmethod
813 def _boundsForPolynomial(initialPars, lowers=[], uppers=[]):
814 if not len(lowers):
815 lowers = [np.NINF for p in initialPars]
816 if not len(uppers):
817 uppers = [np.inf for p in initialPars]
818 lowers[1] = 0 # no negative gains
819 return (lowers, uppers)
821 @staticmethod
822 def _boundsForAstier(initialPars, lowers=[], uppers=[]):
823 if not len(lowers):
824 lowers = [np.NINF for p in initialPars]
825 if not len(uppers):
826 uppers = [np.inf for p in initialPars]
827 return (lowers, uppers)
829 @staticmethod
830 def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative):
831 """Return a boolean array to mask bad points.
833 Parameters
834 ----------
835 means : `numpy.array`
836 Input array with mean signal values.
838 variances : `numpy.array`
839 Input array with variances at each mean value.
841 maxDeviationPositive : `float`
842 Maximum deviation from being constant for the variance/mean
843 ratio, in the positive direction.
845 maxDeviationNegative : `float`
846 Maximum deviation from being constant for the variance/mean
847 ratio, in the negative direction.
849 Return
850 ------
851 goodPoints : `numpy.array` [`bool`]
852 Boolean array to select good (`True`) and bad (`False`)
853 points.
855 Notes
856 -----
857 A linear function has a constant ratio, so find the median
858 value of the ratios, and exclude the points that deviate
859 from that by more than a factor of maxDeviationPositive/negative.
860 Asymmetric deviations are supported as we expect the PTC to turn
861 down as the flux increases, but sometimes it anomalously turns
862 upwards just before turning over, which ruins the fits, so it
863 is wise to be stricter about restricting positive outliers than
864 negative ones.
866 Too high and points that are so bad that fit will fail will be included
867 Too low and the non-linear points will be excluded, biasing the NL fit.
869 This function also masks points after the variance starts decreasing.
870 """
872 assert(len(means) == len(variances))
873 ratios = [b/a for (a, b) in zip(means, variances)]
874 medianRatio = np.nanmedian(ratios)
875 ratioDeviations = [(r/medianRatio)-1 for r in ratios]
877 # so that it doesn't matter if the deviation is expressed as positive or negative
878 maxDeviationPositive = abs(maxDeviationPositive)
879 maxDeviationNegative = -1. * abs(maxDeviationNegative)
881 goodPoints = np.array([True if (r < maxDeviationPositive and r > maxDeviationNegative)
882 else False for r in ratioDeviations])
884 # Discard points when variance starts decreasing
885 pivot = np.where(np.array(np.diff(variances)) < 0)[0]
886 if len(pivot) == 0:
887 return goodPoints
888 else:
889 pivot = np.min(pivot)
890 goodPoints[pivot:len(goodPoints)] = False
891 return goodPoints
893 def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9):
894 """"""
895 nBad = Counter(array)[0]
896 if nBad == 0:
897 return array
899 if warn:
900 msg = f"Found {nBad} zeros in array at elements {[x for x in np.where(array==0)[0]]}"
901 self.log.warn(msg)
903 array[array == 0] = substituteValue
904 return array
906 def fitPtc(self, dataset, ptcFitType):
907 """Fit the photon transfer curve to a polynimial or to Astier+19 approximation.
909 Fit the photon transfer curve with either a polynomial of the order
910 specified in the task config, or using the Astier approximation.
912 Sigma clipping is performed iteratively for the fit, as well as an
913 initial clipping of data points that are more than
914 config.initialNonLinearityExclusionThreshold away from lying on a
915 straight line. This other step is necessary because the photon transfer
916 curve turns over catastrophically at very high flux (because saturation
917 drops the variance to ~0) and these far outliers cause the initial fit
918 to fail, meaning the sigma cannot be calculated to perform the
919 sigma-clipping.
921 Parameters
922 ----------
923 dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
924 The dataset containing the means, variances and exposure times
926 ptcFitType : `str`
927 Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or
928 'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC
930 Returns
931 -------
932 dataset: `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
933 This is the same dataset as the input paramter, however, it has been modified
934 to include information such as the fit vectors and the fit parameters. See
935 the class `PhotonTransferCurveDatase`.
936 """
938 matrixSide = self.config.maximumRangeCovariancesAstier
939 nanMatrix = np.empty((matrixSide, matrixSide))
940 nanMatrix[:] = np.nan
942 for amp in dataset.ampNames:
943 lenInputTimes = len(dataset.rawExpTimes[amp])
944 listNanMatrix = np.empty((lenInputTimes, matrixSide, matrixSide))
945 listNanMatrix[:] = np.nan
947 dataset.covariances[amp] = listNanMatrix
948 dataset.covariancesModel[amp] = listNanMatrix
949 dataset.covariancesSqrtWeights[amp] = listNanMatrix
950 dataset.aMatrix[amp] = nanMatrix
951 dataset.bMatrix[amp] = nanMatrix
952 dataset.covariancesNoB[amp] = listNanMatrix
953 dataset.covariancesModelNoB[amp] = listNanMatrix
954 dataset.covariancesSqrtWeightsNoB[amp] = listNanMatrix
955 dataset.aMatrixNoB[amp] = nanMatrix
957 def errFunc(p, x, y):
958 return ptcFunc(p, x) - y
960 sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers
961 maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers
963 for i, ampName in enumerate(dataset.ampNames):
964 timeVecOriginal = np.array(dataset.rawExpTimes[ampName])
965 meanVecOriginal = np.array(dataset.rawMeans[ampName])
966 varVecOriginal = np.array(dataset.rawVars[ampName])
967 varVecOriginal = self._makeZeroSafe(varVecOriginal)
969 goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal,
970 self.config.initialNonLinearityExclusionThresholdPositive,
971 self.config.initialNonLinearityExclusionThresholdNegative)
972 if not (goodPoints.any()):
973 msg = (f"\nSERIOUS: All points in goodPoints: {goodPoints} are bad."
974 f"Setting {ampName} to BAD.")
975 self.log.warn(msg)
976 # The first and second parameters of initial fit are discarded (bias and gain)
977 # for the final NL coefficients
978 dataset.badAmps.append(ampName)
979 dataset.expIdMask[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
980 dataset.gain[ampName] = np.nan
981 dataset.gainErr[ampName] = np.nan
982 dataset.noise[ampName] = np.nan
983 dataset.noiseErr[ampName] = np.nan
984 dataset.ptcFitPars[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
985 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
986 dataset.ptcFitParsError[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
987 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
988 dataset.ptcFitChiSq[ampName] = np.nan
989 dataset.finalVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
990 dataset.finalModelVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
991 dataset.finalMeans[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
992 continue
994 mask = goodPoints
996 if ptcFitType == 'EXPAPPROXIMATION':
997 ptcFunc = funcAstier
998 parsIniPtc = [-1e-9, 1.0, 10.] # a00, gain, noise
999 # lowers and uppers obtained from studies by C. Lage (UC Davis, 11/2020).
1000 bounds = self._boundsForAstier(parsIniPtc, lowers=[-1e-4, 0.5, -100],
1001 uppers=[1e-4, 2.5, 100])
1002 if ptcFitType == 'POLYNOMIAL':
1003 ptcFunc = funcPolynomial
1004 parsIniPtc = self._initialParsForPolynomial(self.config.polynomialFitDegree + 1)
1005 bounds = self._boundsForPolynomial(parsIniPtc)
1007 # Before bootstrap fit, do an iterative fit to get rid of outliers
1008 count = 1
1009 while count <= maxIterationsPtcOutliers:
1010 # Note that application of the mask actually shrinks the array
1011 # to size rather than setting elements to zero (as we want) so
1012 # always update mask itself and re-apply to the original data
1013 meanTempVec = meanVecOriginal[mask]
1014 varTempVec = varVecOriginal[mask]
1015 res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec))
1016 pars = res.x
1018 # change this to the original from the temp because the masks are ANDed
1019 # meaning once a point is masked it's always masked, and the masks must
1020 # always be the same length for broadcasting
1021 sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal)
1022 newMask = np.array([True if np.abs(r) < sigmaCutPtcOutliers else False for r in sigResids])
1023 mask = mask & newMask
1024 if not (mask.any() and newMask.any()):
1025 msg = (f"\nSERIOUS: All points in either mask: {mask} or newMask: {newMask} are bad. "
1026 f"Setting {ampName} to BAD.")
1027 self.log.warn(msg)
1028 # The first and second parameters of initial fit are discarded (bias and gain)
1029 # for the final NL coefficients
1030 dataset.badAmps.append(ampName)
1031 dataset.expIdMask[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1032 dataset.gain[ampName] = np.nan
1033 dataset.gainErr[ampName] = np.nan
1034 dataset.noise[ampName] = np.nan
1035 dataset.noiseErr[ampName] = np.nan
1036 dataset.ptcFitPars[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1)
1037 if ptcFitType in ["POLYNOMIAL", ] else
1038 np.repeat(np.nan, 3))
1039 dataset.ptcFitParsError[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1)
1040 if ptcFitType in ["POLYNOMIAL", ] else
1041 np.repeat(np.nan, 3))
1042 dataset.ptcFitChiSq[ampName] = np.nan
1043 dataset.finalVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1044 dataset.finalModelVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1045 dataset.finalMeans[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1046 break
1047 nDroppedTotal = Counter(mask)[False]
1048 self.log.debug(f"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}")
1049 count += 1
1050 # objects should never shrink
1051 assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal))
1053 if not (mask.any() and newMask.any()):
1054 continue
1055 dataset.expIdMask[ampName] = mask # store the final mask
1056 parsIniPtc = pars
1057 meanVecFinal = meanVecOriginal[mask]
1058 varVecFinal = varVecOriginal[mask]
1060 if Counter(mask)[False] > 0:
1061 self.log.info((f"Number of points discarded in PTC of amplifier {ampName}:" +
1062 f" {Counter(mask)[False]} out of {len(meanVecOriginal)}"))
1064 if (len(meanVecFinal) < len(parsIniPtc)):
1065 msg = (f"\nSERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of"
1066 f"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.")
1067 self.log.warn(msg)
1068 # The first and second parameters of initial fit are discarded (bias and gain)
1069 # for the final NL coefficients
1070 dataset.badAmps.append(ampName)
1071 dataset.expIdMask[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1072 dataset.gain[ampName] = np.nan
1073 dataset.gainErr[ampName] = np.nan
1074 dataset.noise[ampName] = np.nan
1075 dataset.noiseErr[ampName] = np.nan
1076 dataset.ptcFitPars[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
1077 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
1078 dataset.ptcFitParsError[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
1079 ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
1080 dataset.ptcFitChiSq[ampName] = np.nan
1081 dataset.finalVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1082 dataset.finalModelVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1083 dataset.finalMeans[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
1084 continue
1086 # Fit the PTC
1087 if self.config.doFitBootstrap:
1088 parsFit, parsFitErr, reducedChiSqPtc = fitBootstrap(parsIniPtc, meanVecFinal,
1089 varVecFinal, ptcFunc,
1090 weightsY=1./np.sqrt(varVecFinal))
1091 else:
1092 parsFit, parsFitErr, reducedChiSqPtc = fitLeastSq(parsIniPtc, meanVecFinal,
1093 varVecFinal, ptcFunc,
1094 weightsY=1./np.sqrt(varVecFinal))
1095 dataset.ptcFitPars[ampName] = parsFit
1096 dataset.ptcFitParsError[ampName] = parsFitErr
1097 dataset.ptcFitChiSq[ampName] = reducedChiSqPtc
1098 # Masked variances (measured and modeled) and means. Need to pad the array so astropy.Table does
1099 # not crash (the mask may vary per amp).
1100 padLength = len(dataset.rawExpTimes[ampName]) - len(varVecFinal)
1101 dataset.finalVars[ampName] = np.pad(varVecFinal, (0, padLength), 'constant',
1102 constant_values=np.nan)
1103 dataset.finalModelVars[ampName] = np.pad(ptcFunc(parsFit, meanVecFinal), (0, padLength),
1104 'constant', constant_values=np.nan)
1105 dataset.finalMeans[ampName] = np.pad(meanVecFinal, (0, padLength), 'constant',
1106 constant_values=np.nan)
1108 if ptcFitType == 'EXPAPPROXIMATION':
1109 ptcGain = parsFit[1]
1110 ptcGainErr = parsFitErr[1]
1111 ptcNoise = np.sqrt(np.fabs(parsFit[2]))
1112 ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2]))
1113 if ptcFitType == 'POLYNOMIAL':
1114 ptcGain = 1./parsFit[1]
1115 ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1])
1116 ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain
1117 ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain
1118 dataset.gain[ampName] = ptcGain
1119 dataset.gainErr[ampName] = ptcGainErr
1120 dataset.noise[ampName] = ptcNoise
1121 dataset.noiseErr[ampName] = ptcNoiseErr
1122 if not len(dataset.ptcFitType) == 0:
1123 dataset.ptcFitType = ptcFitType
1124 if len(dataset.badAmps) == 0:
1125 dataset.badAmps = np.repeat(np.nan, len(list(dataset.rawExpTimes.values())[0]))
1127 return dataset