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