lsst.cp.pipe  20.0.0-5-g5e0f401
ptc.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 
23 __all__ = ['MeasurePhotonTransferCurveTask',
24  'MeasurePhotonTransferCurveTaskConfig',
25  'PhotonTransferCurveDataset']
26 
27 import numpy as np
28 import matplotlib.pyplot as plt
29 from sqlite3 import OperationalError
30 from collections import Counter
31 from dataclasses import dataclass
32 
33 import lsst.afw.math as afwMath
34 import lsst.pex.config as pexConfig
35 import lsst.pipe.base as pipeBase
36 from .utils import (NonexistentDatasetTaskDataIdContainer, PairedVisitListTaskRunner,
37  checkExpLengthEqual, fitLeastSq, fitBootstrap, funcPolynomial, funcAstier)
38 from scipy.optimize import least_squares
39 
40 from lsst.ip.isr.linearize import Linearizer
41 import datetime
42 
43 from .astierCovPtcUtils import (fftSize, CovFft, computeCovDirect, fitData)
44 
45 
46 class MeasurePhotonTransferCurveTaskConfig(pexConfig.Config):
47  """Config class for photon transfer curve measurement task"""
48  ccdKey = pexConfig.Field(
49  dtype=str,
50  doc="The key by which to pull a detector from a dataId, e.g. 'ccd' or 'detector'.",
51  default='ccd',
52  )
53  ptcFitType = pexConfig.ChoiceField(
54  dtype=str,
55  doc="Fit PTC to approximation in Astier+19 (Equation 16) or to a polynomial.",
56  default="POLYNOMIAL",
57  allowed={
58  "POLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegree' to set 'n').",
59  "EXPAPPROXIMATION": "Approximation in Astier+19 (Eq. 16).",
60  "FULLCOVARIANCE": "Full covariances model in Astier+19 (Eq. 20)"
61  }
62  )
63  maximumRangeCovariancesAstier = pexConfig.Field(
64  dtype=int,
65  doc="Maximum range of covariances as in Astier+19",
66  default=8,
67  )
68  covAstierRealSpace = pexConfig.Field(
69  dtype=bool,
70  doc="Calculate covariances in real space or via FFT? (see appendix A of Astier+19).",
71  default=False,
72  )
73  polynomialFitDegree = pexConfig.Field(
74  dtype=int,
75  doc="Degree of polynomial to fit the PTC, when 'ptcFitType'=POLYNOMIAL.",
76  default=3,
77  )
78  doCreateLinearizer = pexConfig.Field(
79  dtype=bool,
80  doc="Calculate non-linearity and persist linearizer?",
81  default=False,
82  )
83  linearizerType = pexConfig.ChoiceField(
84  dtype=str,
85  doc="Linearizer type, if doCreateLinearizer=True",
86  default="LINEARIZEPOLYNOMIAL",
87  allowed={
88  "LINEARIZEPOLYNOMIAL": "n-degree polynomial (use 'polynomialFitDegreeNonLinearity' to set 'n').",
89  "LINEARIZESQUARED": "c0 quadratic coefficient derived from coefficients of polynomiual fit",
90  "LOOKUPTABLE": "Loouk table formed from linear part of polynomial fit."
91  }
92  )
93  polynomialFitDegreeNonLinearity = pexConfig.Field(
94  dtype=int,
95  doc="If doCreateLinearizer, degree of polynomial to fit the meanSignal vs exposureTime" +
96  " curve to produce the table for LinearizeLookupTable.",
97  default=3,
98  )
99  binSize = pexConfig.Field(
100  dtype=int,
101  doc="Bin the image by this factor in both dimensions.",
102  default=1,
103  )
104  minMeanSignal = pexConfig.Field(
105  dtype=float,
106  doc="Minimum value (inclusive) of mean signal (in DN) above which to consider.",
107  default=0,
108  )
109  maxMeanSignal = pexConfig.Field(
110  dtype=float,
111  doc="Maximum value (inclusive) of mean signal (in DN) below which to consider.",
112  default=9e6,
113  )
114  initialNonLinearityExclusionThresholdPositive = pexConfig.RangeField(
115  dtype=float,
116  doc="Initially exclude data points with a variance that are more than a factor of this from being"
117  " linear in the positive direction, from the PTC fit. Note that these points will also be"
118  " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
119  " to allow an accurate determination of the sigmas for said iterative fit.",
120  default=0.12,
121  min=0.0,
122  max=1.0,
123  )
124  initialNonLinearityExclusionThresholdNegative = pexConfig.RangeField(
125  dtype=float,
126  doc="Initially exclude data points with a variance that are more than a factor of this from being"
127  " linear in the negative direction, from the PTC fit. Note that these points will also be"
128  " excluded from the non-linearity fit. This is done before the iterative outlier rejection,"
129  " to allow an accurate determination of the sigmas for said iterative fit.",
130  default=0.25,
131  min=0.0,
132  max=1.0,
133  )
134  sigmaCutPtcOutliers = pexConfig.Field(
135  dtype=float,
136  doc="Sigma cut for outlier rejection in PTC.",
137  default=5.0,
138  )
139  maskNameList = pexConfig.ListField(
140  dtype=str,
141  doc="Mask list to exclude from statistics calculations.",
142  default=['SUSPECT', 'BAD', 'NO_DATA'],
143  )
144  nSigmaClipPtc = pexConfig.Field(
145  dtype=float,
146  doc="Sigma cut for afwMath.StatisticsControl()",
147  default=5.5,
148  )
149  nIterSigmaClipPtc = pexConfig.Field(
150  dtype=int,
151  doc="Number of sigma-clipping iterations for afwMath.StatisticsControl()",
152  default=1,
153  )
154  maxIterationsPtcOutliers = pexConfig.Field(
155  dtype=int,
156  doc="Maximum number of iterations for outlier rejection in PTC.",
157  default=2,
158  )
159  doFitBootstrap = pexConfig.Field(
160  dtype=bool,
161  doc="Use bootstrap for the PTC fit parameters and errors?.",
162  default=False,
163  )
164  maxAduForLookupTableLinearizer = pexConfig.Field(
165  dtype=int,
166  doc="Maximum DN value for the LookupTable linearizer.",
167  default=2**18,
168  )
169  instrumentName = pexConfig.Field(
170  dtype=str,
171  doc="Instrument name.",
172  default='',
173  )
174 
175 
176 @dataclass
178  """A simple class to hold the output from the
179  `calculateLinearityResidualAndLinearizers` function.
180  """
181  # Normalized coefficients for polynomial NL correction
182  polynomialLinearizerCoefficients: list
183  # Normalized coefficient for quadratic polynomial NL correction (c0)
184  quadraticPolynomialLinearizerCoefficient: float
185  # LUT array row for the amplifier at hand
186  linearizerTableRow: list
187  meanSignalVsTimePolyFitPars: list
188  meanSignalVsTimePolyFitParsErr: list
189  meanSignalVsTimePolyFitReducedChiSq: float
190 
191 
193  """A simple class to hold the output data from the PTC task.
194 
195  The dataset is made up of a dictionary for each item, keyed by the
196  amplifiers' names, which much be supplied at construction time.
197 
198  New items cannot be added to the class to save accidentally saving to the
199  wrong property, and the class can be frozen if desired.
200 
201  inputVisitPairs records the visits used to produce the data.
202  When fitPtc() or fitCovariancesAstier() is run, a mask is built up, which is by definition
203  always the same length as inputVisitPairs, rawExpTimes, rawMeans
204  and rawVars, and is a list of bools, which are incrementally set to False
205  as points are discarded from the fits.
206 
207  PTC fit parameters for polynomials are stored in a list in ascending order
208  of polynomial term, i.e. par[0]*x^0 + par[1]*x + par[2]*x^2 etc
209  with the length of the list corresponding to the order of the polynomial
210  plus one.
211 
212  Parameters
213  ----------
214  ampNames : `list`
215  List with the names of the amplifiers of the detector at hand.
216 
217  ptcFitType : `str`
218  Type of model fitted to the PTC: "POLYNOMIAL", "EXPAPPROXIMATION", or "FULLCOVARIANCE".
219 
220  Returns
221  -------
222  `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
223  Output dataset from MeasurePhotonTransferCurveTask.
224  """
225 
226  def __init__(self, ampNames, ptcFitType):
227  # add items to __dict__ directly because __setattr__ is overridden
228 
229  # instance variables
230  self.__dict__["ptcFitType"] = ptcFitType
231  self.__dict__["ampNames"] = ampNames
232  self.__dict__["badAmps"] = []
233 
234  # raw data variables
235  # visitMask is the mask produced after outlier rejection. The mask produced by "FULLCOVARIANCE"
236  # may differ from the one produced in the other two PTC fit types.
237  self.__dict__["inputVisitPairs"] = {ampName: [] for ampName in ampNames}
238  self.__dict__["visitMask"] = {ampName: [] for ampName in ampNames}
239  self.__dict__["rawExpTimes"] = {ampName: [] for ampName in ampNames}
240  self.__dict__["rawMeans"] = {ampName: [] for ampName in ampNames}
241  self.__dict__["rawVars"] = {ampName: [] for ampName in ampNames}
242 
243  # Gain and noise
244  self.__dict__["gain"] = {ampName: -1. for ampName in ampNames}
245  self.__dict__["gainErr"] = {ampName: -1. for ampName in ampNames}
246  self.__dict__["noise"] = {ampName: -1. for ampName in ampNames}
247  self.__dict__["noiseErr"] = {ampName: -1. for ampName in ampNames}
248 
249  # if ptcFitTye in ["POLYNOMIAL", "EXPAPPROXIMATION"]
250  # fit information
251  self.__dict__["ptcFitPars"] = {ampName: [] for ampName in ampNames}
252  self.__dict__["ptcFitParsError"] = {ampName: [] for ampName in ampNames}
253  self.__dict__["ptcFitReducedChiSquared"] = {ampName: [] for ampName in ampNames}
254 
255  # if ptcFitTye in ["FULLCOVARIANCE"]
256  # "covariancesTuple" is a numpy recarray with entries of the form
257  # ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName']
258  # "covariancesFits" has CovFit objects that fit the measured covariances to Eq. 20 of Astier+19.
259  # In "covariancesFitsWithNoB", "b"=0 in the model described by Eq. 20 of Astier+19.
260  self.__dict__["covariancesTuple"] = {ampName: [] for ampName in ampNames}
261  self.__dict__["covariancesFitsWithNoB"] = {ampName: [] for ampName in ampNames}
262  self.__dict__["covariancesFits"] = {ampName: [] for ampName in ampNames}
263  self.__dict__["aMatrix"] = {ampName: [] for ampName in ampNames}
264  self.__dict__["bMatrix"] = {ampName: [] for ampName in ampNames}
265 
266  # "final" means that the "raw" vectors above had "visitMask" applied.
267  self.__dict__["finalVars"] = {ampName: [] for ampName in ampNames}
268  self.__dict__["finalModelVars"] = {ampName: [] for ampName in ampNames}
269  self.__dict__["finalMeans"] = {ampName: [] for ampName in ampNames}
270 
271  def __setattr__(self, attribute, value):
272  """Protect class attributes"""
273  if attribute not in self.__dict__:
274  raise AttributeError(f"{attribute} is not already a member of PhotonTransferCurveDataset, which"
275  " does not support setting of new attributes.")
276  else:
277  self.__dict__[attribute] = value
278 
279  def getVisitsUsed(self, ampName):
280  """Get the visits used, i.e. not discarded, for a given amp.
281 
282  If no mask has been created yet, all visits are returned.
283  """
284  if len(self.visitMask[ampName]) == 0:
285  return self.inputVisitPairs[ampName]
286 
287  # if the mask exists it had better be the same length as the visitPairs
288  assert len(self.visitMask[ampName]) == len(self.inputVisitPairs[ampName])
289 
290  pairs = self.inputVisitPairs[ampName]
291  mask = self.visitMask[ampName]
292  # cast to bool required because numpy
293  return [(v1, v2) for ((v1, v2), m) in zip(pairs, mask) if bool(m) is True]
294 
295  def getGoodAmps(self):
296  return [amp for amp in self.ampNames if amp not in self.badAmps]
297 
298 
299 class MeasurePhotonTransferCurveTask(pipeBase.CmdLineTask):
300  """A class to calculate, fit, and plot a PTC from a set of flat pairs.
301 
302  The Photon Transfer Curve (var(signal) vs mean(signal)) is a standard tool
303  used in astronomical detectors characterization (e.g., Janesick 2001,
304  Janesick 2007). If ptcFitType is "EXPAPPROXIMATION" or "POLYNOMIAL", this task calculates the
305  PTC from a series of pairs of flat-field images; each pair taken at identical exposure
306  times. The difference image of each pair is formed to eliminate fixed pattern noise,
307  and then the variance of the difference image and the mean of the average image
308  are used to produce the PTC. An n-degree polynomial or the approximation in Equation
309  16 of Astier+19 ("The Shape of the Photon Transfer Curve of CCD sensors",
310  arXiv:1905.08677) can be fitted to the PTC curve. These models include
311  parameters such as the gain (e/DN) and readout noise.
312 
313  Linearizers to correct for signal-chain non-linearity are also calculated.
314  The `Linearizer` class, in general, can support per-amp linearizers, but in this
315  task this is not supported.
316 
317  If ptcFitType is "FULLCOVARIANCE", the covariances of the difference images are calculated via the
318  DFT methods described in Astier+19 and the variances for the PTC are given by the cov[0,0] elements
319  at each signal level. The full model in Equation 20 of Astier+19 is fit to the PTC to get the gain
320  and the noise.
321 
322  Parameters
323  ----------
324 
325  *args: `list`
326  Positional arguments passed to the Task constructor. None used at this
327  time.
328  **kwargs: `dict`
329  Keyword arguments passed on to the Task constructor. None used at this
330  time.
331 
332  """
333 
334  RunnerClass = PairedVisitListTaskRunner
335  ConfigClass = MeasurePhotonTransferCurveTaskConfig
336  _DefaultName = "measurePhotonTransferCurve"
337 
338  def __init__(self, *args, **kwargs):
339  pipeBase.CmdLineTask.__init__(self, *args, **kwargs)
340  plt.interactive(False) # stop windows popping up when plotting. When headless, use 'agg' backend too
341  self.config.validate()
342  self.config.freeze()
343 
344  @classmethod
345  def _makeArgumentParser(cls):
346  """Augment argument parser for the MeasurePhotonTransferCurveTask."""
347  parser = pipeBase.ArgumentParser(name=cls._DefaultName)
348  parser.add_argument("--visit-pairs", dest="visitPairs", nargs="*",
349  help="Visit pairs to use. Each pair must be of the form INT,INT e.g. 123,456")
350  parser.add_id_argument("--id", datasetType="photonTransferCurveDataset",
351  ContainerClass=NonexistentDatasetTaskDataIdContainer,
352  help="The ccds to use, e.g. --id ccd=0..100")
353  return parser
354 
355  @pipeBase.timeMethod
356  def runDataRef(self, dataRef, visitPairs):
357  """Run the Photon Transfer Curve (PTC) measurement task.
358 
359  For a dataRef (which is each detector here),
360  and given a list of visit pairs (postISR) at different exposure times,
361  measure the PTC.
362 
363  Parameters
364  ----------
365  dataRef : list of lsst.daf.persistence.ButlerDataRef
366  dataRef for the detector for the visits to be fit.
367 
368  visitPairs : `iterable` of `tuple` of `int`
369  Pairs of visit numbers to be processed together
370  """
371 
372  # setup necessary objects
373  detNum = dataRef.dataId[self.config.ccdKey]
374  detector = dataRef.get('camera')[dataRef.dataId[self.config.ccdKey]]
375  # expand some missing fields that we need for lsstCam. This is a work-around
376  # for Gen2 problems that I (RHL) don't feel like solving. The calibs pipelines
377  # (which inherit from CalibTask) use addMissingKeys() to do basically the same thing
378  #
379  # Basically, the butler's trying to look up the fields in `raw_visit` which won't work
380  for name in dataRef.getButler().getKeys('bias'):
381  if name not in dataRef.dataId:
382  try:
383  dataRef.dataId[name] = \
384  dataRef.getButler().queryMetadata('raw', [name], detector=detNum)[0]
385  except OperationalError:
386  pass
387 
388  amps = detector.getAmplifiers()
389  ampNames = [amp.getName() for amp in amps]
390  datasetPtc = PhotonTransferCurveDataset(ampNames, self.config.ptcFitType)
391  self.log.info('Measuring PTC using %s visits for detector %s' % (visitPairs, detector.getId()))
392 
393  tupleRecords = []
394  allTags = []
395  for (v1, v2) in visitPairs:
396  # Get postISR exposures.
397  dataRef.dataId['expId'] = v1
398  exp1 = dataRef.get("postISRCCD", immediate=True)
399  dataRef.dataId['expId'] = v2
400  exp2 = dataRef.get("postISRCCD", immediate=True)
401  del dataRef.dataId['expId']
402 
403  checkExpLengthEqual(exp1, exp2, v1, v2, raiseWithMessage=True)
404  expTime = exp1.getInfo().getVisitInfo().getExposureTime()
405 
406  tupleRows = []
407  nAmpsNan = 0
408  for ampNumber, amp in enumerate(detector):
409  ampName = amp.getName()
410  # covAstier: (i, j, var (cov[0,0]), cov, npix)
411  doRealSpace = self.config.covAstierRealSpace
412  muDiff, varDiff, covAstier = self.measureMeanVarCov(exp1, exp2, region=amp.getBBox(),
413  covAstierRealSpace=doRealSpace)
414  if np.isnan(muDiff) or np.isnan(varDiff) or np.isnan(covAstier):
415  msg = f"NaN mean or cov in amp {ampNumber} in visit pair {v1}, {v2} of detector {detNum}."
416  self.log.warn(msg)
417  nAmpsNan += 1
418  continue
419 
420  datasetPtc.rawExpTimes[ampName].append(expTime)
421  datasetPtc.rawMeans[ampName].append(muDiff)
422  datasetPtc.rawVars[ampName].append(varDiff)
423  datasetPtc.inputVisitPairs[ampName].append((v1, v2))
424 
425  tupleRows += [(muDiff, ) + covRow + (ampNumber, expTime, ampName) for covRow in covAstier]
426  tags = ['mu', 'i', 'j', 'var', 'cov', 'npix', 'ext', 'expTime', 'ampName']
427  if nAmpsNan == len(ampNames):
428  msg = f"NaN mean in all amps of visit pair {v1}, {v2} of detector {detNum}."
429  self.log.warn(msg)
430  continue
431  allTags += tags
432  tupleRecords += tupleRows
433  covariancesWithTags = np.core.records.fromrecords(tupleRecords, names=allTags)
434 
435  if self.config.ptcFitType in ["FULLCOVARIANCE", ]:
436  # Calculate covariances and fit them, including the PTC, to Astier+19 full model (Eq. 20)
437  datasetPtc = self.fitCovariancesAstier(datasetPtc, covariancesWithTags)
438  elif self.config.ptcFitType in ["EXPAPPROXIMATION", "POLYNOMIAL"]:
439  # Fit the PTC to a polynomial or to Astier+19 exponential approximation (Eq. 16)
440  # Fill up PhotonTransferCurveDataset object.
441  datasetPtc = self.fitPtc(datasetPtc, self.config.ptcFitType)
442 
443  # Fit a poynomial to calculate non-linearity and persist linearizer.
444  if self.config.doCreateLinearizer:
445  numberAmps = len(amps)
446  numberAduValues = self.config.maxAduForLookupTableLinearizer
447  lookupTableArray = np.zeros((numberAmps, numberAduValues), dtype=np.float32)
448 
449  # Fit (non)linearity of signal vs time curve.
450  # Fill up PhotonTransferCurveDataset object.
451  # Fill up array for LUT linearizer (tableArray).
452  # Produce coefficients for Polynomial ans Squared linearizers.
453  # Build linearizer objects.
454  linearizer = self.fitNonLinearityAndBuildLinearizers(datasetPtc, detector,
455  tableArray=lookupTableArray,
456  log=self.log)
457 
458  if self.config.linearizerType == "LINEARIZEPOLYNOMIAL":
459  linDataType = 'linearizePolynomial'
460  linMsg = "polynomial (coefficients for a polynomial correction)."
461  elif self.config.linearizerType == "LINEARIZESQUARED":
462  linDataType = 'linearizePolynomial'
463  linMsg = "squared (c0, derived from k_i coefficients of a polynomial fit)."
464  elif self.config.linearizerType == "LOOKUPTABLE":
465  linDataType = 'linearizePolynomial'
466  linMsg = "lookup table (linear component of polynomial fit)."
467  else:
468  raise RuntimeError("Invalid config.linearizerType {selg.config.linearizerType}. "
469  "Supported: 'LOOKUPTABLE', 'LINEARIZESQUARED', or 'LINEARIZEPOLYNOMIAL'")
470 
471  butler = dataRef.getButler()
472  self.log.info(f"Writing linearizer: \n {linMsg}")
473 
474  detName = detector.getName()
475  now = datetime.datetime.utcnow()
476  calibDate = now.strftime("%Y-%m-%d")
477 
478  butler.put(linearizer, datasetType=linDataType, dataId={'detector': detNum,
479  'detectorName': detName, 'calibDate': calibDate})
480 
481  self.log.info(f"Writing PTC data to {dataRef.getUri(write=True)}")
482  dataRef.put(datasetPtc, datasetType="photonTransferCurveDataset")
483 
484  return pipeBase.Struct(exitStatus=0)
485 
486  def fitCovariancesAstier(self, dataset, covariancesWithTagsArray):
487  """Fit measured flat covariances to full model in Astier+19.
488 
489  Parameters
490  ----------
491  dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
492  The dataset containing information such as the means, variances and exposure times.
493 
494  covariancesWithTagsArray : `numpy.recarray`
495  Tuple with at least (mu, cov, var, i, j, npix), where:
496  mu : 0.5*(m1 + m2), where:
497  mu1: mean value of flat1
498  mu2: mean value of flat2
499  cov: covariance value at lag(i, j)
500  var: variance(covariance value at lag(0, 0))
501  i: lag dimension
502  j: lag dimension
503  npix: number of pixels used for covariance calculation.
504 
505  Returns
506  -------
507  dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
508  This is the same dataset as the input paramter, however, it has been modified
509  to include information such as the fit vectors and the fit parameters. See
510  the class `PhotonTransferCurveDatase`.
511  """
512 
513  covFits, covFitsNoB = fitData(covariancesWithTagsArray, maxMu=self.config.maxMeanSignal,
514  r=self.config.maximumRangeCovariancesAstier)
515 
516  dataset.covariancesTuple = covariancesWithTagsArray
517  dataset.covariancesFits = covFits
518  dataset.covariancesFitsWithNoB = covFitsNoB
519  dataset = self.getOutputPtcDataCovAstier(dataset, covFits)
520 
521  return dataset
522 
523  def getOutputPtcDataCovAstier(self, dataset, covFits):
524  """Get output data for PhotonTransferCurveCovAstierDataset from CovFit objects.
525 
526  Parameters
527  ----------
528  dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
529  The dataset containing information such as the means, variances and exposure times.
530 
531  covFits: `dict`
532  Dictionary of CovFit objects, with amp names as keys.
533 
534  Returns
535  -------
536  dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
537  This is the same dataset as the input paramter, however, it has been modified
538  to include extra information such as the mask 1D array, gains, reoudout noise, measured signal,
539  measured variance, modeled variance, a, and b coefficient matrices (see Astier+19) per amplifier.
540  See the class `PhotonTransferCurveDatase`.
541  """
542 
543  for i, amp in enumerate(covFits):
544  fit = covFits[amp]
545  meanVecFinal, varVecFinal, varVecModel, wc = fit.getNormalizedFitData(0, 0, divideByMu=False)
546  gain = fit.getGain()
547  dataset.visitMask[amp] = fit.getMaskVar()
548  dataset.gain[amp] = gain
549  dataset.gainErr[amp] = fit.getGainErr()
550  dataset.noise[amp] = np.sqrt(np.fabs(fit.getRon()))
551  dataset.noiseErr[amp] = fit.getRonErr()
552  dataset.finalVars[amp].append(varVecFinal/(gain**2))
553  dataset.finalModelVars[amp].append(varVecModel/(gain**2))
554  dataset.finalMeans[amp].append(meanVecFinal/gain)
555  dataset.aMatrix[amp].append(fit.getA())
556  dataset.bMatrix[amp].append(fit.getB())
557 
558  return dataset
559 
560  def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False):
561  """Calculate the mean of each of two exposures and the variance and covariance of their difference.
562 
563  The variance is calculated via afwMath, and the covariance via the methods in Astier+19 (appendix A).
564  In theory, var = covariance[0,0]. This should be validated, and in the future, we may decide to just
565  keep one (covariance).
566 
567  Parameters
568  ----------
569  exposure1 : `lsst.afw.image.exposure.exposure.ExposureF`
570  First exposure of flat field pair.
571 
572  exposure2 : `lsst.afw.image.exposure.exposure.ExposureF`
573  Second exposure of flat field pair.
574 
575  region : `lsst.geom.Box2I`, optional
576  Region of each exposure where to perform the calculations (e.g, an amplifier).
577 
578  covAstierRealSpace : `bool`, optional
579  Should the covariannces in Astier+19 be calculated in real space or via FFT?
580  See Appendix A of Astier+19.
581 
582  Returns
583  -------
584  mu : `float` or `NaN`
585  0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in
586  both exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
587 
588  varDiff : `float` or `NaN`
589  Half of the clipped variance of the difference of the regions inthe two input
590  exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
591 
592  covDiffAstier : `list` or `NaN`
593  List with tuples of the form (dx, dy, var, cov, npix), where:
594  dx : `int`
595  Lag in x
596  dy : `int`
597  Lag in y
598  var : `float`
599  Variance at (dx, dy).
600  cov : `float`
601  Covariance at (dx, dy).
602  nPix : `int`
603  Number of pixel pairs used to evaluate var and cov.
604  If either mu1 or m2 are NaN's, the returned value is NaN.
605  """
606 
607  if region is not None:
608  im1Area = exposure1.maskedImage[region]
609  im2Area = exposure2.maskedImage[region]
610  else:
611  im1Area = exposure1.maskedImage
612  im2Area = exposure2.maskedImage
613 
614  im1Area = afwMath.binImage(im1Area, self.config.binSize)
615  im2Area = afwMath.binImage(im2Area, self.config.binSize)
616 
617  im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList)
618  im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
619  self.config.nIterSigmaClipPtc,
620  im1MaskVal)
621  im1StatsCtrl.setNanSafe(True)
622  im1StatsCtrl.setAndMask(im1MaskVal)
623 
624  im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList)
625  im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
626  self.config.nIterSigmaClipPtc,
627  im2MaskVal)
628  im2StatsCtrl.setNanSafe(True)
629  im2StatsCtrl.setAndMask(im2MaskVal)
630 
631  # Clipped mean of images; then average of mean.
632  mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue()
633  mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue()
634  if np.isnan(mu1) or np.isnan(mu2):
635  return np.nan, np.nan, np.nan
636  mu = 0.5*(mu1 + mu2)
637 
638  # Take difference of pairs
639  # symmetric formula: diff = (mu2*im1-mu1*im2)/(0.5*(mu1+mu2))
640  temp = im2Area.clone()
641  temp *= mu1
642  diffIm = im1Area.clone()
643  diffIm *= mu2
644  diffIm -= temp
645  diffIm /= mu
646 
647  diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList)
648  diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
649  self.config.nIterSigmaClipPtc,
650  diffImMaskVal)
651  diffImStatsCtrl.setNanSafe(True)
652  diffImStatsCtrl.setAndMask(diffImMaskVal)
653 
654  varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue())
655 
656  # Get the mask and identify good pixels as '1', and the rest as '0'.
657  w1 = np.where(im1Area.getMask().getArray() == 0, 1, 0)
658  w2 = np.where(im2Area.getMask().getArray() == 0, 1, 0)
659 
660  w12 = w1*w2
661  wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0)
662  w = w12*wDiff
663 
664  maxRangeCov = self.config.maximumRangeCovariancesAstier
665  if covAstierRealSpace:
666  covDiffAstier = computeCovDirect(diffIm.getImage().getArray(), w, maxRangeCov)
667  else:
668  shapeDiff = diffIm.getImage().getArray().shape
669  fftShape = (fftSize(shapeDiff[0] + maxRangeCov), fftSize(shapeDiff[1]+maxRangeCov))
670  c = CovFft(diffIm.getImage().getArray(), w, fftShape, maxRangeCov)
671  covDiffAstier = c.reportCovFft(maxRangeCov)
672 
673  return mu, varDiff, covDiffAstier
674 
675  def computeCovDirect(self, diffImage, weightImage, maxRange):
676  """Compute covariances of diffImage in real space.
677 
678  For lags larger than ~25, it is slower than the FFT way.
679  Taken from https://github.com/PierreAstier/bfptc/
680 
681  Parameters
682  ----------
683  diffImage : `numpy.array`
684  Image to compute the covariance of.
685 
686  weightImage : `numpy.array`
687  Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
688 
689  maxRange : `int`
690  Last index of the covariance to be computed.
691 
692  Returns
693  -------
694  outList : `list`
695  List with tuples of the form (dx, dy, var, cov, npix), where:
696  dx : `int`
697  Lag in x
698  dy : `int`
699  Lag in y
700  var : `float`
701  Variance at (dx, dy).
702  cov : `float`
703  Covariance at (dx, dy).
704  nPix : `int`
705  Number of pixel pairs used to evaluate var and cov.
706  """
707  outList = []
708  var = 0
709  # (dy,dx) = (0,0) has to be first
710  for dy in range(maxRange + 1):
711  for dx in range(0, maxRange + 1):
712  if (dx*dy > 0):
713  cov1, nPix1 = self.covDirectValue(diffImage, weightImage, dx, dy)
714  cov2, nPix2 = self.covDirectValue(diffImage, weightImage, dx, -dy)
715  cov = 0.5*(cov1 + cov2)
716  nPix = nPix1 + nPix2
717  else:
718  cov, nPix = self.covDirectValue(diffImage, weightImage, dx, dy)
719  if (dx == 0 and dy == 0):
720  var = cov
721  outList.append((dx, dy, var, cov, nPix))
722 
723  return outList
724 
725  def covDirectValue(self, diffImage, weightImage, dx, dy):
726  """Compute covariances of diffImage in real space at lag (dx, dy).
727 
728  Taken from https://github.com/PierreAstier/bfptc/ (c.f., appendix of Astier+19).
729 
730  Parameters
731  ----------
732  diffImage : `numpy.array`
733  Image to compute the covariance of.
734 
735  weightImage : `numpy.array`
736  Weight image of diffImage (1's and 0's for good and bad pixels, respectively).
737 
738  dx : `int`
739  Lag in x.
740 
741  dy : `int`
742  Lag in y.
743 
744  Returns
745  -------
746  cov : `float`
747  Covariance at (dx, dy)
748 
749  nPix : `int`
750  Number of pixel pairs used to evaluate var and cov.
751  """
752  (nCols, nRows) = diffImage.shape
753  # switching both signs does not change anything:
754  # it just swaps im1 and im2 below
755  if (dx < 0):
756  (dx, dy) = (-dx, -dy)
757  # now, we have dx >0. We have to distinguish two cases
758  # depending on the sign of dy
759  if dy >= 0:
760  im1 = diffImage[dy:, dx:]
761  w1 = weightImage[dy:, dx:]
762  im2 = diffImage[:nCols - dy, :nRows - dx]
763  w2 = weightImage[:nCols - dy, :nRows - dx]
764  else:
765  im1 = diffImage[:nCols + dy, dx:]
766  w1 = weightImage[:nCols + dy, dx:]
767  im2 = diffImage[-dy:, :nRows - dx]
768  w2 = weightImage[-dy:, :nRows - dx]
769  # use the same mask for all 3 calculations
770  wAll = w1*w2
771  # do not use mean() because weightImage=0 pixels would then count
772  nPix = wAll.sum()
773  im1TimesW = im1*wAll
774  s1 = im1TimesW.sum()/nPix
775  s2 = (im2*wAll).sum()/nPix
776  p = (im1TimesW*im2).sum()/nPix
777  cov = p - s1*s2
778 
779  return cov, nPix
780 
781  def fitNonLinearityAndBuildLinearizers(self, datasetPtc, detector, tableArray=None, log=None):
782  """Fit non-linearity function and build linearizer objects.
783 
784  Parameters
785  ----------
786  datasePtct : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
787  The dataset containing information such as the means, variances and exposure times.
788  nLinearity
789 
790  detector : `lsst.afw.cameraGeom.Detector`
791  Detector object.
792 
793  tableArray : `np.array`, optional
794  Optional. Look-up table array with size rows=nAmps and columns=DN values.
795  It will be modified in-place if supplied.
796 
797  log : `lsst.log.Log`, optional
798  Logger to handle messages.
799 
800  Returns
801  -------
802  linearizer : `lsst.ip.isr.Linearizer`
803  Linearizer object
804  """
805 
806  # Fit NonLinearity
807  datasetNonLinearity = self.fitNonLinearity(datasetPtc, tableArray=tableArray)
808 
809  # Produce linearizer
810  now = datetime.datetime.utcnow()
811  calibDate = now.strftime("%Y-%m-%d")
812  linType = self.config.linearizerType
813 
814  if linType == "LOOKUPTABLE":
815  tableArray = tableArray
816  else:
817  tableArray = None
818 
819  linearizer = self.buildLinearizerObject(datasetNonLinearity, detector, calibDate, linType,
820  instruName=self.config.instrumentName,
821  tableArray=tableArray,
822  log=log)
823 
824  return linearizer
825 
826  def fitNonLinearity(self, datasetPtc, tableArray=None):
827  """Fit a polynomial to signal vs effective time curve to calculate linearity and residuals.
828 
829  Parameters
830  ----------
831  datasetPtc : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
832  The dataset containing the means, variances and exposure times.
833 
834  tableArray : `np.array`
835  Optional. Look-up table array with size rows=nAmps and columns=DN values.
836  It will be modified in-place if supplied.
837 
838  Returns
839  -------
840  datasetNonLinearity : `dict`
841  Dictionary of `lsst.cp.pipe.ptc.LinearityResidualsAndLinearizersDataset`
842  dataclasses. Each one holds the output of `calculateLinearityResidualAndLinearizers` per
843  amplifier.
844  """
845  datasetNonLinearity = {ampName: [] for ampName in datasetPtc.ampNames}
846  for i, ampName in enumerate(datasetPtc.ampNames):
847  # If a mask is not found, use all points.
848  if (len(datasetPtc.visitMask[ampName]) == 0):
849  self.log.warn(f"Mask not found for {ampName} in non-linearity fit. Using all points.")
850  mask = np.repeat(True, len(datasetPtc.rawExpTimes[ampName]))
851  else:
852  mask = datasetPtc.visitMask[ampName]
853 
854  timeVecFinal = np.array(datasetPtc.rawExpTimes[ampName])[mask]
855  meanVecFinal = np.array(datasetPtc.rawMeans[ampName])[mask]
856 
857  # Non-linearity residuals (NL of mean vs time curve): percentage, and fit to a quadratic function
858  # In this case, len(parsIniNonLinearity) = 3 indicates that we want a quadratic fit
859  datasetLinRes = self.calculateLinearityResidualAndLinearizers(timeVecFinal, meanVecFinal)
860 
861  # LinearizerLookupTable
862  if tableArray is not None:
863  tableArray[i, :] = datasetLinRes.linearizerTableRow
864 
865  datasetNonLinearity[ampName] = datasetLinRes
866 
867  return datasetNonLinearity
868 
869  def calculateLinearityResidualAndLinearizers(self, exposureTimeVector, meanSignalVector):
870  """Calculate linearity residual and fit an n-order polynomial to the mean vs time curve
871  to produce corrections (deviation from linear part of polynomial) for a particular amplifier
872  to populate LinearizeLookupTable.
873  Use the coefficients of this fit to calculate the correction coefficients for LinearizePolynomial
874  and LinearizeSquared."
875 
876  Parameters
877  ---------
878 
879  exposureTimeVector: `list` of `float`
880  List of exposure times for each flat pair
881 
882  meanSignalVector: `list` of `float`
883  List of mean signal from diference image of flat pairs
884 
885  Returns
886  -------
887  dataset : `lsst.cp.pipe.ptc.LinearityResidualsAndLinearizersDataset`
888  The dataset containing the fit parameters, the NL correction coefficients, and the
889  LUT row for the amplifier at hand.
890 
891  Notes
892  -----
893  datase members:
894 
895  dataset.polynomialLinearizerCoefficients : `list` of `float`
896  Coefficients for LinearizePolynomial, where corrImage = uncorrImage + sum_i c_i uncorrImage^(2 +
897  i).
898  c_(j-2) = -k_j/(k_1^j) with units DN^(1-j) (c.f., Eq. 37 of 2003.05978). The units of k_j are
899  DN/t^j, and they are fit from meanSignalVector = k0 + k1*exposureTimeVector +
900  k2*exposureTimeVector^2 + ... + kn*exposureTimeVector^n, with
901  n = "polynomialFitDegreeNonLinearity". k_0 and k_1 and degenerate with bias level and gain,
902  and are not used by the non-linearity correction. Therefore, j = 2...n in the above expression
903  (see `LinearizePolynomial` class in `linearize.py`.)
904 
905  dataset.quadraticPolynomialLinearizerCoefficient : `float`
906  Coefficient for LinearizeSquared, where corrImage = uncorrImage + c0*uncorrImage^2.
907  c0 = -k2/(k1^2), where k1 and k2 are fit from
908  meanSignalVector = k0 + k1*exposureTimeVector + k2*exposureTimeVector^2 +...
909  + kn*exposureTimeVector^n, with n = "polynomialFitDegreeNonLinearity".
910 
911  dataset.linearizerTableRow : `list` of `float`
912  One dimensional array with deviation from linear part of n-order polynomial fit
913  to mean vs time curve. This array will be one row (for the particular amplifier at hand)
914  of the table array for LinearizeLookupTable.
915 
916  dataset.meanSignalVsTimePolyFitPars : `list` of `float`
917  Parameters from n-order polynomial fit to meanSignalVector vs exposureTimeVector.
918 
919  dataset.meanSignalVsTimePolyFitParsErr : `list` of `float`
920  Parameters from n-order polynomial fit to meanSignalVector vs exposureTimeVector.
921 
922  dataset.meanSignalVsTimePolyFitReducedChiSq : `float`
923  Reduced unweighted chi squared from polynomial fit to meanSignalVector vs exposureTimeVector.
924  """
925 
926  # Lookup table linearizer
927  parsIniNonLinearity = self._initialParsForPolynomial(self.config.polynomialFitDegreeNonLinearity + 1)
928  if self.config.doFitBootstrap:
929  parsFit, parsFitErr, reducedChiSquaredNonLinearityFit = fitBootstrap(parsIniNonLinearity,
930  exposureTimeVector,
931  meanSignalVector,
932  funcPolynomial)
933  else:
934  parsFit, parsFitErr, reducedChiSquaredNonLinearityFit = fitLeastSq(parsIniNonLinearity,
935  exposureTimeVector,
936  meanSignalVector,
937  funcPolynomial)
938 
939  # LinearizeLookupTable:
940  # Use linear part to get time at wich signal is maxAduForLookupTableLinearizer DN
941  tMax = (self.config.maxAduForLookupTableLinearizer - parsFit[0])/parsFit[1]
942  timeRange = np.linspace(0, tMax, self.config.maxAduForLookupTableLinearizer)
943  signalIdeal = parsFit[0] + parsFit[1]*timeRange
944  signalUncorrected = funcPolynomial(parsFit, timeRange)
945  linearizerTableRow = signalIdeal - signalUncorrected # LinearizerLookupTable has corrections
946  # LinearizePolynomial and LinearizeSquared:
947  # Check that magnitude of higher order (>= 3) coefficents of the polyFit are small,
948  # i.e., less than threshold = 1e-10 (typical quadratic and cubic coefficents are ~1e-6
949  # and ~1e-12).
950  k1 = parsFit[1]
951  polynomialLinearizerCoefficients = []
952  for i, coefficient in enumerate(parsFit):
953  c = -coefficient/(k1**i)
954  polynomialLinearizerCoefficients.append(c)
955  if np.fabs(c) > 1e-10:
956  msg = f"Coefficient {c} in polynomial fit larger than threshold 1e-10."
957  self.log.warn(msg)
958  # Coefficient for LinearizedSquared. Called "c0" in linearize.py
959  c0 = polynomialLinearizerCoefficients[2]
960 
961  dataset = LinearityResidualsAndLinearizersDataset([], None, [], [], [], None)
962  dataset.polynomialLinearizerCoefficients = polynomialLinearizerCoefficients
963  dataset.quadraticPolynomialLinearizerCoefficient = c0
964  dataset.linearizerTableRow = linearizerTableRow
965  dataset.meanSignalVsTimePolyFitPars = parsFit
966  dataset.meanSignalVsTimePolyFitParsErr = parsFitErr
967  dataset.meanSignalVsTimePolyFitReducedChiSq = reducedChiSquaredNonLinearityFit
968 
969  return dataset
970 
971  def buildLinearizerObject(self, datasetNonLinearity, detector, calibDate, linearizerType, instruName='',
972  tableArray=None, log=None):
973  """Build linearizer object to persist.
974 
975  Parameters
976  ----------
977  datasetNonLinearity : `dict`
978  Dictionary of `lsst.cp.pipe.ptc.LinearityResidualsAndLinearizersDataset` objects.
979 
980  detector : `lsst.afw.cameraGeom.Detector`
981  Detector object
982 
983  calibDate : `datetime.datetime`
984  Calibration date
985 
986  linearizerType : `str`
987  'LOOKUPTABLE', 'LINEARIZESQUARED', or 'LINEARIZEPOLYNOMIAL'
988 
989  instruName : `str`, optional
990  Instrument name
991 
992  tableArray : `np.array`, optional
993  Look-up table array with size rows=nAmps and columns=DN values
994 
995  log : `lsst.log.Log`, optional
996  Logger to handle messages
997 
998  Returns
999  -------
1000  linearizer : `lsst.ip.isr.Linearizer`
1001  Linearizer object
1002  """
1003  detName = detector.getName()
1004  detNum = detector.getId()
1005  if linearizerType == "LOOKUPTABLE":
1006  if tableArray is not None:
1007  linearizer = Linearizer(detector=detector, table=tableArray, log=log)
1008  else:
1009  raise RuntimeError("tableArray must be provided when creating a LookupTable linearizer")
1010  elif linearizerType in ("LINEARIZESQUARED", "LINEARIZEPOLYNOMIAL"):
1011  linearizer = Linearizer(log=log)
1012  else:
1013  raise RuntimeError("Invalid linearizerType {linearizerType} to build a Linearizer object. "
1014  "Supported: 'LOOKUPTABLE', 'LINEARIZESQUARED', or 'LINEARIZEPOLYNOMIAL'")
1015  for i, amp in enumerate(detector.getAmplifiers()):
1016  ampName = amp.getName()
1017  datasetNonLinAmp = datasetNonLinearity[ampName]
1018  if linearizerType == "LOOKUPTABLE":
1019  linearizer.linearityCoeffs[ampName] = [i, 0]
1020  linearizer.linearityType[ampName] = "LookupTable"
1021  elif linearizerType == "LINEARIZESQUARED":
1022  linearizer.fitParams[ampName] = datasetNonLinAmp.meanSignalVsTimePolyFitPars
1023  linearizer.fitParamsErr[ampName] = datasetNonLinAmp.meanSignalVsTimePolyFitParsErr
1024  linearizer.linearityFitReducedChiSquared[ampName] = (
1025  datasetNonLinAmp.meanSignalVsTimePolyFitReducedChiSq)
1026  linearizer.linearityCoeffs[ampName] = [
1027  datasetNonLinAmp.quadraticPolynomialLinearizerCoefficient]
1028  linearizer.linearityType[ampName] = "Squared"
1029  elif linearizerType == "LINEARIZEPOLYNOMIAL":
1030  linearizer.fitParams[ampName] = datasetNonLinAmp.meanSignalVsTimePolyFitPars
1031  linearizer.fitParamsErr[ampName] = datasetNonLinAmp.meanSignalVsTimePolyFitParsErr
1032  linearizer.linearityFitReducedChiSquared[ampName] = (
1033  datasetNonLinAmp.meanSignalVsTimePolyFitReducedChiSq)
1034  # Slice correction coefficients (starting at 2) for polynomial linearizer
1035  # (and squared linearizer above). The first and second are reduntant with
1036  # the bias and gain, respectively, and are not used by LinearizerPolynomial.
1037  polyLinCoeffs = np.array(datasetNonLinAmp.polynomialLinearizerCoefficients[2:])
1038  linearizer.linearityCoeffs[ampName] = polyLinCoeffs
1039  linearizer.linearityType[ampName] = "Polynomial"
1040  linearizer.linearityBBox[ampName] = amp.getBBox()
1041  linearizer.validate()
1042  calibId = f"detectorName={detName} detector={detNum} calibDate={calibDate} ccd={detNum} filter=NONE"
1043 
1044  try:
1045  raftName = detName.split("_")[0]
1046  calibId += f" raftName={raftName}"
1047  except Exception:
1048  raftname = "NONE"
1049  calibId += f" raftName={raftname}"
1050 
1051  serial = detector.getSerial()
1052  linearizer.updateMetadata(instrumentName=instruName, detectorId=f"{detNum}",
1053  calibId=calibId, serial=serial, detectorName=f"{detName}")
1054 
1055  return linearizer
1056 
1057  @staticmethod
1058  def _initialParsForPolynomial(order):
1059  assert(order >= 2)
1060  pars = np.zeros(order, dtype=np.float)
1061  pars[0] = 10
1062  pars[1] = 1
1063  pars[2:] = 0.0001
1064  return pars
1065 
1066  @staticmethod
1067  def _boundsForPolynomial(initialPars):
1068  lowers = [np.NINF for p in initialPars]
1069  uppers = [np.inf for p in initialPars]
1070  lowers[1] = 0 # no negative gains
1071  return (lowers, uppers)
1072 
1073  @staticmethod
1074  def _boundsForAstier(initialPars):
1075  lowers = [np.NINF for p in initialPars]
1076  uppers = [np.inf for p in initialPars]
1077  return (lowers, uppers)
1078 
1079  @staticmethod
1080  def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative):
1081  """Return a boolean array to mask bad points.
1082 
1083  A linear function has a constant ratio, so find the median
1084  value of the ratios, and exclude the points that deviate
1085  from that by more than a factor of maxDeviationPositive/negative.
1086  Asymmetric deviations are supported as we expect the PTC to turn
1087  down as the flux increases, but sometimes it anomalously turns
1088  upwards just before turning over, which ruins the fits, so it
1089  is wise to be stricter about restricting positive outliers than
1090  negative ones.
1091 
1092  Too high and points that are so bad that fit will fail will be included
1093  Too low and the non-linear points will be excluded, biasing the NL fit."""
1094  ratios = [b/a for (a, b) in zip(means, variances)]
1095  medianRatio = np.median(ratios)
1096  ratioDeviations = [(r/medianRatio)-1 for r in ratios]
1097 
1098  # so that it doesn't matter if the deviation is expressed as positive or negative
1099  maxDeviationPositive = abs(maxDeviationPositive)
1100  maxDeviationNegative = -1. * abs(maxDeviationNegative)
1101 
1102  goodPoints = np.array([True if (r < maxDeviationPositive and r > maxDeviationNegative)
1103  else False for r in ratioDeviations])
1104  return goodPoints
1105 
1106  def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9):
1107  """"""
1108  nBad = Counter(array)[0]
1109  if nBad == 0:
1110  return array
1111 
1112  if warn:
1113  msg = f"Found {nBad} zeros in array at elements {[x for x in np.where(array==0)[0]]}"
1114  self.log.warn(msg)
1115 
1116  array[array == 0] = substituteValue
1117  return array
1118 
1119  def fitPtc(self, dataset, ptcFitType):
1120  """Fit the photon transfer curve to a polynimial or to Astier+19 approximation.
1121 
1122  Fit the photon transfer curve with either a polynomial of the order
1123  specified in the task config, or using the Astier approximation.
1124 
1125  Sigma clipping is performed iteratively for the fit, as well as an
1126  initial clipping of data points that are more than
1127  config.initialNonLinearityExclusionThreshold away from lying on a
1128  straight line. This other step is necessary because the photon transfer
1129  curve turns over catastrophically at very high flux (because saturation
1130  drops the variance to ~0) and these far outliers cause the initial fit
1131  to fail, meaning the sigma cannot be calculated to perform the
1132  sigma-clipping.
1133 
1134  Parameters
1135  ----------
1136  dataset : `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
1137  The dataset containing the means, variances and exposure times
1138 
1139  ptcFitType : `str`
1140  Fit a 'POLYNOMIAL' (degree: 'polynomialFitDegree') or
1141  'EXPAPPROXIMATION' (Eq. 16 of Astier+19) to the PTC
1142 
1143  Returns
1144  -------
1145  dataset: `lsst.cp.pipe.ptc.PhotonTransferCurveDataset`
1146  This is the same dataset as the input paramter, however, it has been modified
1147  to include information such as the fit vectors and the fit parameters. See
1148  the class `PhotonTransferCurveDatase`.
1149  """
1150 
1151  def errFunc(p, x, y):
1152  return ptcFunc(p, x) - y
1153 
1154  sigmaCutPtcOutliers = self.config.sigmaCutPtcOutliers
1155  maxIterationsPtcOutliers = self.config.maxIterationsPtcOutliers
1156 
1157  for i, ampName in enumerate(dataset.ampNames):
1158  timeVecOriginal = np.array(dataset.rawExpTimes[ampName])
1159  meanVecOriginal = np.array(dataset.rawMeans[ampName])
1160  varVecOriginal = np.array(dataset.rawVars[ampName])
1161  varVecOriginal = self._makeZeroSafe(varVecOriginal)
1162 
1163  mask = ((meanVecOriginal >= self.config.minMeanSignal) &
1164  (meanVecOriginal <= self.config.maxMeanSignal))
1165 
1166  goodPoints = self._getInitialGoodPoints(meanVecOriginal, varVecOriginal,
1167  self.config.initialNonLinearityExclusionThresholdPositive,
1168  self.config.initialNonLinearityExclusionThresholdNegative)
1169  mask = mask & goodPoints
1170 
1171  if ptcFitType == 'EXPAPPROXIMATION':
1172  ptcFunc = funcAstier
1173  parsIniPtc = [-1e-9, 1.0, 10.] # a00, gain, noise
1174  bounds = self._boundsForAstier(parsIniPtc)
1175  if ptcFitType == 'POLYNOMIAL':
1176  ptcFunc = funcPolynomial
1177  parsIniPtc = self._initialParsForPolynomial(self.config.polynomialFitDegree + 1)
1178  bounds = self._boundsForPolynomial(parsIniPtc)
1179 
1180  # Before bootstrap fit, do an iterative fit to get rid of outliers
1181  count = 1
1182  while count <= maxIterationsPtcOutliers:
1183  # Note that application of the mask actually shrinks the array
1184  # to size rather than setting elements to zero (as we want) so
1185  # always update mask itself and re-apply to the original data
1186  meanTempVec = meanVecOriginal[mask]
1187  varTempVec = varVecOriginal[mask]
1188  res = least_squares(errFunc, parsIniPtc, bounds=bounds, args=(meanTempVec, varTempVec))
1189  pars = res.x
1190 
1191  # change this to the original from the temp because the masks are ANDed
1192  # meaning once a point is masked it's always masked, and the masks must
1193  # always be the same length for broadcasting
1194  sigResids = (varVecOriginal - ptcFunc(pars, meanVecOriginal))/np.sqrt(varVecOriginal)
1195  newMask = np.array([True if np.abs(r) < sigmaCutPtcOutliers else False for r in sigResids])
1196  mask = mask & newMask
1197 
1198  nDroppedTotal = Counter(mask)[False]
1199  self.log.debug(f"Iteration {count}: discarded {nDroppedTotal} points in total for {ampName}")
1200  count += 1
1201  # objects should never shrink
1202  assert (len(mask) == len(timeVecOriginal) == len(meanVecOriginal) == len(varVecOriginal))
1203 
1204  dataset.visitMask[ampName] = mask # store the final mask
1205  parsIniPtc = pars
1206  meanVecFinal = meanVecOriginal[mask]
1207  varVecFinal = varVecOriginal[mask]
1208 
1209  if Counter(mask)[False] > 0:
1210  self.log.info((f"Number of points discarded in PTC of amplifier {ampName}:" +
1211  f" {Counter(mask)[False]} out of {len(meanVecOriginal)}"))
1212 
1213  if (len(meanVecFinal) < len(parsIniPtc)):
1214  msg = (f"\nSERIOUS: Not enough data points ({len(meanVecFinal)}) compared to the number of"
1215  f"parameters of the PTC model({len(parsIniPtc)}). Setting {ampName} to BAD.")
1216  self.log.warn(msg)
1217  # The first and second parameters of initial fit are discarded (bias and gain)
1218  # for the final NL coefficients
1219  dataset.badAmps.append(ampName)
1220  dataset.gain[ampName] = np.nan
1221  dataset.gainErr[ampName] = np.nan
1222  dataset.noise[ampName] = np.nan
1223  dataset.noiseErr[ampName] = np.nan
1224  dataset.ptcFitPars[ampName] = np.nan
1225  dataset.ptcFitParsError[ampName] = np.nan
1226  dataset.ptcFitReducedChiSquared[ampName] = np.nan
1227  continue
1228 
1229  # Fit the PTC
1230  if self.config.doFitBootstrap:
1231  parsFit, parsFitErr, reducedChiSqPtc = fitBootstrap(parsIniPtc, meanVecFinal,
1232  varVecFinal, ptcFunc)
1233  else:
1234  parsFit, parsFitErr, reducedChiSqPtc = fitLeastSq(parsIniPtc, meanVecFinal,
1235  varVecFinal, ptcFunc)
1236  dataset.ptcFitPars[ampName] = parsFit
1237  dataset.ptcFitParsError[ampName] = parsFitErr
1238  dataset.ptcFitReducedChiSquared[ampName] = reducedChiSqPtc
1239 
1240  if ptcFitType == 'EXPAPPROXIMATION':
1241  ptcGain = parsFit[1]
1242  ptcGainErr = parsFitErr[1]
1243  ptcNoise = np.sqrt(np.fabs(parsFit[2]))
1244  ptcNoiseErr = 0.5*(parsFitErr[2]/np.fabs(parsFit[2]))*np.sqrt(np.fabs(parsFit[2]))
1245  if ptcFitType == 'POLYNOMIAL':
1246  ptcGain = 1./parsFit[1]
1247  ptcGainErr = np.fabs(1./parsFit[1])*(parsFitErr[1]/parsFit[1])
1248  ptcNoise = np.sqrt(np.fabs(parsFit[0]))*ptcGain
1249  ptcNoiseErr = (0.5*(parsFitErr[0]/np.fabs(parsFit[0]))*(np.sqrt(np.fabs(parsFit[0]))))*ptcGain
1250  dataset.gain[ampName] = ptcGain
1251  dataset.gainErr[ampName] = ptcGainErr
1252  dataset.noise[ampName] = ptcNoise
1253  dataset.noiseErr[ampName] = ptcNoiseErr
1254  if not len(dataset.ptcFitType) == 0:
1255  dataset.ptcFitType = ptcFitType
1256 
1257  return dataset
lsst::ip::isr::linearize::Linearizer
lsst.cp.pipe.ptc.PhotonTransferCurveDataset.__setattr__
def __setattr__(self, attribute, value)
Definition: ptc.py:271
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.fitCovariancesAstier
def fitCovariancesAstier(self, dataset, covariancesWithTagsArray)
Definition: ptc.py:486
lsst.cp.pipe.ptc.PhotonTransferCurveDataset
Definition: ptc.py:192
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.fitPtc
def fitPtc(self, dataset, ptcFitType)
Definition: ptc.py:1119
lsst.cp.pipe.utils.fitBootstrap
def fitBootstrap(initialParams, dataX, dataY, function, confidenceSigma=1.)
Definition: utils.py:222
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.fitNonLinearityAndBuildLinearizers
def fitNonLinearityAndBuildLinearizers(self, datasetPtc, detector, tableArray=None, log=None)
Definition: ptc.py:781
lsst.cp.pipe.astierCovPtcUtils.fftSize
def fftSize(s)
Definition: astierCovPtcUtils.py:122
lsst.cp.pipe.utils.funcPolynomial
def funcPolynomial(pars, x)
Definition: utils.py:288
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._boundsForAstier
def _boundsForAstier(initialPars)
Definition: ptc.py:1074
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.measureMeanVarCov
def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False)
Definition: ptc.py:560
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.__init__
def __init__(self, *args, **kwargs)
Definition: ptc.py:338
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._boundsForPolynomial
def _boundsForPolynomial(initialPars)
Definition: ptc.py:1067
lsst.cp.pipe.utils.fitLeastSq
def fitLeastSq(initialParams, dataX, dataY, function)
Definition: utils.py:164
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.getOutputPtcDataCovAstier
def getOutputPtcDataCovAstier(self, dataset, covFits)
Definition: ptc.py:523
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTaskConfig
Definition: ptc.py:46
lsst.cp.pipe.ptc.PhotonTransferCurveDataset.getGoodAmps
def getGoodAmps(self)
Definition: ptc.py:295
lsst.cp.pipe.ptc.PhotonTransferCurveDataset.getVisitsUsed
def getVisitsUsed(self, ampName)
Definition: ptc.py:279
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.covDirectValue
def covDirectValue(self, diffImage, weightImage, dx, dy)
Definition: ptc.py:725
lsst.cp.pipe.astierCovPtcUtils.fitData
def fitData(tupleName, maxMu=1e9, r=8)
Definition: astierCovPtcUtils.py:323
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._DefaultName
string _DefaultName
Definition: ptc.py:336
lsst.cp.pipe.astierCovPtcUtils.CovFft
Definition: astierCovPtcUtils.py:28
lsst.cp.pipe.utils.checkExpLengthEqual
def checkExpLengthEqual(exp1, exp2, v1=None, v2=None, raiseWithMessage=False)
Definition: utils.py:324
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask
Definition: ptc.py:299
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._makeZeroSafe
def _makeZeroSafe(self, array, warn=True, substituteValue=1e-9)
Definition: ptc.py:1106
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._getInitialGoodPoints
def _getInitialGoodPoints(means, variances, maxDeviationPositive, maxDeviationNegative)
Definition: ptc.py:1080
lsst::ip::isr::linearize
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.runDataRef
def runDataRef(self, dataRef, visitPairs)
Definition: ptc.py:356
lsst::afw::math
lsst.cp.pipe.ptc.LinearityResidualsAndLinearizersDataset
Definition: ptc.py:177
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask._initialParsForPolynomial
def _initialParsForPolynomial(order)
Definition: ptc.py:1058
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.computeCovDirect
def computeCovDirect(self, diffImage, weightImage, maxRange)
Definition: ptc.py:675
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.fitNonLinearity
def fitNonLinearity(self, datasetPtc, tableArray=None)
Definition: ptc.py:826
lsst::pipe::base
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.calculateLinearityResidualAndLinearizers
def calculateLinearityResidualAndLinearizers(self, exposureTimeVector, meanSignalVector)
Definition: ptc.py:869
lsst.cp.pipe.ptc.MeasurePhotonTransferCurveTask.buildLinearizerObject
def buildLinearizerObject(self, datasetNonLinearity, detector, calibDate, linearizerType, instruName='', tableArray=None, log=None)
Definition: ptc.py:971
lsst.cp.pipe.ptc.PhotonTransferCurveDataset.__init__
def __init__(self, ampNames, ptcFitType)
Definition: ptc.py:226