lsst.cp.pipe  20.0.0-23-g10eeb28+96244339af
cpFlatNormTask.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 <http://www.gnu.org/licenses/>.
21 import numpy as np
22 from collections import defaultdict
23 
24 import lsst.afw.math as afwMath
25 import lsst.daf.base as dafBase
26 import lsst.pex.config as pexConfig
27 import lsst.pipe.base as pipeBase
29 from lsst.cp.pipe.cpCombine import VignetteExposure
30 
31 from ._lookupStaticCalibration import lookupStaticCalibration
32 
33 __all__ = ["CpFlatMeasureTask", "CpFlatMeasureTaskConfig",
34  "CpFlatNormalizationTask", "CpFlatNormalizationTaskConfig"]
35 
36 
37 class CpFlatMeasureConnections(pipeBase.PipelineTaskConnections,
38  dimensions=("instrument", "exposure", "detector")):
39  inputExp = cT.Input(
40  name="postISRCCD",
41  doc="Input exposure to measure statistics from.",
42  storageClass="Exposure",
43  dimensions=("instrument", "exposure", "detector"),
44  )
45  outputStats = cT.Output(
46  name="flatStats",
47  doc="Output statistics to write.",
48  storageClass="PropertyList",
49  dimensions=("instrument", "exposure", "detector"),
50  )
51 
52 
53 class CpFlatMeasureTaskConfig(pipeBase.PipelineTaskConfig,
54  pipelineConnections=CpFlatMeasureConnections):
55  maskNameList = pexConfig.ListField(
56  dtype=str,
57  doc="Mask list to exclude from statistics calculations.",
58  default=['DETECTED', 'BAD', 'NO_DATA'],
59  )
60  doVignette = pexConfig.Field(
61  dtype=bool,
62  doc="Mask vignetted regions?",
63  default=True,
64  )
65  numSigmaClip = pexConfig.Field(
66  dtype=float,
67  doc="Rejection threshold (sigma) for statistics clipping.",
68  default=3.0,
69  )
70  clipMaxIter = pexConfig.Field(
71  dtype=int,
72  doc="Max number of clipping iterations to apply.",
73  default=3,
74  )
75 
76 
77 class CpFlatMeasureTask(pipeBase.PipelineTask,
78  pipeBase.CmdLineTask):
79  """Apply extra masking and measure image statistics.
80  """
81  ConfigClass = CpFlatMeasureTaskConfig
82  _DefaultName = "cpFlatMeasure"
83 
84  def run(self, inputExp):
85  """Mask ISR processed FLAT exposures to ensure consistent statistics.
86 
87  Parameters
88  ----------
89  inputExp : `lsst.afw.image.Exposure`
90  Post-ISR processed exposure to measure.
91 
92  Returns
93  -------
94  outputStats : `lsst.daf.base.PropertyList`
95  List containing the statistics.
96  """
97  if self.config.doVignette:
98  VignetteExposure(inputExp, doUpdateMask=True, maskPlane='BAD',
99  doSetValue=False, log=self.log)
100  mask = inputExp.getMask()
101  maskVal = mask.getPlaneBitMask(self.config.maskNameList)
102  statsControl = afwMath.StatisticsControl(self.config.numSigmaClip,
103  self.config.clipMaxIter,
104  maskVal)
105  statsControl.setAndMask(maskVal)
106 
107  outputStats = dafBase.PropertyList()
108 
109  # Detector level:
110  stats = afwMath.makeStatistics(inputExp.getMaskedImage(),
111  afwMath.MEANCLIP | afwMath.STDEVCLIP | afwMath.NPOINT,
112  statsControl)
113  outputStats['DETECTOR_MEDIAN'] = stats.getValue(afwMath.MEANCLIP)
114  outputStats['DETECTOR_SIGMA'] = stats.getValue(afwMath.STDEVCLIP)
115  outputStats['DETECTOR_N'] = stats.getValue(afwMath.NPOINT)
116  self.log.info("Stats: median=%f sigma=%f n=%d",
117  outputStats['DETECTOR_MEDIAN'],
118  outputStats['DETECTOR_SIGMA'],
119  outputStats['DETECTOR_N'])
120 
121  # AMP LEVEL:
122  for ampIdx, amp in enumerate(inputExp.getDetector()):
123  ampName = amp.getName()
124  ampExp = inputExp.Factory(inputExp, amp.getBBox())
125  stats = afwMath.makeStatistics(ampExp.getMaskedImage(),
126  afwMath.MEANCLIP | afwMath.STDEVCLIP | afwMath.NPOINT,
127  statsControl)
128  outputStats[f'AMP_NAME_{ampIdx}'] = ampName
129  outputStats[f'AMP_MEDIAN_{ampIdx}'] = stats.getValue(afwMath.MEANCLIP)
130  outputStats[f'AMP_SIGMA_{ampIdx}'] = stats.getValue(afwMath.STDEVCLIP)
131  outputStats[f'AMP_N_{ampIdx}'] = stats.getValue(afwMath.NPOINT)
132 
133  return pipeBase.Struct(
134  outputStats=outputStats
135  )
136 
137 
138 class CpFlatNormalizationConnections(pipeBase.PipelineTaskConnections,
139  dimensions=("instrument", "physical_filter")):
140  inputMDs = cT.Input(
141  name="cpFlatProc_metadata",
142  doc="Input metadata for each visit/detector in input set.",
143  storageClass="PropertyList",
144  dimensions=("instrument", "physical_filter", "detector", "exposure"),
145  multiple=True,
146  )
147  camera = cT.PrerequisiteInput(
148  name="camera",
149  doc="Input camera to use for gain lookup.",
150  storageClass="Camera",
151  dimensions=("instrument",),
152  lookupFunction=lookupStaticCalibration,
153  isCalibration=True,
154  )
155 
156  outputScales = cT.Output(
157  name="cpFlatNormScales",
158  doc="Output combined proposed calibration.",
159  storageClass="StructuredDataDict",
160  dimensions=("instrument", "physical_filter"),
161  )
162 
163 
164 class CpFlatNormalizationTaskConfig(pipeBase.PipelineTaskConfig,
165  pipelineConnections=CpFlatNormalizationConnections):
166  level = pexConfig.ChoiceField(
167  dtype=str,
168  doc="Which level to apply normalizations.",
169  default='DETECTOR',
170  allowed={
171  'DETECTOR': "Correct using full detector statistics.",
172  'AMP': "Correct using individual amplifiers.",
173  },
174  )
175  scaleMaxIter = pexConfig.Field(
176  dtype=int,
177  doc="Max number of iterations to use in scale solver.",
178  default=10,
179  )
180 
181 
182 class CpFlatNormalizationTask(pipeBase.PipelineTask,
183  pipeBase.CmdLineTask):
184  """Rescale merged flat frames to remove unequal screen illumination.
185  """
186  ConfigClass = CpFlatNormalizationTaskConfig
187  _DefaultName = "cpFlatNorm"
188 
189  def runQuantum(self, butlerQC, inputRefs, outputRefs):
190  inputs = butlerQC.get(inputRefs)
191 
192  # Use the dimensions of the inputs for generating
193  # output scales.
194  dimensions = [exp.dataId.byName() for exp in inputRefs.inputMDs]
195  inputs['inputDims'] = dimensions
196 
197  outputs = self.run(**inputs)
198  butlerQC.put(outputs, outputRefs)
199 
200  def run(self, inputMDs, inputDims, camera):
201  """Normalize FLAT exposures to a consistent level.
202 
203  Parameters
204  ----------
205  inputMDs : `list` [`lsst.daf.base.PropertyList`]
206  Amplifier-level metadata used to construct scales.
207  inputDims : `list` [`dict`]
208  List of dictionaries of input data dimensions/values.
209  Each list entry should contain:
210 
211  ``"exposure"``
212  exposure id value (`int`)
213  ``"detector"``
214  detector id value (`int`)
215 
216  Returns
217  -------
218  outputScales : `dict` [`dict` [`dict` [`float`]]]
219  Dictionary of scales, indexed by detector (`int`),
220  amplifier (`int`), and exposure (`int`).
221 
222  Raises
223  ------
224  KeyError
225  Raised if the input dimensions do not contain detector and
226  exposure, or if the metadata does not contain the expected
227  statistic entry.
228  """
229  expSet = sorted(set([d['exposure'] for d in inputDims]))
230  detSet = sorted(set([d['detector'] for d in inputDims]))
231 
232  expMap = {exposureId: idx for idx, exposureId in enumerate(expSet)}
233  detMap = {detectorId: idx for idx, detectorId in enumerate(detSet)}
234 
235  nExp = len(expSet)
236  nDet = len(detSet)
237  if self.config.level == 'DETECTOR':
238  bgMatrix = np.zeros((nDet, nExp))
239  bgCounts = np.ones((nDet, nExp))
240  elif self.config.level == 'AMP':
241  nAmp = len(camera[0])
242  bgMatrix = np.zeros((nDet * nAmp, nExp))
243  bgCounts = np.ones((nDet * nAmp, nExp))
244 
245  for inMetadata, inDimensions in zip(inputMDs, inputDims):
246  try:
247  exposureId = inDimensions['exposure']
248  detectorId = inDimensions['detector']
249  except Exception as e:
250  raise KeyError("Cannot find expected dimensions in %s" % (inDimensions, )) from e
251 
252  if self.config.level == 'DETECTOR':
253  detId = detMap[detectorId]
254  expId = expMap[exposureId]
255  try:
256  value = inMetadata.get('DETECTOR_MEDIAN')
257  count = inMetadata.get('DETECTOR_N')
258  except Exception as e:
259  raise KeyError("Cannot read expected metadata string.") from e
260 
261  if np.isfinite(value):
262  bgMatrix[detId][expId] = value
263  bgCounts[detId][expId] = count
264  else:
265  bgMatrix[detId][expId] = np.nan
266  bgCounts[detId][expId] = 1
267 
268  elif self.config.level == 'AMP':
269  detector = camera[detectorId]
270  nAmp = len(detector)
271 
272  detId = detMap[detectorId] * nAmp
273  expId = expMap[exposureId]
274 
275  for ampIdx, amp in enumerate(detector):
276  try:
277  value = inMetadata.get(f'AMP_MEDIAN_{ampIdx}')
278  count = inMetadata.get(f'AMP_N_{ampIdx}')
279  except Exception as e:
280  raise KeyError("cannot read expected metadata string.") from e
281 
282  detAmpId = detId + ampIdx
283  if np.isfinite(value):
284  bgMatrix[detAmpId][expId] = value
285  bgCounts[detAmpId][expId] = count
286  else:
287  bgMatrix[detAmpId][expId] = np.nan
288  bgMatrix[detAmpId][expId] = 1
289 
290  scaleResult = self.measureScales(bgMatrix, bgCounts, iterations=self.config.scaleMaxIter)
291  expScales = scaleResult.expScales
292  detScales = scaleResult.detScales
293 
294  outputScales = defaultdict(lambda: defaultdict(lambda: defaultdict(lambda: defaultdict(float))))
295 
296  if self.config.level == 'DETECTOR':
297  for detId, det in enumerate(detSet):
298  for amp in camera[detId]:
299  for expId, exp in enumerate(expSet):
300  outputScales['expScale'][det][amp.getName()][exp] = expScales[expId].tolist()
301  outputScales['detScale'][det] = detScales[detId].tolist()
302  elif self.config.level == 'AMP':
303  for detId, det in enumerate(detSet):
304  for ampIdx, amp in enumerate(camera[detId]):
305  for expId, exp in enumerate(expSet):
306  outputScales['expScale'][det][amp.getName()][exp] = expScales[expId].tolist()
307  detAmpId = detId + ampIdx
308  outputScales['detScale'][det][amp.getName()] = detScales[detAmpId].tolist()
309 
310  return pipeBase.Struct(
311  outputScales=outputScales,
312  )
313 
314  def measureScales(self, bgMatrix, bgCounts=None, iterations=10):
315  """Convert backgrounds to exposure and detector components.
316 
317  Parameters
318  ----------
319  bgMatrix : `np.ndarray`, (nDetectors, nExposures)
320  Input backgrounds indexed by exposure (axis=0) and
321  detector (axis=1).
322  bgCounts : `np.ndarray`, (nDetectors, nExposures), optional
323  Input pixel counts used to in measuring bgMatrix, indexed
324  identically.
325  iterations : `int`, optional
326  Number of iterations to use in decomposition.
327 
328  Returns
329  -------
330  scaleResult : `lsst.pipe.base.Struct`
331  Result struct containing fields:
332 
333  ``vectorE``
334  Output E vector of exposure level scalings
335  (`np.array`, (nExposures)).
336  ``vectorG``
337  Output G vector of detector level scalings
338  (`np.array`, (nExposures)).
339  ``bgModel``
340  Expected model bgMatrix values, calculated from E and G
341  (`np.ndarray`, (nDetectors, nExposures)).
342 
343  Notes
344  -----
345 
346  The set of background measurements B[exposure, detector] of
347  flat frame data should be defined by a "Cartesian" product of
348  two vectors, E[exposure] and G[detector]. The E vector
349  represents the total flux incident on the focal plane. In a
350  perfect camera, this is simply the sum along the columns of B
351  (np.sum(B, axis=0)).
352 
353  However, this simple model ignores differences in detector
354  gains, the vignetting of the detectors, and the illumination
355  pattern of the source lamp. The G vector describes these
356  detector dependent differences, which should be identical over
357  different exposures. For a perfect lamp of unit total
358  intensity, this is simply the sum along the rows of B
359  (np.sum(B, axis=1)). This algorithm divides G by the total
360  flux level, to provide the relative (not absolute) scales
361  between detectors.
362 
363  The algorithm here, from pipe_drivers/constructCalibs.py and
364  from there from Eugene Magnier/PanSTARRS [1]_, attempts to
365  iteratively solve this decomposition from initial "perfect" E
366  and G vectors. The operation is performed in log space to
367  reduce the multiply and divides to linear additions and
368  subtractions.
369 
370  References
371  ----------
372  .. [1] https://svn.pan-starrs.ifa.hawaii.edu/trac/ipp/browser/trunk/psModules/src/detrend/pmFlatNormalize.c # noqa: E501
373 
374  """
375  numExps = bgMatrix.shape[1]
376  numChips = bgMatrix.shape[0]
377  if bgCounts is None:
378  bgCounts = np.ones_like(bgMatrix)
379 
380  logMeas = np.log(bgMatrix)
381  logMeas = np.ma.masked_array(logMeas, ~np.isfinite(logMeas))
382  logG = np.zeros(numChips)
383  logE = np.array([np.average(logMeas[:, iexp] - logG,
384  weights=bgCounts[:, iexp]) for iexp in range(numExps)])
385 
386  for iter in range(iterations):
387  logG = np.array([np.average(logMeas[ichip, :] - logE,
388  weights=bgCounts[ichip, :]) for ichip in range(numChips)])
389 
390  bad = np.isnan(logG)
391  if np.any(bad):
392  logG[bad] = logG[~bad].mean()
393 
394  logE = np.array([np.average(logMeas[:, iexp] - logG,
395  weights=bgCounts[:, iexp]) for iexp in range(numExps)])
396  fluxLevel = np.average(np.exp(logG), weights=np.sum(bgCounts, axis=1))
397 
398  logG -= np.log(fluxLevel)
399  self.log.debug(f"ITER {iter}: Flux: {fluxLevel}")
400  self.log.debug(f"Exps: {np.exp(logE)}")
401  self.log.debug(f"{np.mean(logG)}")
402 
403  logE = np.array([np.average(logMeas[:, iexp] - logG,
404  weights=bgCounts[:, iexp]) for iexp in range(numExps)])
405 
406  bgModel = np.exp(logE[np.newaxis, :] - logG[:, np.newaxis])
407  return pipeBase.Struct(
408  expScales=np.exp(logE),
409  detScales=np.exp(logG),
410  bgModel=bgModel,
411  )
lsst.cp.pipe.cpFlatNormTask.CpFlatMeasureTask.run
def run(self, inputExp)
Definition: cpFlatNormTask.py:84
lsst.cp.pipe.cpCombine.VignetteExposure
def VignetteExposure(exposure, polygon=None, doUpdateMask=True, maskPlane='BAD', doSetValue=False, vignetteValue=0.0, log=None)
Definition: cpCombine.py:528
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationTask.runQuantum
def runQuantum(self, butlerQC, inputRefs, outputRefs)
Definition: cpFlatNormTask.py:189
lsst.cp.pipe.cpFlatNormTask.CpFlatMeasureTaskConfig
Definition: cpFlatNormTask.py:54
lsst.cp.pipe.cpCombine
Definition: cpCombine.py:1
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationTask.run
def run(self, inputMDs, inputDims, camera)
Definition: cpFlatNormTask.py:200
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationConnections
Definition: cpFlatNormTask.py:139
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationTaskConfig
Definition: cpFlatNormTask.py:165
lsst.cp.pipe.cpFlatNormTask.CpFlatMeasureTask
Definition: cpFlatNormTask.py:78
lsst::pex::config
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationTask
Definition: cpFlatNormTask.py:183
lsst::daf::base
lsst::afw::math
lsst.cp.pipe.cpFlatNormTask.CpFlatNormalizationTask.measureScales
def measureScales(self, bgMatrix, bgCounts=None, iterations=10)
Definition: cpFlatNormTask.py:314
lsst::pipe::base
lsst::pipe::base::connectionTypes
lsst.cp.pipe.cpFlatNormTask.CpFlatMeasureConnections
Definition: cpFlatNormTask.py:38