LSST Applications  21.0.0-172-gfb10e10a+18fedfabac,22.0.0+297cba6710,22.0.0+80564b0ff1,22.0.0+8d77f4f51a,22.0.0+a28f4c53b1,22.0.0+dcf3732eb2,22.0.1-1-g7d6de66+2a20fdde0d,22.0.1-1-g8e32f31+297cba6710,22.0.1-1-geca5380+7fa3b7d9b6,22.0.1-12-g44dc1dc+2a20fdde0d,22.0.1-15-g6a90155+515f58c32b,22.0.1-16-g9282f48+790f5f2caa,22.0.1-2-g92698f7+dcf3732eb2,22.0.1-2-ga9b0f51+7fa3b7d9b6,22.0.1-2-gd1925c9+bf4f0e694f,22.0.1-24-g1ad7a390+a9625a72a8,22.0.1-25-g5bf6245+3ad8ecd50b,22.0.1-25-gb120d7b+8b5510f75f,22.0.1-27-g97737f7+2a20fdde0d,22.0.1-32-gf62ce7b1+aa4237961e,22.0.1-4-g0b3f228+2a20fdde0d,22.0.1-4-g243d05b+871c1b8305,22.0.1-4-g3a563be+32dcf1063f,22.0.1-4-g44f2e3d+9e4ab0f4fa,22.0.1-42-gca6935d93+ba5e5ca3eb,22.0.1-5-g15c806e+85460ae5f3,22.0.1-5-g58711c4+611d128589,22.0.1-5-g75bb458+99c117b92f,22.0.1-6-g1c63a23+7fa3b7d9b6,22.0.1-6-g50866e6+84ff5a128b,22.0.1-6-g8d3140d+720564cf76,22.0.1-6-gd805d02+cc5644f571,22.0.1-8-ge5750ce+85460ae5f3,master-g6e05de7fdc+babf819c66,master-g99da0e417a+8d77f4f51a,w.2021.48
LSST Data Management Base Package
crosstalk.py
Go to the documentation of this file.
1 #
2 # This file is part of obs_decam.
3 #
4 # Developed for the LSST Data Management System.
5 # This product includes software developed by the LSST Project
6 # (http://www.lsst.org).
7 # Copyright 2018 LSST Corporation.
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 #
22 
23 """Apply crosstalk corrections to DECam images.
24 
25 Some of this code is based on DECam_crosstalk.py, used for the THELI pipeline
26 and written by Thomas Erben (private communication, 2018)
27 """
28 
29 
30 from collections import defaultdict
31 from astropy.io import fits
32 import datetime as dt
33 import logging
34 import os
35 import argparse
36 import re
37 import numpy as np
38 
39 import lsst.pipe.base as pipeBase
40 from lsst.ip.isr import CrosstalkConfig, CrosstalkTask, IsrTask
41 from lsst.obs.decam import DecamMapper
42 from lsst.utils import getPackageDir
43 
44 
45 __all__ = ["DecamCrosstalkTask"]
46 
47 
49  """Configuration for DECam crosstalk removal.
50  """
51 
52  def getSourcesAndCoeffsFile(self, filename='DECam_xtalk_20130606.txt'):
53  """File containing DECam crosstalk coefficients.
54 
55  This text file is provided by NOAO in a particular format with
56  information about the DECam crosstalk coefficients. It is available at
57  http://www.ctio.noao.edu/noao/content/DECam-Calibration-Files
58 
59  Parameters
60  ----------
61  filename : `str`, optional
62  File containing the decam crosstalk coefficients, from NOAO.
63 
64  Returns
65  -------
66  result : `str`
67  Full path to filename.
68  """
69  mapper = DecamMapper()
70  packageName = mapper.getPackageName()
71  packageDir = getPackageDir(packageName)
72  return os.path.join(packageDir, 'decam', filename)
73 
75  """Read crosstalk sources and coefficients from DECam-specific file.
76 
77  Returns
78  -------
79  sources : `defaultdict`
80  Sources of crosstalk read from crosstalk_file.
81  coeffs : `dict`
82  Linear crosstalk coefficients corresponding to sources.
83  """
84  crosstalk_file = self.getSourcesAndCoeffsFilegetSourcesAndCoeffsFile()
85  sources = defaultdict(list)
86  coeffs = {}
87  log = logging.getLogger('lsst.obs.decam.DecamCrosstalkConfig')
88  log.info('Reading crosstalk coefficient data')
89  with open(crosstalk_file) as f:
90  for line in f:
91  li = line.strip()
92  if not li.startswith('#'):
93  elem = li.split()
94  # The xtalk file has image areas like 'ccd01A'; preserve only '01A'
95  sources[elem[0][3:]].append(elem[1][3:])
96  # The linear crosstalk coefficients
97  coeffs[(elem[0][3:], elem[1][3:])] = float(elem[2])
98  return sources, coeffs
99 
100 
102  """Remove crosstalk from a raw DECam image.
103 
104  This must be done after the overscan is corrected but before it is trimmed.
105  This version of crosstalk removal assumes you want to do it as part of a
106  regular stack-ISR workflow.
107  """
108  ConfigClass = DecamCrosstalkConfig
109  _DefaultName = 'decamCrosstalk'
110 
111  def prepCrosstalk(self, dataRef, crosstalk=None):
112  """Retrieve source image data and coefficients to prepare for crosstalk correction.
113 
114  This is called by readIsrData. The crosstalkSources land in isrData.
115 
116  Parameters
117  ----------
118  dataRef : `lsst.daf.persistence.butlerSubset.ButlerDataRef`
119  DataRef of exposure to correct which must include a dataId with
120  at least one visit and a ccdnum corresponding to the detector id.
121  crosstalk : `lsst.ip.isr.CrosstalkCalib`
122  Crosstalk calibration that will be used.
123 
124  Returns
125  -------
126  crosstalkSources : `defaultdict`
127  Contains image data and corresponding crosstalk coefficients for
128  each crosstalk source of the given dataRef victim.
129  """
130  # Define a butler to prepare for accessing crosstalk 'source' image data
131  try:
132  butler = dataRef.getButler()
133  except RuntimeError:
134  self.log.fatal('Cannot get a Butler from the dataRef provided')
135 
136  # If there is no interChip CT, then we're done.
137  if crosstalk.interChip is None or len(crosstalk.interChip) == 0:
138  return None
139 
140  # Retrieve visit and ccdnum from 'victim' so we can look up 'sources'
141  visit = dataRef.dataId['visit']
142  crosstalkSources = defaultdict(list)
143 
144  decamisr = IsrTask() # needed for source_exposure overscan correction
145  for sourceDet in crosstalk.interChip:
146  try:
147  source_exposure = butler.get('raw', dataId={'visit': visit, 'ccd': sourceDet})
148  except RuntimeError:
149  self.log.warning("Cannot access source %s, SKIPPING IT", sourceDet)
150  else:
151  self.log.info("Correcting victim %s from crosstalk source %s.",
152  crosstalk._detectorName, sourceDet)
153 
154  for amp in source_exposure.getDetector():
155  decamisr.overscanCorrection(source_exposure, amp)
156 
157  # Do not assemble here, as DECam CT is corrected before assembly.
158  crosstalkSources[sourceDet] = source_exposure
159  return crosstalkSources
160 
161 
162 class DecamCrosstalkIO(pipeBase.Task):
163  """Remove crosstalk from a raw DECam Multi-Extension FITS file.
164 
165  This must be done after the overscan is corrected but before it is trimmed.
166  NOTE: This version of crosstalk removal assumes you want to do it as an
167  standalone file operation, not as part of a regular stack-ISR workflow.
168  """
169  ConfigClass = DecamCrosstalkConfig
170  _DefaultName = 'decamCrosstalkIO'
171 
173  """Parse command-line arguments to run the crosstalk correction alone.
174 
175  Examples
176  --------
177  runCrosstalkAlone.py -f decam0411673.fits.fz -o decam0411673_xtalk_corrected.fits.fz
178  """
179  parser = argparse.ArgumentParser()
180  parser.add_argument('-f', '--filename',
181  help='Location on disk of the DECam MEF file you wish to crosstalk correct')
182  parser.add_argument('-o', '--outfile',
183  help='Location on disk where crosstalk-corrected DECam MEF file will be saved')
184  args = parser.parse_args()
185  parsed = {'filename': args.filename, 'outfile': args.outfile}
186  return parsed
187 
188  def runCrosstalkAlone(self):
189  """Utility for crosstalk correction directly on a DECam MEF file.
190  """
191  # Override the logger attached to this Task.
192  log = logging.getLogger('lsst.obs.decam.DecamCrosstalkIO')
193  parsed = self.parseCrosstalkIOArgsparseCrosstalkIOArgs()
194  filename = parsed['filename']
195  outfile = parsed['outfile']
196  mef = fits.open(filename)
197  sources, coeffs = self.config.getSourcesAndCoeffs()
198  corrected = subtractCrosstalkIO(mef, sources, coeffs)
199 
200  if outfile:
201  corrected.writeto(outfile, overwrite=True)
202  else:
203  log.warning('No outfile specified, not writing corrected image to disk')
204 
205 
206 def subtractCrosstalkIO(mef, sources, coeffs):
207  """Remove crosstalk from a DECam Multi-Extension FITS raw image.
208 
209  This version uses astropy I/O rather than the Butler due to the cross-amp
210  crosstalk and ease of accessing Multi-Extension FITS data. It does not
211  integrate with isr and does not use the Butler.
212 
213  WARNING: this will run fine on a non-overscan-corrected raw image, but the
214  crosstalk coefficients were computed post-overscan-correction. Therefore,
215  **the input raw image should be overscan corrected before this is run.**
216 
217  Parameters
218  ----------
219  mef : `astropy.io.fits.hdu.hdulist.HDUList`
220  One DECam Multi-Extension FITS image.
221  sources : `defaultdict`
222  Crosstalk source areas affecting flux in a certain victim area.
223  coeffs : `list`
224  Linear crosstalk coefficients.
225 
226  Returns
227  -------
228  mef : `astropy.io.fits.hdu.hdulist.HDUList`
229  New MEF image with crosstalk corrected data and updated header.
230  """
231  log = logging.getLogger('lsst.obs.decam.subtractCrosstalkIO')
232  lowx = {} # Lower x pixel index of given chip region
233  highx = {} # Upper x pixel index of given chip region
234  extension = {} # FITS extension of given chip position from CCDNUM header
235  corrected = [0] * 63 # placeholder for corrected data
236 
237  for i in range(1, 63):
238  ccd = '%02d' % mef[i].header['CCDNUM']
239  extension[ccd] = i
240  for section in 'AB':
241  ccdsec_str = mef[i].header['DATASEC' + section]
242  ccdsec_list = re.findall(r"[\w']+", ccdsec_str)
243  lowx[ccd + section] = int(ccdsec_list[0]) - 1
244  highx[ccd + section] = int(ccdsec_list[1])
245  corrected[i] = mef[i].data.astype(np.float32)
246  for ccd in ['%02d' % i for i in range(1, 63)]:
247  for section in 'AB':
248  victim = ccd + section
249  victim_data = corrected[extension[ccd]][:, lowx[victim]:highx[victim]]
250  for source in sources[victim]:
251  source_data = mef[extension[source[:2]]].data[:, lowx[source]:highx[source]]
252  if lowx[victim] != lowx[source]:
253  # If regions were read out on different sides of the CCD
254  # (i.e., an A and B mismatch) flip data in the x-direction
255  source_data = np.fliplr(source_data)
256  # Perform linear crosstalk correction
257  victim_data[:, :] = victim_data - coeffs[(victim, source)] * source_data
258  log.info('Correcting victim %s from crosstalk source %s for HDU %s', victim, source, ccd)
259 
260  for i in range(1, 63):
261  mef[i].header['HISTORY'] = 'Crosstalk corrected on {0}'.format(
262  dt.datetime.now().strftime("%Y-%m-%dT%H:%M:%S"))
263  mef[i] = fits.ImageHDU(header=mef[i].header, data=corrected[i])
264 
265  return mef
def getSourcesAndCoeffsFile(self, filename='DECam_xtalk_20130606.txt')
Definition: crosstalk.py:52
def prepCrosstalk(self, dataRef, crosstalk=None)
Definition: crosstalk.py:111
std::shared_ptr< FrameSet > append(FrameSet const &first, FrameSet const &second)
Construct a FrameSet that performs two transformations in series.
Definition: functional.cc:33
std::string getPackageDir(std::string const &packageName)
return the root directory of a setup package
Definition: packaging.cc:33
def subtractCrosstalkIO(mef, sources, coeffs)
Definition: crosstalk.py:206
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174