lsst.cp.pipe  21.0.0-7-gfd72ab2+fa9972ccfc
cpSolvePtcTask.py
Go to the documentation of this file.
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 #
22 import numpy as np
23 from collections import Counter
24 
25 import lsst.pex.config as pexConfig
26 import lsst.pipe.base as pipeBase
27 from lsst.cp.pipe.utils import (fitLeastSq, fitBootstrap, funcPolynomial, funcAstier)
28 
29 from scipy.optimize import least_squares
30 
32 
33 from .astierCovPtcUtils import fitData
34 
35 from lsst.ip.isr import PhotonTransferCurveDataset
36 
37 from lsst.cp.pipe._lookupStaticCalibration import lookupStaticCalibration
38 
39 import copy
40 
41 
42 __all__ = ['PhotonTransferCurveSolveConfig', 'PhotonTransferCurveSolveTask']
43 
44 
45 class PhotonTransferCurveSolveConnections(pipeBase.PipelineTaskConnections,
46  dimensions=("instrument", "detector")):
47  inputCovariances = cT.Input(
48  name="ptcCovariances",
49  doc="Tuple with measured covariances from flats.",
50  storageClass="PhotonTransferCurveDataset",
51  dimensions=("instrument", "exposure", "detector"),
52  multiple=True,
53  )
54  camera = cT.PrerequisiteInput(
55  name="camera",
56  doc="Camera the input data comes from.",
57  storageClass="Camera",
58  dimensions=("instrument",),
59  isCalibration=True,
60  lookupFunction=lookupStaticCalibration,
61  )
62  outputPtcDataset = cT.Output(
63  name="ptcDatsetProposal",
64  doc="Output proposed ptc dataset.",
65  storageClass="PhotonTransferCurveDataset",
66  dimensions=("instrument", "detector"),
67  multiple=False,
68  isCalibration=True,
69  )
70 
71 
72 class PhotonTransferCurveSolveConfig(pipeBase.PipelineTaskConfig,
73  pipelineConnections=PhotonTransferCurveSolveConnections):
74  """Configuration for fitting measured covariances.
75  """
76  ptcFitType = pexConfig.ChoiceField(
77  dtype=str,
78  doc="Fit PTC to Eq. 16, Eq. 20 in Astier+19, or to a polynomial.",
79  default="POLYNOMIAL",
80  allowed={
81  "POLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegree' to set 'n').",
82  "EXPAPPROXIMATION": "Approximation in Astier+19 (Eq. 16).",
83  "FULLCOVARIANCE": "Full covariances model in Astier+19 (Eq. 20)"
84  }
85  )
86  maximumRangeCovariancesAstier = pexConfig.Field(
87  dtype=int,
88  doc="Maximum range of covariances as in Astier+19",
89  default=8,
90  )
91  sigmaClipFullFitCovariancesAstier = pexConfig.Field(
92  dtype=float,
93  doc="sigma clip for full model fit for FULLCOVARIANCE ptcFitType ",
94  default=5.0,
95  )
96  maxIterFullFitCovariancesAstier = pexConfig.Field(
97  dtype=int,
98  doc="Maximum number of iterations in full model fit for FULLCOVARIANCE ptcFitType",
99  default=3,
100  )
101  polynomialFitDegree = pexConfig.Field(
102  dtype=int,
103  doc="Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.",
104  default=3,
105  )
106  sigmaCutPtcOutliers = pexConfig.Field(
107  dtype=float,
108  doc="Sigma cut for outlier rejection in PTC.",
109  default=5.0,
110  )
111  maxIterationsPtcOutliers = pexConfig.Field(
112  dtype=int,
113  doc="Maximum number of iterations for outlier rejection in PTC.",
114  default=2,
115  )
116  initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField(
117  dtype=float,
118  doc="Initially exclude data points with a variance that are more than a factor of this from being"
119  " linear in the positive direction, from the PTC fit. Note that these points will also be"
120  " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
121  " to allow an accurate determination of the sigmas for said iterative fit.",
122  default=0.05,
123  min=0.0,
124  max=1.0,
125  )
126  initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField(
127  dtype=float,
128  doc="Initially exclude data points with a variance that are more than a factor of this from being"
129  " linear in the negative direction, from the PTC fit. Note that these points will also be"
130  " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
131  " to allow an accurate determination of the sigmas for said iterative fit.",
132  default=0.25,
133  min=0.0,
134  max=1.0,
135  )
136  minMeanRatioTest = pexConfig.Field(
137  dtype=float,
138  doc="In the initial test to screen out bad points with a ratio test, points with low"
139  " flux can get inadvertantly screened. This test only screens out points with flux"
140  " above this value.",
141  default=20000,
142  )
143  minVarPivotSearch = pexConfig.Field(
144  dtype=float,
145  doc="The code looks for a pivot signal point after which the variance starts decreasing at high-flux"
146  " to exclude then form the PTC model fit. However, sometimes at low fluxes, the variance"
147  " decreases slightly. Set this variable for the variance value, in ADU^2, after which the pivot "
148  " should be sought.",
149  default=10000,
150  )
151  doFitBootstrap = pexConfig.Field(
152  dtype=bool,
153  doc="Use bootstrap for the PTC fit parameters and errors?.",
154  default=False,
155  )
156 
157 
158 class PhotonTransferCurveSolveTask(pipeBase.PipelineTask,
159  pipeBase.CmdLineTask):
160  """Task to fit the PTC from flat covariances.
161  This task assembles the list of individual PTC datasets produced
162  by `PhotonTransferCurveSolveTask` into one single final PTC dataset.
163  The task fits the measured (co)variances to a polynomial model or to
164  the models described in equations 16 and 20 of Astier+19
165  (referred to as `POLYNOMIAL`, `EXPAPPROXIMATION`, and `FULLCOVARIANCE`
166  in the configuration options of the task, respectively). Parameters
167  of interest such as tghe gain and noise are derived from the fits.
168 
169  Astier+19: "The Shape of the Photon Transfer Curve
170  of CCD sensors", arXiv:1905.08677
171  """
172  ConfigClass = PhotonTransferCurveSolveConfig
173  _DefaultName = 'cpPhotonTransferCurveSolve'
174 
175  def runQuantum(self, butlerQC, inputRefs, outputRefs):
176  """Ensure that the input and output dimensions are passed along.
177 
178  Parameters
179  ----------
180  butlerQC : `~lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext`
181  Butler to operate on.
182  inputRefs : `~lsst.pipe.base.connections.InputQuantizedConnection`
183  Input data refs to load.
184  ouptutRefs : `~lsst.pipe.base.connections.OutputQuantizedConnection`
185  Output data refs to persist.
186  """
187  inputs = butlerQC.get(inputRefs)
188  outputs = self.runrun(inputCovariances=inputs['inputCovariances'], camera=inputs['camera'])
189  butlerQC.put(outputs, outputRefs)
190 
191  def run(self, inputCovariances, camera=None, inputExpList=None):
192  """Fit measure covariances to different models.
193 
194  Parameters
195  ----------
196  inputCovariances : `list` [`lsst.ip.isr.PhotonTransferCurveDataset`]
197  List of lsst.ip.isr.PhotonTransferCurveDataset datasets.
198 
199  camera : `lsst.afw.cameraGeom.Camera`, optional
200  Input camera.
201 
202  inputExpList : `list` [`~lsst.afw.image.exposure.exposure.ExposureF`], optional
203  List of exposures.
204 
205  Returns
206  -------
207  results : `lsst.pipe.base.Struct`
208  The results struct containing:
209  ``outputPtcDatset`` : `lsst.ip.isr.PhotonTransferCurveDataset`
210  Final PTC dataset, containing information such as the means, variances,
211  and exposure times.
212  """
213  # Assemble partial PTC datasets into a single dataset.
214  ampNames = np.unique(inputCovariances[0].ampNames)
215  datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType,
216  self.config.maximumRangeCovariancesAstier)
217  for partialPtcDataset in inputCovariances:
218  if partialPtcDataset.ptcFitType == 'DUMMY':
219  continue
220  for ampName in ampNames:
221  datasetPtc.inputExpIdPairs[ampName].append(partialPtcDataset.inputExpIdPairs[ampName])
222  if type(partialPtcDataset.rawExpTimes[ampName]) is list:
223  datasetPtc.rawExpTimes[ampName].append(partialPtcDataset.rawExpTimes[ampName][0])
224  else:
225  datasetPtc.rawExpTimes[ampName].append(partialPtcDataset.rawExpTimes[ampName])
226  if type(partialPtcDataset.rawMeans[ampName]) is list:
227  datasetPtc.rawMeans[ampName].append(partialPtcDataset.rawMeans[ampName][0])
228  else:
229  datasetPtc.rawMeans[ampName].append(partialPtcDataset.rawMeans[ampName])
230  if type(partialPtcDataset.rawVars[ampName]) is list:
231  datasetPtc.rawVars[ampName].append(partialPtcDataset.rawVars[ampName][0])
232  else:
233  datasetPtc.rawVars[ampName].append(partialPtcDataset.rawVars[ampName])
234  datasetPtc.covariances[ampName].append(np.array(partialPtcDataset.covariances[ampName][0]))
235  datasetPtc.covariancesSqrtWeights[ampName].append(
236  np.array(partialPtcDataset.covariancesSqrtWeights[ampName][0]))
237  # Sort arrays that are filled so far in the final dataset by rawMeans index
238  for ampName in ampNames:
239  index = np.argsort(np.ravel(np.array(datasetPtc.rawMeans[ampName])))
240  datasetPtc.inputExpIdPairs[ampName] = np.array(datasetPtc.inputExpIdPairs[ampName])[index]
241  datasetPtc.rawExpTimes[ampName] = np.array(datasetPtc.rawExpTimes[ampName])[index]
242  datasetPtc.rawMeans[ampName] = np.array(datasetPtc.rawMeans[ampName])[index]
243  datasetPtc.rawVars[ampName] = np.array(datasetPtc.rawVars[ampName])[index]
244  datasetPtc.covariances[ampName] = np.array(datasetPtc.covariances[ampName])[index]
245  datasetPtc.covariancesSqrtWeights[ampName] = np.array(
246  datasetPtc.covariancesSqrtWeights[ampName])[index]
247 
248  if self.config.ptcFitType == "FULLCOVARIANCE":
249  # Calculate covariances and fit them, including the PTC, to Astier+19 full model (Eq. 20)
250  # First, fit get the flat pairs that are masked, fitting C_00 vs mu to
251  # the EXPAPPROXIMATION model (Eq. 16 in Astier+19).
252  # The points at these fluxes will also be masked when calculating the other covariances, C_ij)
253  tempDatasetPtc = copy.copy(datasetPtc)
254  tempDatasetPtc.ptcFitType = "EXPAPPROXIMATION"
255  tempDatasetPtc = self.fitPtcfitPtc(tempDatasetPtc)
256  for ampName in datasetPtc.ampNames:
257  datasetPtc.expIdMask[ampName] = tempDatasetPtc.expIdMask[ampName]
258  datasetPtc.fitType = "FULLCOVARIANCE"
259  datasetPtc = self.fitCovariancesAstierfitCovariancesAstier(datasetPtc)
260  # The other options are: self.config.ptcFitType in ("EXPAPPROXIMATION", "POLYNOMIAL")
261  else:
262  # Fit the PTC to a polynomial or to Astier+19 exponential approximation (Eq. 16).
263  # Fill up PhotonTransferCurveDataset object.
264  datasetPtc = self.fitPtcfitPtc(datasetPtc)
265  if inputExpList is not None:
266  # It should be a list of exposures, to get the detector.
267  detector = inputExpList[0].getDetector()
268  else:
269  detector = None
270  datasetPtc.updateMetadata(setDate=True, camera=camera, detector=detector)
271 
272  return pipeBase.Struct(
273  outputPtcDataset=datasetPtc,
274  )
275 
276  def fitCovariancesAstier(self, dataset):
277  """Fit measured flat covariances to full model in Astier+19.
278 
279  Parameters
280  ----------
281  dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
282  The dataset containing information such as the means, (co)variances,
283  and exposure times.
284 
285  Returns
286  -------
287  dataset: `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
288  This is the same dataset as the input paramter, however, it has been modified
289  to include information such as the fit vectors and the fit parameters. See
290  the class `PhotonTransferCurveDatase`.
291  """
292 
293  covFits, covFitsNoB = fitData(dataset,
294  r=self.config.maximumRangeCovariancesAstier)
295 
296  dataset = self.getOutputPtcDataCovAstiergetOutputPtcDataCovAstier(dataset, covFits, covFitsNoB)
297 
298  return dataset
299 
300  def getOutputPtcDataCovAstier(self, dataset, covFits, covFitsNoB):
301  """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects.
302 
303  Parameters
304  ----------
305  dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
306  The dataset containing information such as the means, variances and exposure times.
307  covFits: `dict`
308  Dictionary of CovFit objects, with amp names as keys.
309  covFitsNoB : `dict`
310  Dictionary of CovFit objects, with amp names as keys, and 'b=0' in Eq. 20 of Astier+19.
311 
312  Returns
313  -------
314  dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
315  This is the same dataset as the input paramter, however, it has been modified
316  to include extra information such as the mask 1D array, gains, reoudout noise, measured signal,
317  measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier.
318  See the class `PhotonTransferCurveDatase`.
319  """
320  assert(len(covFits) == len(covFitsNoB))
321 
322  for i, amp in enumerate(dataset.ampNames):
323  lenInputTimes = len(dataset.rawExpTimes[amp])
324  # Not used when ptcFitType is 'FULLCOVARIANCE'
325  dataset.ptcFitPars[amp] = [np.nan]
326  dataset.ptcFitParsError[amp] = [np.nan]
327  dataset.ptcFitChiSq[amp] = np.nan
328  if amp in covFits:
329  fit = covFits[amp]
330  fitNoB = covFitsNoB[amp]
331  # Save full covariances, covariances models, and their weights
332  # dataset.expIdMask is already full
333  dataset.covariances[amp] = fit.cov
334  dataset.covariancesModel[amp] = fit.evalCovModel()
335  dataset.covariancesSqrtWeights[amp] = fit.sqrtW
336  dataset.aMatrix[amp] = fit.getA()
337  dataset.bMatrix[amp] = fit.getB()
338  dataset.covariancesModelNoB[amp] = fitNoB.evalCovModel()
339  dataset.aMatrixNoB[amp] = fitNoB.getA()
340 
341  (meanVecFinal, varVecFinal, varVecModel,
342  wc, varMask) = fit.getFitData(0, 0, divideByMu=False)
343  gain = fit.getGain()
344 
345  dataset.gain[amp] = gain
346  dataset.gainErr[amp] = fit.getGainErr()
347  dataset.noise[amp] = np.sqrt(fit.getRon())
348  dataset.noiseErr[amp] = fit.getRonErr()
349  dataset.finalVars[amp] = varVecFinal
350  dataset.finalModelVars[amp] = varVecModel
351  dataset.finalMeans[amp] = meanVecFinal
352 
353  else:
354  # Bad amp
355  # Entries need to have proper dimensions so read/write with astropy.Table works.
356  matrixSide = self.config.maximumRangeCovariancesAstier
357  nanMatrix = np.full((matrixSide, matrixSide), np.nan)
358  listNanMatrix = np.full((lenInputTimes, matrixSide, matrixSide), np.nan)
359 
360  dataset.covariances[amp] = listNanMatrix
361  dataset.covariancesModel[amp] = listNanMatrix
362  dataset.covariancesSqrtWeights[amp] = listNanMatrix
363  dataset.aMatrix[amp] = nanMatrix
364  dataset.bMatrix[amp] = nanMatrix
365  dataset.covariancesModelNoB[amp] = listNanMatrix
366  dataset.aMatrixNoB[amp] = nanMatrix
367 
368  dataset.expIdMask[amp] = np.repeat(np.nan, lenInputTimes)
369  dataset.gain[amp] = np.nan
370  dataset.gainErr[amp] = np.nan
371  dataset.noise[amp] = np.nan
372  dataset.noiseErr[amp] = np.nan
373  dataset.finalVars[amp] = np.repeat(np.nan, lenInputTimes)
374  dataset.finalModelVars[amp] = np.repeat(np.nan, lenInputTimes)
375  dataset.finalMeans[amp] = np.repeat(np.nan, lenInputTimes)
376 
377  return dataset
378 
379  @staticmethod
380  def _initialParsForPolynomial(order):
381  assert(order >= 2)
382  pars = np.zeros(order, dtype=np.float)
383  pars[0] = 10
384  pars[1] = 1
385  pars[2:] = 0.0001
386  return pars
387 
388  @staticmethod
389  def _boundsForPolynomial(initialPars, lowers=[], uppers=[]):
390  if not len(lowers):
391  lowers = [np.NINF for p in initialPars]
392  if not len(uppers):
393  uppers = [np.inf for p in initialPars]
394  lowers[1] = 0 # no negative gains
395  return (lowers, uppers)
396 
397  @staticmethod
398  def _boundsForAstier(initialPars, lowers=[], uppers=[]):
399  if not len(lowers):
400  lowers = [np.NINF for p in initialPars]
401  if not len(uppers):
402  uppers = [np.inf for p in initialPars]
403  return (lowers, uppers)
404 
405  @staticmethod
406  def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative,
407  minMeanRatioTest, minVarPivotSearch):
408  """Return a boolean array to mask bad points.
409 
410  Parameters
411  ----------
412  means : `numpy.array`
413  Input array with mean signal values.
414  variances : `numpy.array`
415  Input array with variances at each mean value.
416  maxDeviationPositive : `float`
417  Maximum deviation from being constant for the variance/mean
418  ratio, in the positive direction.
419  maxDeviationNegative : `float`
420  Maximum deviation from being constant for the variance/mean
421  ratio, in the negative direction.
422  minMeanRatioTest : `float`
423  Minimum signal value (in ADU) after which to start examining
424  the ratios var/mean.
425  minVarPivotSearch : `float`
426  Minimum variance point (in ADU^2) after which the pivot point
427  wher the variance starts decreasing should be sought.
428 
429  Returns
430  ------
431  goodPoints : `numpy.array` [`bool`]
432  Boolean array to select good (`True`) and bad (`False`)
433  points.
434 
435  Notes
436  -----
437  A linear function has a constant ratio, so find the median
438  value of the ratios, and exclude the points that deviate
439  from that by more than a factor of maxDeviationPositive/negative.
440  Asymmetric deviations are supported as we expect the PTC to turn
441  down as the flux increases, but sometimes it anomalously turns
442  upwards just before turning over, which ruins the fits, so it
443  is wise to be stricter about restricting positive outliers than
444  negative ones.
445  Too high and points that are so bad that fit will fail will be included
446  Too low and the non-linear points will be excluded, biasing the NL fit.
447  This function also masks points after the variance starts decreasing.
448  """
449 
450  assert(len(means) == len(variances))
451  ratios = [b/a for (a, b) in zip(means, variances)]
452  medianRatio = np.nanmedian(ratios)
453  ratioDeviations = [0.0 if a < minMeanRatioTest else (r/medianRatio)-1
454  for (a, r) in zip(means, ratios)]
455 
456  # so that it doesn't matter if the deviation is expressed as positive or negative
457  maxDeviationPositive = abs(maxDeviationPositive)
458  maxDeviationNegative = -1. * abs(maxDeviationNegative)
459 
460  goodPoints = np.array([True if (r < maxDeviationPositive and r > maxDeviationNegative)
461  else False for r in ratioDeviations])
462 
463  # Eliminate points beyond which the variance decreases
464  pivot = np.where(np.array(np.diff(variances)) < 0)[0]
465  if len(pivot) > 0:
466  # For small values, sometimes the variance decreases slightly
467  # Only look when var > self.config.minVarPivotSearch
468  pivot = [p for p in pivot if variances[p] > minVarPivotSearch]
469  if len(pivot) > 0:
470  pivot = np.min(pivot)
471  goodPoints[pivot+1:len(goodPoints)] = False
472 
473  return goodPoints
474 
475  def _makeZeroSafe(self, array, substituteValue=1e-9):
476  """"""
477  array = np.array(array)
478  nBad = Counter(np.ravel(array))[0]
479  if nBad == 0:
480  return array
481 
482  index, = np.where(array == 0)
483  if len(index):
484  msg = f"Found {nBad} zeros in array at elements {index}"
485  self.log.warn(msg)
486 
487  array[index] = substituteValue
488 
489  return array
490 
491  def fitPtc(self, dataset):
492  """Fit the photon transfer curve to a polynomial or to Astier+19 approximation.
493 
494  Fit the photon transfer curve with either a polynomial of the order
495  specified in the task config, or using the exponential approximation
496  in Astier+19 (Eq. 16).
497 
498  Sigma clipping is performed iteratively for the fit, as well as an
499  initial clipping of data points that are more than
500  config.initialNonLinearityExclusionThreshold away from lying on a
501  straight line. This other step is necessary because the photon transfer
502  curve turns over catastrophically at very high flux (because saturation
503  drops the variance to ~0) and these far outliers cause the initial fit
504  to fail, meaning the sigma cannot be calculated to perform the
505  sigma-clipping.
506 
507  Parameters
508  ----------
509  dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
510  The dataset containing the means, variances and exposure times.
511 
512  Returns
513  -------
514  dataset: `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
515  This is the same dataset as the input parameter, however, it has been modified
516  to include information such as the fit vectors and the fit parameters. See
517  the class `PhotonTransferCurveDatase`.
518 
519  Raises
520  ------
521  RuntimeError:
522  Raises if dataset.ptcFitType is None or empty.
523  """
524  if dataset.ptcFitType:
525  ptcFitType = dataset.ptcFitType
526  else:
527  raise RuntimeError("ptcFitType is None of empty in PTC dataset.")
528  matrixSide = self.config.maximumRangeCovariancesAstier
529  nanMatrix = np.empty((matrixSide, matrixSide))
530  nanMatrix[:] = np.nan
531 
532  for amp in dataset.ampNames:
533  lenInputTimes = len(dataset.rawExpTimes[amp])
534  listNanMatrix = np.empty((lenInputTimes, matrixSide, matrixSide))
535  listNanMatrix[:] = np.nan
536 
537  dataset.covariancesModel[amp] = listNanMatrix
538  dataset.aMatrix[amp] = nanMatrix
539  dataset.bMatrix[amp] = nanMatrix
540  dataset.covariancesModelNoB[amp] = listNanMatrix
541  dataset.aMatrixNoB[amp] = nanMatrix
542 
543  def errFunc(p, x, y):
544  return ptcFunc(p, x) - y
545 
546  sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers
547  maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers
548 
549  for i, ampName in enumerate(dataset.ampNames):
550  timeVecOriginal = np.ravel(np.array(dataset.rawExpTimes[ampName]))
551  meanVecOriginal = np.ravel(np.array(dataset.rawMeans[ampName]))
552  varVecOriginal = np.ravel(np.array(dataset.rawVars[ampName]))
553  varVecOriginal = self._makeZeroSafe_makeZeroSafe(varVecOriginal)
554 
555  goodPoints = self._getInitialGoodPoints_getInitialGoodPoints(meanVecOriginal, varVecOriginal,
556  self.config.initialNonLinearityExclusionThresholdPositive,
557  self.config.initialNonLinearityExclusionThresholdNegative,
558  self.config.minMeanRatioTest,
559  self.config.minVarPivotSearch)
560  if not (goodPoints.any()):
561  msg = (f"SERIOUS: All points in goodPoints: {goodPoints} are bad."
562  f"Setting {ampName} to BAD.")
563  self.log.warn(msg)
564  # Fill entries with NaNs
565  self.fillBadAmpfillBadAmp(dataset, ptcFitType, ampName)
566  continue
567 
568  mask = goodPoints
569 
570  if ptcFitType == 'EXPAPPROXIMATION':
571  ptcFunc = funcAstier
572  parsIniPtc = [-1e-9, 1.0, 10.] # a00, gain, noisei^2
573  # lowers and uppers obtained from BOT data studies by C. Lage (UC Davis, 11/2020).
574  bounds = self._boundsForAstier_boundsForAstier(parsIniPtc, lowers=[-1e-4, 0.5, -2000],
575  uppers=[1e-4, 2.5, 2000])
576  if ptcFitType == 'POLYNOMIAL':
577  ptcFunc = funcPolynomial
578  parsIniPtc = self._initialParsForPolynomial_initialParsForPolynomial(self.config.polynomialFitDegree + 1)
579  bounds = self._boundsForPolynomial_boundsForPolynomial(parsIniPtc)
580 
581  # Before bootstrap fit, do an iterative fit to get rid of outliers
582  count = 1
583  while count <= maxIterationsPtcOutliers:
584  # Note that application of the mask actually shrinks the array
585  # to size rather than setting elements to zero (as we want) so
586  # always update mask itself and re-apply to the original data
587  meanTempVec = meanVecOriginal[mask]
588  varTempVec = varVecOriginal[mask]
589  res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec))
590  pars = res.x
591 
592  # change this to the original from the temp because the masks are ANDed
593  # meaning once a point is masked it's always masked, and the masks must
594  # always be the same length for broadcasting
595  sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal)
596  newMask = np.array([True if np.abs(r) < sigmaCutPtcOutliers else False for r in sigResids])
597  mask = mask & newMask
598  if not (mask.any() and newMask.any()):
599  msg = (f"SERIOUS: All points in either mask: {mask} or newMask: {newMask} are bad. "
600  f"Setting {ampName} to BAD.")
601  self.log.warn(msg)
602  # Fill entries with NaNs
603  self.fillBadAmpfillBadAmp(dataset, ptcFitType, ampName)
604  break
605  nDroppedTotal = Counter(mask)[False]
606  self.log.debug(f"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}")
607  count += 1
608  # objects should never shrink
609  assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal))
610  if not (mask.any() and newMask.any()):
611  continue
612  dataset.expIdMask[ampName] = mask # store the final mask
613  parsIniPtc = pars
614  meanVecFinal = meanVecOriginal[mask]
615  varVecFinal = varVecOriginal[mask]
616 
617  if Counter(mask)[False] > 0:
618  self.log.info((f"Number of points discarded in PTC of amplifier {ampName}:"
619  f" {Counter(mask)[False]} out of {len(meanVecOriginal)}"))
620 
621  if (len(meanVecFinal) < len(parsIniPtc)):
622  msg = (f"SERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of "
623  f"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.")
624  self.log.warn(msg)
625  # Fill entries with NaNs
626  self.fillBadAmpfillBadAmp(dataset, ptcFitType, ampName)
627  continue
628  # Fit the PTC
629  if self.config.doFitBootstrap:
630  parsFit, parsFitErr, reducedChiSqPtc = fitBootstrap(parsIniPtc, meanVecFinal,
631  varVecFinal, ptcFunc,
632  weightsY=1./np.sqrt(varVecFinal))
633  else:
634  parsFit, parsFitErr, reducedChiSqPtc = fitLeastSq(parsIniPtc, meanVecFinal,
635  varVecFinal, ptcFunc,
636  weightsY=1./np.sqrt(varVecFinal))
637  dataset.ptcFitPars[ampName] = parsFit
638  dataset.ptcFitParsError[ampName] = parsFitErr
639  dataset.ptcFitChiSq[ampName] = reducedChiSqPtc
640  # Masked variances (measured and modeled) and means. Need to pad the array so astropy.Table does
641  # not crash (the mask may vary per amp).
642  padLength = len(dataset.rawExpTimes[ampName]) - len(varVecFinal)
643  dataset.finalVars[ampName] = np.pad(varVecFinal, (0, padLength), 'constant',
644  constant_values=np.nan)
645  dataset.finalModelVars[ampName] = np.pad(ptcFunc(parsFit, meanVecFinal), (0, padLength),
646  'constant', constant_values=np.nan)
647  dataset.finalMeans[ampName] = np.pad(meanVecFinal, (0, padLength), 'constant',
648  constant_values=np.nan)
649  if ptcFitType == 'EXPAPPROXIMATION':
650  ptcGain = parsFit[1]
651  ptcGainErr = parsFitErr[1]
652  ptcNoise = np.sqrt(np.fabs(parsFit[2]))
653  ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2]))
654  if ptcFitType == 'POLYNOMIAL':
655  ptcGain = 1./parsFit[1]
656  ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1])
657  ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain
658  ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain
659  dataset.gain[ampName] = ptcGain
660  dataset.gainErr[ampName] = ptcGainErr
661  dataset.noise[ampName] = ptcNoise
662  dataset.noiseErr[ampName] = ptcNoiseErr
663 
664  if not len(dataset.ptcFitType) == 0:
665  dataset.ptcFitType = ptcFitType
666  if len(dataset.badAmps) == 0:
667  dataset.badAmps = np.repeat(np.nan, len(list(dataset.rawExpTimes.values())[0]))
668 
669  return dataset
670 
671  def fillBadAmp(self, dataset, ptcFitType, ampName):
672  """Fill the dataset with NaNs if there are not enough good points.
673 
674  Parameters
675  ----------
676  dataset : `lsst.ip.isr.ptcDataset.PhotonTransferCurveDataset`
677  The dataset containing the means, variances and exposure times.
678  ptcFitType : `str`
679  Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or
680  'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC.
681  ampName : `str`
682  Amplifier name.
683  """
684  dataset.badAmps.append(ampName)
685  dataset.expIdMask[ampName] = np.repeat(False, len(dataset.rawExpTimes[ampName]))
686  dataset.gain[ampName] = np.nan
687  dataset.gainErr[ampName] = np.nan
688  dataset.noise[ampName] = np.nan
689  dataset.noiseErr[ampName] = np.nan
690  dataset.ptcFitPars[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
691  ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
692  dataset.ptcFitParsError[ampName] = (np.repeat(np.nan, self.config.polynomialFitDegree + 1) if
693  ptcFitType in ["POLYNOMIAL", ] else np.repeat(np.nan, 3))
694  dataset.ptcFitChiSq[ampName] = np.nan
695  dataset.finalVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
696  dataset.finalModelVars[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
697  dataset.finalMeans[ampName] = np.repeat(np.nan, len(dataset.rawExpTimes[ampName]))
698 
699  return
def _boundsForPolynomial(initialPars, lowers=[], uppers=[])
def runQuantum(self, butlerQC, inputRefs, outputRefs)
def run(self, inputCovariances, camera=None, inputExpList=None)
def _makeZeroSafe(self, array, substituteValue=1e-9)
def _boundsForAstier(initialPars, lowers=[], uppers=[])
def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative, minMeanRatioTest, minVarPivotSearch)
def getOutputPtcDataCovAstier(self, dataset, covFits, covFitsNoB)
def fillBadAmp(self, dataset, ptcFitType, ampName)
def fitBootstrap(initialParams, dataX, dataY, function, weightsY=None, confidenceSigma=1.)
Definition: utils.py:411
def fitLeastSq(initialParams, dataX, dataY, function, weightsY=None)
Definition: utils.py:344