LSST Applications  21.0.0+04719a4bac,21.0.0-1-ga51b5d4+f5e6047307,21.0.0-11-g2b59f77+a9c1acf22d,21.0.0-11-ga42c5b2+86977b0b17,21.0.0-12-gf4ce030+76814010d2,21.0.0-13-g1721dae+760e7a6536,21.0.0-13-g3a573fe+768d78a30a,21.0.0-15-g5a7caf0+f21cbc5713,21.0.0-16-g0fb55c1+b60e2d390c,21.0.0-19-g4cded4ca+71a93a33c0,21.0.0-2-g103fe59+bb20972958,21.0.0-2-g45278ab+04719a4bac,21.0.0-2-g5242d73+3ad5d60fb1,21.0.0-2-g7f82c8f+8babb168e8,21.0.0-2-g8f08a60+06509c8b61,21.0.0-2-g8faa9b5+616205b9df,21.0.0-2-ga326454+8babb168e8,21.0.0-2-gde069b7+5e4aea9c2f,21.0.0-2-gecfae73+1d3a86e577,21.0.0-2-gfc62afb+3ad5d60fb1,21.0.0-25-g1d57be3cd+e73869a214,21.0.0-3-g357aad2+ed88757d29,21.0.0-3-g4a4ce7f+3ad5d60fb1,21.0.0-3-g4be5c26+3ad5d60fb1,21.0.0-3-g65f322c+e0b24896a3,21.0.0-3-g7d9da8d+616205b9df,21.0.0-3-ge02ed75+a9c1acf22d,21.0.0-4-g591bb35+a9c1acf22d,21.0.0-4-g65b4814+b60e2d390c,21.0.0-4-gccdca77+0de219a2bc,21.0.0-4-ge8a399c+6c55c39e83,21.0.0-5-gd00fb1e+05fce91b99,21.0.0-6-gc675373+3ad5d60fb1,21.0.0-64-g1122c245+4fb2b8f86e,21.0.0-7-g04766d7+cd19d05db2,21.0.0-7-gdf92d54+04719a4bac,21.0.0-8-g5674e7b+d1bd76f71f,master-gac4afde19b+a9c1acf22d,w.2021.13
LSST Data Management Base Package
Classes | Functions
lsst.obs.decam.crosstalk Namespace Reference

Classes

class  DecamCrosstalkConfig
 
class  DecamCrosstalkTask
 
class  DecamCrosstalkIO
 

Functions

def subtractCrosstalkIO (mef, sources, coeffs)
 

Function Documentation

◆ subtractCrosstalkIO()

def lsst.obs.decam.crosstalk.subtractCrosstalkIO (   mef,
  sources,
  coeffs 
)
Remove crosstalk from a DECam Multi-Extension FITS raw image.

This version uses astropy I/O rather than the Butler due to the cross-amp
crosstalk and ease of accessing Multi-Extension FITS data. It does not
integrate with isr and does not use the Butler.

WARNING: this will run fine on a non-overscan-corrected raw image, but the
crosstalk coefficients were computed post-overscan-correction. Therefore,
**the input raw image should be overscan corrected before this is run.**

Parameters
----------
mef : `astropy.io.fits.hdu.hdulist.HDUList`
    One DECam Multi-Extension FITS image.
sources : `defaultdict`
    Crosstalk source areas affecting flux in a certain victim area.
coeffs : `list`
    Linear crosstalk coefficients.

Returns
-------
mef : `astropy.io.fits.hdu.hdulist.HDUList`
    New MEF image with crosstalk corrected data and updated header.

Definition at line 205 of file crosstalk.py.

205 def subtractCrosstalkIO(mef, sources, coeffs):
206  """Remove crosstalk from a DECam Multi-Extension FITS raw image.
207 
208  This version uses astropy I/O rather than the Butler due to the cross-amp
209  crosstalk and ease of accessing Multi-Extension FITS data. It does not
210  integrate with isr and does not use the Butler.
211 
212  WARNING: this will run fine on a non-overscan-corrected raw image, but the
213  crosstalk coefficients were computed post-overscan-correction. Therefore,
214  **the input raw image should be overscan corrected before this is run.**
215 
216  Parameters
217  ----------
218  mef : `astropy.io.fits.hdu.hdulist.HDUList`
219  One DECam Multi-Extension FITS image.
220  sources : `defaultdict`
221  Crosstalk source areas affecting flux in a certain victim area.
222  coeffs : `list`
223  Linear crosstalk coefficients.
224 
225  Returns
226  -------
227  mef : `astropy.io.fits.hdu.hdulist.HDUList`
228  New MEF image with crosstalk corrected data and updated header.
229  """
230  log = lsst.log.Log.getLogger('obs.decam.subtractCrosstalkIO')
231  lowx = {} # Lower x pixel index of given chip region
232  highx = {} # Upper x pixel index of given chip region
233  extension = {} # FITS extension of given chip position from CCDNUM header
234  corrected = [0] * 63 # placeholder for corrected data
235 
236  for i in range(1, 63):
237  ccd = '%02d' % mef[i].header['CCDNUM']
238  extension[ccd] = i
239  for section in 'AB':
240  ccdsec_str = mef[i].header['DATASEC' + section]
241  ccdsec_list = re.findall(r"[\w']+", ccdsec_str)
242  lowx[ccd + section] = int(ccdsec_list[0]) - 1
243  highx[ccd + section] = int(ccdsec_list[1])
244  corrected[i] = mef[i].data.astype(np.float32)
245  for ccd in ['%02d' % i for i in range(1, 63)]:
246  for section in 'AB':
247  victim = ccd + section
248  victim_data = corrected[extension[ccd]][:, lowx[victim]:highx[victim]]
249  for source in sources[victim]:
250  source_data = mef[extension[source[:2]]].data[:, lowx[source]:highx[source]]
251  if lowx[victim] != lowx[source]:
252  # If regions were read out on different sides of the CCD
253  # (i.e., an A and B mismatch) flip data in the x-direction
254  source_data = np.fliplr(source_data)
255  # Perform linear crosstalk correction
256  victim_data[:, :] = victim_data - coeffs[(victim, source)] * source_data
257  log.info('Correcting victim %s from crosstalk source %s for HDU %s' % (victim, source, ccd))
258 
259  for i in range(1, 63):
260  mef[i].header['HISTORY'] = 'Crosstalk corrected on {0}'.format(
261  dt.datetime.now().strftime("%Y-%m-%dT%H:%M:%S"))
262  mef[i] = fits.ImageHDU(header=mef[i].header, data=corrected[i])
263 
264  return mef
static Log getLogger(Log const &logger)
Definition: Log.h:760
def subtractCrosstalkIO(mef, sources, coeffs)
Definition: crosstalk.py:205
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174