LSSTApplications  17.0+124,17.0+14,17.0+73,18.0.0+37,18.0.0+80,18.0.0-4-g68ffd23+4,18.1.0-1-g0001055+12,18.1.0-1-g03d53ef+5,18.1.0-1-g1349e88+55,18.1.0-1-g2505f39+44,18.1.0-1-g5315e5e+4,18.1.0-1-g5e4b7ea+14,18.1.0-1-g7e8fceb+4,18.1.0-1-g85f8cd4+48,18.1.0-1-g8ff0b9f+4,18.1.0-1-ga2c679d+1,18.1.0-1-gd55f500+35,18.1.0-10-gb58edde+2,18.1.0-11-g0997b02+4,18.1.0-13-gfe4edf0b+12,18.1.0-14-g259bd21+21,18.1.0-19-gdb69f3f+2,18.1.0-2-g5f9922c+24,18.1.0-2-gd3b74e5+11,18.1.0-2-gfbf3545+32,18.1.0-26-g728bddb4+5,18.1.0-27-g6ff7ca9+2,18.1.0-3-g52aa583+25,18.1.0-3-g8ea57af+9,18.1.0-3-gb69f684+42,18.1.0-3-gfcaddf3+6,18.1.0-32-gd8786685a,18.1.0-4-gf3f9b77+6,18.1.0-5-g1dd662b+2,18.1.0-5-g6dbcb01+41,18.1.0-6-gae77429+3,18.1.0-7-g9d75d83+9,18.1.0-7-gae09a6d+30,18.1.0-9-gc381ef5+4,w.2019.45
LSSTDataManagementBasePackage
Classes | Functions
lsst.obs.decam.crosstalk Namespace Reference

Classes

class  DecamCrosstalkConfig
 
class  DecamCrosstalkIO
 
class  DecamCrosstalkTask
 

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 288 of file crosstalk.py.

288 def subtractCrosstalkIO(mef, sources, coeffs):
289  """Remove crosstalk from a DECam Multi-Extension FITS raw image.
290 
291  This version uses astropy I/O rather than the Butler due to the cross-amp
292  crosstalk and ease of accessing Multi-Extension FITS data. It does not
293  integrate with isr and does not use the Butler.
294 
295  WARNING: this will run fine on a non-overscan-corrected raw image, but the
296  crosstalk coefficients were computed post-overscan-correction. Therefore,
297  **the input raw image should be overscan corrected before this is run.**
298 
299  Parameters
300  ----------
301  mef : `astropy.io.fits.hdu.hdulist.HDUList`
302  One DECam Multi-Extension FITS image.
303  sources : `defaultdict`
304  Crosstalk source areas affecting flux in a certain victim area.
305  coeffs : `list`
306  Linear crosstalk coefficients.
307 
308  Returns
309  -------
310  mef : `astropy.io.fits.hdu.hdulist.HDUList`
311  New MEF image with crosstalk corrected data and updated header.
312  """
313  log = lsst.log.Log.getLogger('obs.decam.subtractCrosstalkIO')
314  lowx = {} # Lower x pixel index of given chip region
315  highx = {} # Upper x pixel index of given chip region
316  extension = {} # FITS extension of given chip position from CCDNUM header
317  corrected = [0] * 63 # placeholder for corrected data
318 
319  for i in range(1, 63):
320  ccd = '%02d' % mef[i].header['CCDNUM']
321  extension[ccd] = i
322  for section in 'AB':
323  ccdsec_str = mef[i].header['DATASEC' + section]
324  ccdsec_list = re.findall(r"[\w']+", ccdsec_str)
325  lowx[ccd + section] = int(ccdsec_list[0]) - 1
326  highx[ccd + section] = int(ccdsec_list[1])
327  corrected[i] = mef[i].data.astype(np.float32)
328  for ccd in ['%02d' % i for i in range(1, 63)]:
329  for section in 'AB':
330  victim = ccd + section
331  victim_data = corrected[extension[ccd]][:, lowx[victim]:highx[victim]]
332  for source in sources[victim]:
333  source_data = mef[extension[source[:2]]].data[:, lowx[source]:highx[source]]
334  if lowx[victim] != lowx[source]:
335  # If regions were read out on different sides of the CCD
336  # (i.e., an A and B mismatch) flip data in the x-direction
337  source_data = np.fliplr(source_data)
338  # Perform linear crosstalk correction
339  victim_data[:, :] = victim_data - coeffs[(victim, source)] * source_data
340  log.info('Correcting victim %s from crosstalk source %s for HDU %s' % (victim, source, ccd))
341 
342  for i in range(1, 63):
343  mef[i].header['HISTORY'] = 'Crosstalk corrected on {0}'.format(
344  dt.datetime.now().strftime("%Y-%m-%dT%H:%M:%S"))
345  mef[i] = fits.ImageHDU(header=mef[i].header, data=corrected[i])
346 
347  return mef
348 
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174
def subtractCrosstalkIO(mef, sources, coeffs)
Definition: crosstalk.py:288
static Log getLogger(Log const &logger)
Definition: Log.h:760