LSST Applications  21.0.0-147-g0e635eb1+1acddb5be5,22.0.0+052faf71bd,22.0.0+1ea9a8b2b2,22.0.0+6312710a6c,22.0.0+729191ecac,22.0.0+7589c3a021,22.0.0+9f079a9461,22.0.1-1-g7d6de66+b8044ec9de,22.0.1-1-g87000a6+536b1ee016,22.0.1-1-g8e32f31+6312710a6c,22.0.1-10-gd060f87+016f7cdc03,22.0.1-12-g9c3108e+df145f6f68,22.0.1-16-g314fa6d+c825727ab8,22.0.1-19-g93a5c75+d23f2fb6d8,22.0.1-19-gb93eaa13+aab3ef7709,22.0.1-2-g8ef0a89+b8044ec9de,22.0.1-2-g92698f7+9f079a9461,22.0.1-2-ga9b0f51+052faf71bd,22.0.1-2-gac51dbf+052faf71bd,22.0.1-2-gb66926d+6312710a6c,22.0.1-2-gcb770ba+09e3807989,22.0.1-20-g32debb5+b8044ec9de,22.0.1-23-gc2439a9a+fb0756638e,22.0.1-3-g496fd5d+09117f784f,22.0.1-3-g59f966b+1e6ba2c031,22.0.1-3-g849a1b8+f8b568069f,22.0.1-3-gaaec9c0+c5c846a8b1,22.0.1-32-g5ddfab5d3+60ce4897b0,22.0.1-4-g037fbe1+64e601228d,22.0.1-4-g8623105+b8044ec9de,22.0.1-5-g096abc9+d18c45d440,22.0.1-5-g15c806e+57f5c03693,22.0.1-7-gba73697+57f5c03693,master-g6e05de7fdc+c1283a92b8,master-g72cdda8301+729191ecac,w.2021.39
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 206 of file crosstalk.py.

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('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 subtractCrosstalkIO(mef, sources, coeffs)
Definition: crosstalk.py:206
def format(config, name=None, writeSourceLine=True, prefix="", verbose=False)
Definition: history.py:174