lsst.cp.pipe  21.0.0-16-g434b3eb+fcd236b602
cpExtractPtcTask.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 
24 import lsst.afw.math as afwMath
25 import lsst.pex.config as pexConfig
26 import lsst.pipe.base as pipeBase
27 from lsst.cp.pipe.utils import (arrangeFlatsByExpTime, arrangeFlatsByExpId,
28  sigmaClipCorrection)
29 
31 
32 from .astierCovPtcUtils import (CovFastFourierTransform, computeCovDirect)
33 from .astierCovPtcFit import makeCovArray
34 
35 from lsst.ip.isr import PhotonTransferCurveDataset
36 
37 
38 __all__ = ['PhotonTransferCurveExtractConfig', 'PhotonTransferCurveExtractTask']
39 
40 
41 class PhotonTransferCurveExtractConnections(pipeBase.PipelineTaskConnections,
42  dimensions=("instrument", "detector")):
43 
44  inputExp = cT.Input(
45  name="ptcInputExposurePairs",
46  doc="Input post-ISR processed exposure pairs (flats) to"
47  "measure covariances from.",
48  storageClass="Exposure",
49  dimensions=("instrument", "exposure", "detector"),
50  multiple=True,
51  deferLoad=False,
52  )
53 
54  outputCovariances = cT.Output(
55  name="ptcCovariances",
56  doc="Extracted flat (co)variances.",
57  storageClass="PhotonTransferCurveDataset",
58  dimensions=("instrument", "exposure", "detector"),
59  multiple=True,
60  )
61 
62 
63 class PhotonTransferCurveExtractConfig(pipeBase.PipelineTaskConfig,
64  pipelineConnections=PhotonTransferCurveExtractConnections):
65  """Configuration for the measurement of covariances from flats.
66  """
67  matchByExposureId = pexConfig.Field(
68  dtype=bool,
69  doc="Should exposures by matched by ID rather than exposure time?",
70  default=False,
71  )
72  maximumRangeCovariancesAstier = pexConfig.Field(
73  dtype=int,
74  doc="Maximum range of covariances as in Astier+19",
75  default=8,
76  )
77  covAstierRealSpace = pexConfig.Field(
78  dtype=bool,
79  doc="Calculate covariances in real space or via FFT? (see appendix A of Astier+19).",
80  default=False,
81  )
82  binSize = pexConfig.Field(
83  dtype=int,
84  doc="Bin the image by this factor in both dimensions.",
85  default=1,
86  )
87  minMeanSignal = pexConfig.DictField(
88  keytype=str,
89  itemtype=float,
90  doc="Minimum values (inclusive) of mean signal (in ADU) above which to consider, per amp."
91  " The same cut is applied to all amps if this dictionary is of the form"
92  " {'ALL_AMPS': value}",
93  default={'ALL_AMPS': 0.0},
94  )
95  maxMeanSignal = pexConfig.DictField(
96  keytype=str,
97  itemtype=float,
98  doc="Maximum values (inclusive) of mean signal (in ADU) below 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': 1e6},
102  )
103  maskNameList = pexConfig.ListField(
104  dtype=str,
105  doc="Mask list to exclude from statistics calculations.",
106  default=['SUSPECT', 'BAD', 'NO_DATA', 'SAT'],
107  )
108  nSigmaClipPtc = pexConfig.Field(
109  dtype=float,
110  doc="Sigma cut for afwMath.StatisticsControl()",
111  default=5.5,
112  )
113  nIterSigmaClipPtc = pexConfig.Field(
114  dtype=int,
115  doc="Number of sigma-clipping iterations for afwMath.StatisticsControl()",
116  default=3,
117  )
118  minNumberGoodPixelsForCovariance = pexConfig.Field(
119  dtype=int,
120  doc="Minimum number of acceptable good pixels per amp to calculate the covariances (via FFT or"
121  " direclty).",
122  default=10000,
123  )
124  thresholdDiffAfwVarVsCov00 = pexConfig.Field(
125  dtype=float,
126  doc="If the absolute fractional differece between afwMath.VARIANCECLIP and Cov00 "
127  "for a region of a difference image is greater than this threshold (percentage), "
128  "a warning will be issued.",
129  default=1.,
130  )
131  detectorMeasurementRegion = pexConfig.ChoiceField(
132  dtype=str,
133  doc="Region of each exposure where to perform the calculations (amplifier or full image).",
134  default='AMP',
135  allowed={
136  "AMP": "Amplifier of the detector.",
137  "FULL": "Full image."
138  }
139  )
140 
141 
142 class PhotonTransferCurveExtractTask(pipeBase.PipelineTask,
143  pipeBase.CmdLineTask):
144  """Task to measure covariances from flat fields.
145  This task receives as input a list of flat-field images
146  (flats), and sorts these flats in pairs taken at the
147  same time (if there's a different number of flats,
148  those flats are discarded). The mean, variance, and
149  covariances are measured from the difference of the flat
150  pairs at a given time. The variance is calculated
151  via afwMath, and the covariance via the methods in Astier+19
152  (appendix A). In theory, var = covariance[0,0]. This should
153  be validated, and in the future, we may decide to just keep
154  one (covariance).
155 
156  The measured covariances at a particular time (along with
157  other quantities such as the mean) are stored in a PTC dataset
158  object (`PhotonTransferCurveDataset`), which gets partially
159  filled. The number of partially-filled PTC dataset objects
160  will be less than the number of input exposures, but gen3
161  requires/assumes that the number of input dimensions matches
162  bijectively the number of output dimensions. Therefore, a
163  number of "dummy" PTC dataset are inserted in the output list
164  that has the partially-filled PTC datasets with the covariances.
165  This output list will be used as input of
166  `PhotonTransferCurveSolveTask`, which will assemble the multiple
167  `PhotonTransferCurveDataset`s into a single one in order to fit
168  the measured covariances as a function of flux to a particular
169  model.
170 
171  Astier+19: "The Shape of the Photon Transfer Curve of CCD
172  sensors", arXiv:1905.08677.
173  """
174  ConfigClass = PhotonTransferCurveExtractConfig
175  _DefaultName = 'cpPtcExtract'
176 
177  def runQuantum(self, butlerQC, inputRefs, outputRefs):
178  """Ensure that the input and output dimensions are passed along.
179 
180  Parameters
181  ----------
182  butlerQC : `~lsst.daf.butler.butlerQuantumContext.ButlerQuantumContext`
183  Butler to operate on.
184  inputRefs : `~lsst.pipe.base.connections.InputQuantizedConnection`
185  Input data refs to load.
186  ouptutRefs : `~lsst.pipe.base.connections.OutputQuantizedConnection`
187  Output data refs to persist.
188  """
189  inputs = butlerQC.get(inputRefs)
190  # Dictionary, keyed by expTime, with flat exposures
191  if self.config.matchByExposureId:
192  inputs['inputExp'] = arrangeFlatsByExpId(inputs['inputExp'])
193  else:
194  inputs['inputExp'] = arrangeFlatsByExpTime(inputs['inputExp'])
195  # Ids of input list of exposures
196  inputs['inputDims'] = [expId.dataId['exposure'] for expId in inputRefs.inputExp]
197  outputs = self.runrun(**inputs)
198  butlerQC.put(outputs, outputRefs)
199 
200  def run(self, inputExp, inputDims):
201  """Measure covariances from difference of flat pairs
202 
203  Parameters
204  ----------
205  inputExp : `dict` [`float`,
206  (`~lsst.afw.image.exposure.exposure.ExposureF`,
207  `~lsst.afw.image.exposure.exposure.ExposureF`, ...,
208  `~lsst.afw.image.exposure.exposure.ExposureF`)]
209  Dictionary that groups flat-field exposures that have the same
210  exposure time (seconds).
211 
212  inputDims : `list`
213  List of exposure IDs.
214  """
215  # inputExp.values() returns a view, which we turn into a list. We then
216  # access the first exposure to get teh detector.
217  detector = list(inputExp.values())[0][0].getDetector()
218  detNum = detector.getId()
219  amps = detector.getAmplifiers()
220  ampNames = [amp.getName() for amp in amps]
221 
222  # Each amp may have a different min and max ADU signal specified in the config.
223  maxMeanSignalDict = {ampName: 1e6 for ampName in ampNames}
224  minMeanSignalDict = {ampName: 0.0 for ampName in ampNames}
225  for ampName in ampNames:
226  if 'ALL_AMPS' in self.config.maxMeanSignal:
227  maxMeanSignalDict[ampName] = self.config.maxMeanSignal['ALL_AMPS']
228  elif ampName in self.config.maxMeanSignal:
229  maxMeanSignalDict[ampName] = self.config.maxMeanSignal[ampName]
230 
231  if 'ALL_AMPS' in self.config.minMeanSignal:
232  minMeanSignalDict[ampName] = self.config.minMeanSignal['ALL_AMPS']
233  elif ampName in self.config.minMeanSignal:
234  minMeanSignalDict[ampName] = self.config.minMeanSignal[ampName]
235  # These are the column names for `tupleRows` below.
236  tags = [('mu', '<f8'), ('afwVar', '<f8'), ('i', '<i8'), ('j', '<i8'), ('var', '<f8'),
237  ('cov', '<f8'), ('npix', '<i8'), ('ext', '<i8'), ('expTime', '<f8'), ('ampName', '<U3')]
238  # Create a dummy ptcDataset
239  dummyPtcDataset = PhotonTransferCurveDataset(ampNames, 'DUMMY',
240  self.config.maximumRangeCovariancesAstier)
241  # Initialize amps of `dummyPtcDatset`.
242  for ampName in ampNames:
243  dummyPtcDataset.setAmpValues(ampName)
244  # Output list with PTC datasets.
245  partialPtcDatasetList = []
246  # The number of output references needs to match that of input references:
247  # initialize outputlist with dummy PTC datasets.
248  for i in range(len(inputDims)):
249  partialPtcDatasetList.append(dummyPtcDataset)
250 
251  for expTime in inputExp:
252  exposures = inputExp[expTime]
253  if len(exposures) == 1:
254  self.log.warn(f"Only one exposure found at expTime {expTime}. Dropping exposure "
255  f"{exposures[0].getInfo().getVisitInfo().getExposureId()}.")
256  continue
257  else:
258  # Only use the first two exposures at expTime
259  exp1, exp2 = exposures[0], exposures[1]
260  if len(exposures) > 2:
261  self.log.warn(f"Already found 2 exposures at expTime {expTime}. "
262  "Ignoring exposures: "
263  f"{i.getInfo().getVisitInfo().getExposureId() for i in exposures[2:]}")
264  expId1 = exp1.getInfo().getVisitInfo().getExposureId()
265  expId2 = exp2.getInfo().getVisitInfo().getExposureId()
266  nAmpsNan = 0
267  partialPtcDataset = PhotonTransferCurveDataset(ampNames, '',
268  self.config.maximumRangeCovariancesAstier)
269  for ampNumber, amp in enumerate(detector):
270  ampName = amp.getName()
271  # covAstier: [(i, j, var (cov[0,0]), cov, npix) for (i,j) in {maxLag, maxLag}^2]
272  doRealSpace = self.config.covAstierRealSpace
273  if self.config.detectorMeasurementRegion == 'AMP':
274  region = amp.getBBox()
275  elif self.config.detectorMeasurementRegion == 'FULL':
276  region = None
277  # `measureMeanVarCov` is the function that measures the variance and covariances from
278  # the difference image of two flats at the same exposure time.
279  # The variable `covAstier` is of the form: [(i, j, var (cov[0,0]), cov, npix) for (i,j)
280  # in {maxLag, maxLag}^2]
281  muDiff, varDiff, covAstier = self.measureMeanVarCovmeasureMeanVarCov(exp1, exp2, region=region,
282  covAstierRealSpace=doRealSpace)
283  # Correction factor for sigma clipping. Function returns 1/sqrt(varFactor),
284  # so it needs to be squared. varDiff is calculated via afwMath.VARIANCECLIP.
285  varFactor = sigmaClipCorrection(self.config.nSigmaClipPtc)**2
286  varDiff *= varFactor
287 
288  expIdMask = True
289  if np.isnan(muDiff) or np.isnan(varDiff) or (covAstier is None):
290  msg = (f"NaN mean or var, or None cov in amp {ampName} in exposure pair {expId1},"
291  f" {expId2} of detector {detNum}.")
292  self.log.warn(msg)
293  nAmpsNan += 1
294  expIdMask = False
295  covArray = np.full((1, self.config.maximumRangeCovariancesAstier,
296  self.config.maximumRangeCovariancesAstier), np.nan)
297  covSqrtWeights = np.full_like(covArray, np.nan)
298 
299  if (muDiff <= minMeanSignalDict[ampName]) or (muDiff >= maxMeanSignalDict[ampName]):
300  expIdMask = False
301 
302  if covAstier is not None:
303  tupleRows = [(muDiff, varDiff) + covRow + (ampNumber, expTime,
304  ampName) for covRow in covAstier]
305  tempStructArray = np.array(tupleRows, dtype=tags)
306  covArray, vcov, _ = makeCovArray(tempStructArray,
307  self.config.maximumRangeCovariancesAstier)
308  covSqrtWeights = np.nan_to_num(1./np.sqrt(vcov))
309 
310  # Correct covArray for sigma clipping:
311  # 1) Apply varFactor twice for the whole covariance matrix
312  covArray *= varFactor**2
313  # 2) But, only once for the variance element of the matrix, covArray[0,0]
314  covArray[0, 0] /= varFactor
315 
316  partialPtcDataset.setAmpValues(ampName, rawExpTime=[expTime], rawMean=[muDiff],
317  rawVar=[varDiff], inputExpIdPair=[(expId1, expId2)],
318  expIdMask=[expIdMask], covArray=covArray,
319  covSqrtWeights=covSqrtWeights)
320  # Use location of exp1 to save PTC dataset from (exp1, exp2) pair.
321  # expId1 and expId2, as returned by getInfo().getVisitInfo().getExposureId(),
322  # and the exposure IDs stured in inoutDims,
323  # may have the zero-padded detector number appended at
324  # the end (in gen3). A temporary fix is to consider expId//1000 and/or
325  # inputDims//1000.
326  # Below, np.where(expId1 == np.array(inputDims)) (and the other analogous
327  # comparisons) returns a tuple with a single-element array, so [0][0]
328  # is necessary to extract the required index.
329  try:
330  datasetIndex = np.where(expId1 == np.array(inputDims))[0][0]
331  except IndexError:
332  try:
333  datasetIndex = np.where(expId1//1000 == np.array(inputDims))[0][0]
334  except IndexError:
335  datasetIndex = np.where(expId1//1000 == np.array(inputDims)//1000)[0][0]
336  partialPtcDatasetList[datasetIndex] = partialPtcDataset
337  if nAmpsNan == len(ampNames):
338  msg = f"NaN mean in all amps of exposure pair {expId1}, {expId2} of detector {detNum}."
339  self.log.warn(msg)
340  return pipeBase.Struct(
341  outputCovariances=partialPtcDatasetList,
342  )
343 
344  def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False):
345  """Calculate the mean of each of two exposures and the variance
346  and covariance of their difference. The variance is calculated
347  via afwMath, and the covariance via the methods in Astier+19
348  (appendix A). In theory, var = covariance[0,0]. This should
349  be validated, and in the future, we may decide to just keep
350  one (covariance).
351 
352  Parameters
353  ----------
354  exposure1 : `lsst.afw.image.exposure.exposure.ExposureF`
355  First exposure of flat field pair.
356  exposure2 : `lsst.afw.image.exposure.exposure.ExposureF`
357  Second exposure of flat field pair.
358  region : `lsst.geom.Box2I`, optional
359  Region of each exposure where to perform the calculations (e.g, an amplifier).
360  covAstierRealSpace : `bool`, optional
361  Should the covariannces in Astier+19 be calculated in real space or via FFT?
362  See Appendix A of Astier+19.
363 
364  Returns
365  -------
366  mu : `float` or `NaN`
367  0.5*(mu1 + mu2), where mu1, and mu2 are the clipped means of the regions in
368  both exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
369  varDiff : `float` or `NaN`
370  Half of the clipped variance of the difference of the regions inthe two input
371  exposures. If either mu1 or m2 are NaN's, the returned value is NaN.
372  covDiffAstier : `list` or `NaN`
373  List with tuples of the form (dx, dy, var, cov, npix), where:
374  dx : `int`
375  Lag in x
376  dy : `int`
377  Lag in y
378  var : `float`
379  Variance at (dx, dy).
380  cov : `float`
381  Covariance at (dx, dy).
382  nPix : `int`
383  Number of pixel pairs used to evaluate var and cov.
384  If either mu1 or m2 are NaN's, the returned value is NaN.
385  """
386 
387  if region is not None:
388  im1Area = exposure1.maskedImage[region]
389  im2Area = exposure2.maskedImage[region]
390  else:
391  im1Area = exposure1.maskedImage
392  im2Area = exposure2.maskedImage
393 
394  if self.config.binSize > 1:
395  im1Area = afwMath.binImage(im1Area, self.config.binSize)
396  im2Area = afwMath.binImage(im2Area, self.config.binSize)
397 
398  im1MaskVal = exposure1.getMask().getPlaneBitMask(self.config.maskNameList)
399  im1StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
400  self.config.nIterSigmaClipPtc,
401  im1MaskVal)
402  im1StatsCtrl.setNanSafe(True)
403  im1StatsCtrl.setAndMask(im1MaskVal)
404 
405  im2MaskVal = exposure2.getMask().getPlaneBitMask(self.config.maskNameList)
406  im2StatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
407  self.config.nIterSigmaClipPtc,
408  im2MaskVal)
409  im2StatsCtrl.setNanSafe(True)
410  im2StatsCtrl.setAndMask(im2MaskVal)
411 
412  # Clipped mean of images; then average of mean.
413  mu1 = afwMath.makeStatistics(im1Area, afwMath.MEANCLIP, im1StatsCtrl).getValue()
414  mu2 = afwMath.makeStatistics(im2Area, afwMath.MEANCLIP, im2StatsCtrl).getValue()
415  if np.isnan(mu1) or np.isnan(mu2):
416  self.log.warn(f"Mean of amp in image 1 or 2 is NaN: {mu1}, {mu2}.")
417  return np.nan, np.nan, None
418  mu = 0.5*(mu1 + mu2)
419 
420  # Take difference of pairs
421  # symmetric formula: diff = (mu2*im1-mu1*im2)/(0.5*(mu1+mu2))
422  temp = im2Area.clone()
423  temp *= mu1
424  diffIm = im1Area.clone()
425  diffIm *= mu2
426  diffIm -= temp
427  diffIm /= mu
428 
429  diffImMaskVal = diffIm.getMask().getPlaneBitMask(self.config.maskNameList)
430  diffImStatsCtrl = afwMath.StatisticsControl(self.config.nSigmaClipPtc,
431  self.config.nIterSigmaClipPtc,
432  diffImMaskVal)
433  diffImStatsCtrl.setNanSafe(True)
434  diffImStatsCtrl.setAndMask(diffImMaskVal)
435 
436  # Variance calculation via afwMath
437  varDiff = 0.5*(afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue())
438 
439  # Covariances calculations
440  # Get the pixels that were not clipped
441  varClip = afwMath.makeStatistics(diffIm, afwMath.VARIANCECLIP, diffImStatsCtrl).getValue()
442  meanClip = afwMath.makeStatistics(diffIm, afwMath.MEANCLIP, diffImStatsCtrl).getValue()
443  cut = meanClip + self.config.nSigmaClipPtc*np.sqrt(varClip)
444  unmasked = np.where(np.fabs(diffIm.image.array) <= cut, 1, 0)
445 
446  # Get the pixels in the mask planes of teh differenc eimage that were ignored
447  # by the clipping algorithm
448  wDiff = np.where(diffIm.getMask().getArray() == 0, 1, 0)
449  # Combine the two sets of pixels ('1': use; '0': don't use) into a final weight matrix
450  # to be used in the covariance calculations below.
451  w = unmasked*wDiff
452 
453  if np.sum(w) < self.config.minNumberGoodPixelsForCovariance:
454  self.log.warn(f"Number of good points for covariance calculation ({np.sum(w)}) is less "
455  f"(than threshold {self.config.minNumberGoodPixelsForCovariance})")
456  return np.nan, np.nan, None
457 
458  maxRangeCov = self.config.maximumRangeCovariancesAstier
459  if covAstierRealSpace:
460  # Calculate covariances in real space.
461  covDiffAstier = computeCovDirect(diffIm.image.array, w, maxRangeCov)
462  else:
463  # Calculate covariances via FFT (default).
464  shapeDiff = np.array(diffIm.image.array.shape)
465  # Calculate the sizes of FFT dimensions.
466  s = shapeDiff + maxRangeCov
467  tempSize = np.array(np.log(s)/np.log(2.)).astype(int)
468  fftSize = np.array(2**(tempSize+1)).astype(int)
469  fftShape = (fftSize[0], fftSize[1])
470 
471  c = CovFastFourierTransform(diffIm.image.array, w, fftShape, maxRangeCov)
472  covDiffAstier = c.reportCovFastFourierTransform(maxRangeCov)
473 
474  # Compare Cov[0,0] and afwMath.VARIANCECLIP
475  # covDiffAstier[0] is the Cov[0,0] element, [3] is the variance, and there's a factor of 0.5
476  # difference with afwMath.VARIANCECLIP.
477  thresholdPercentage = self.config.thresholdDiffAfwVarVsCov00
478  fractionalDiff = 100*np.fabs(1 - varDiff/(covDiffAstier[0][3]*0.5))
479  if fractionalDiff >= thresholdPercentage:
480  self.log.warn("Absolute fractional difference between afwMatch.VARIANCECLIP and Cov[0,0] "
481  f"is more than {thresholdPercentage}%: {fractionalDiff}")
482 
483  return mu, varDiff, covDiffAstier
def measureMeanVarCov(self, exposure1, exposure2, region=None, covAstierRealSpace=False)
def makeCovArray(inputTuple, maxRangeFromTuple=8)
def computeCovDirect(diffImage, weightImage, maxRange)
def sigmaClipCorrection(nSigClip)
Definition: utils.py:42
def arrangeFlatsByExpId(exposureList)
Definition: utils.py:571
def arrangeFlatsByExpTime(exposureList)
Definition: utils.py:546